Skip to content

Commit 113982e

Browse files
authored
Merge pull request #38 from cval26/dct-types
Add type-based DCT procedure names
2 parents 53fc613 + 4bca750 commit 113982e

File tree

3 files changed

+104
-20
lines changed

3 files changed

+104
-20
lines changed

doc/specs/fftpack.md

Lines changed: 23 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -809,13 +809,16 @@ end program demo_dzfftb
809809

810810
## DCT type-1 (DCT-1)
811811

812-
### Initialize DCT-1: `dcosti`
812+
### Initialize DCT-1: `dcosti` or `dct_t1i`
813813

814814
#### Description
815815

816816
Initializes the array `wsave` which is used in subroutine `dcost`.
817817
The prime factorization of `n` together with a tabulation of the trigonometric functions are computed and stored in `wsave`.
818818

819+
The two procedures are completely equivalent and expect the same arguments.
820+
It is a matter of personal preference which one you choose to use.
821+
819822
#### Status
820823

821824
Experimental
@@ -851,18 +854,21 @@ program demo_dcosti
851854
end program demo_dcosti
852855
```
853856

854-
### Compute DCT-1: `dcost`
857+
### Compute DCT-1: `dcost` or `dct_t1`
855858

856859
#### Description
857860

858861
Computes the DCT-1 of the input real data.
859862
The transform is defined below at output parameter `x`.
860863

864+
The two procedures are completely equivalent and expect the same arguments.
865+
It is a matter of personal preference which one you choose to use.
866+
861867
For real input data `x` of length `n`, the DCT-1 of `x` is equivalent, up to a
862868
scaling factor, to the DFT of the even extension of `x` with length `2*(n-1)`,
863869
where the first and last entries of the original data are not repeated in the
864-
extension. For example, the DCT-1 of input data *abcde* (size \[5\]) is
865-
equivalent to the DFT of data *abcdedcb* (size \[2*4=8\]).
870+
extension. For example, the DCT-1 of input data *abcde* (size \(5\)) is
871+
equivalent to the DFT of data *abcdedcb* (size \(2*4=8\)).
866872

867873
Also, `dcost` is the unnormalized inverse of itself. This means that a call of
868874
`dcost` followed by another call of `dcost` will multiply the input sequence `x`
@@ -932,7 +938,7 @@ end program demo_dcost
932938

933939
## DCT of types 2, 3 (DCT-2, 3), a.k.a "Quarter" cosine transforms
934940

935-
### Initialize DCT-2, 3: `dcosqi`
941+
### Initialize DCT-2, 3: `dcosqi` or `dct_t23i`
936942

937943
#### Description
938944

@@ -941,6 +947,9 @@ The prime factorization of `n` together with
941947
a tabulation of the trigonometric functions are computed and
942948
stored in `wsave`.
943949

950+
The two procedures are completely equivalent and expect the same arguments.
951+
It is a matter of personal preference which one you choose to use.
952+
944953
#### Status
945954

946955
Experimental
@@ -978,13 +987,16 @@ program demo_dcosqi
978987
end program demo_dcosqi
979988
```
980989

981-
### Compute DCT-3: `dcosqf`
990+
### Compute DCT-3: `dcosqf` or `dct_t3`
982991

983992
#### Description
984993

985994
Computes the DCT-3 of the input real data.
986995
The transform is defined below at output parameter `x`.
987996

997+
The two procedures are completely equivalent and expect the same arguments.
998+
It is a matter of personal preference which one you choose to use.
999+
9881000
Also, `dcosqf` (DCT-3) is the unnormalized inverse of `dcosqb` (DCT-2), since a
9891001
call of `dcosqf` followed by a call of `dcosqb` will multiply the input sequence
9901002
`x` by `4*n`.
@@ -1049,13 +1061,16 @@ program demo_dcosqf
10491061
end program demo_dcosqf
10501062
```
10511063

1052-
### Compute DCT-2: `dcosqb`
1064+
### Compute DCT-2: `dcosqb` or `dct_t2`
10531065

10541066
#### Description
10551067

10561068
Computes the DCT-2 of the input real data.
10571069
The transform is defined below at output parameter `x`.
10581070

1071+
The two procedures are completely equivalent and expect the same arguments.
1072+
It is a matter of personal preference which one you choose to use.
1073+
10591074
For real input data `x` of length `n`, the DCT-2 of `x` is equivalent, up to a
10601075
scaling factor, to the DFT of the even extension of `x` with length `4*n`,
10611076
where all the even-frequency entries are zero.
@@ -1076,7 +1091,7 @@ Pure subroutine.
10761091

10771092
#### Syntax
10781093

1079-
`call [[fftpack(module):dcosqf(interface)]](n, x, wsave)`
1094+
`call [[fftpack(module):dcosqb(interface)]](n, x, wsave)`
10801095

10811096
#### Arguments
10821097

src/fftpack.f90

Lines changed: 50 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -17,6 +17,8 @@ module fftpack
1717
public :: dcosqi, dcosqf, dcosqb
1818
public :: dcosti, dcost
1919
public :: dct, idct
20+
public :: dct_t1i, dct_t1
21+
public :: dct_t23i, dct_t2, dct_t3
2022

2123
public :: rk
2224

@@ -125,7 +127,7 @@ end subroutine dzfftb
125127
!> Version: experimental
126128
!>
127129
!> Initialize `dcosqf` and `dcosqb`.
128-
!> ([Specification](../page/specs/fftpack.html#dcosqi))
130+
!> ([Specification](../page/specs/fftpack.html#initialize-dct-2-3-dcosqi-or-dct_t23i))
129131
pure subroutine dcosqi(n, wsave)
130132
import rk
131133
integer, intent(in) :: n
@@ -135,7 +137,7 @@ end subroutine dcosqi
135137
!> Version: experimental
136138
!>
137139
!> Forward transform of quarter wave data.
138-
!> ([Specification](../page/specs/fftpack.html#dcosqf))
140+
!> ([Specification](../page/specs/fftpack.html#compute-dct-3-dcosqf-or-dct_t3))
139141
pure subroutine dcosqf(n, x, wsave)
140142
import rk
141143
integer, intent(in) :: n
@@ -146,7 +148,7 @@ end subroutine dcosqf
146148
!> Version: experimental
147149
!>
148150
!> Unnormalized inverse of `dcosqf`.
149-
!> ([Specification](../page/specs/fftpack.html#dcosqb))
151+
!> ([Specification](../page/specs/fftpack.html#compute-dct-2-dcosqb-or-dct_t2))
150152
pure subroutine dcosqb(n, x, wsave)
151153
import rk
152154
integer, intent(in) :: n
@@ -156,7 +158,8 @@ end subroutine dcosqb
156158

157159
!> Version: experimental
158160
!>
159-
!> Initialize `dcost`. ([Specification](../page/specs/fftpack.html#dcosti))
161+
!> Initialize `dcost`.
162+
!> ([Specification](../page/specs/fftpack.html#initialize-dct-1-dcosti-or-dct_t1i))
160163
pure subroutine dcosti(n, wsave)
161164
import rk
162165
integer, intent(in) :: n
@@ -166,7 +169,7 @@ end subroutine dcosti
166169
!> Version: experimental
167170
!>
168171
!> Discrete fourier cosine transform of an even sequence.
169-
!> ([Specification](../page/specs/fftpack.html#dcost))
172+
!> ([Specification](../page/specs/fftpack.html#compute-dct-1-dcost-or-dct_t1))
170173
pure subroutine dcost(n, x, wsave)
171174
import rk
172175
integer, intent(in) :: n
@@ -245,7 +248,7 @@ end function irfft_rk
245248
!> Version: experimental
246249
!>
247250
!> Dsicrete cosine transforms.
248-
!> ([Specification](../page/specs/fftpack.html#dct))
251+
!> ([Specification](../page/specs/fftpack.html#simplified-dct-of-types-1-2-3-dct))
249252
interface dct
250253
pure module function dct_rk(x, n, type) result(result)
251254
real(kind=rk), intent(in) :: x(:)
@@ -258,7 +261,7 @@ end function dct_rk
258261
!> Version: experimental
259262
!>
260263
!> Inverse discrete cosine transforms.
261-
!> ([Specification](../page/specs/fftpack.html#idct))
264+
!> ([Specification](../page/specs/fftpack.html#simplified-inverse-dct-of-types-1-2-3-idct))
262265
interface idct
263266
pure module function idct_rk(x, n, type) result(result)
264267
real(kind=rk), intent(in) :: x(:)
@@ -268,6 +271,46 @@ pure module function idct_rk(x, n, type) result(result)
268271
end function idct_rk
269272
end interface idct
270273

274+
!> Version: experimental
275+
!>
276+
!> Initialize DCT type-1
277+
!> ([Specification](../page/specs/fftpack.html#initialize-dct-1-dcosti-or-dct_t1i))
278+
interface dct_t1i
279+
procedure :: dcosti
280+
end interface dct_t1i
281+
282+
!> Version: experimental
283+
!>
284+
!> Perform DCT type-1
285+
!> ([Specification](../page/specs/fftpack.html#compute-dct-1-dcost-or-dct_t1))
286+
interface dct_t1
287+
procedure :: dcost
288+
end interface dct_t1
289+
290+
!> Version: experimental
291+
!>
292+
!> Initialize DCT types 2, 3
293+
!> ([Specification](../page/specs/fftpack.html#initialize-dct-2-3-dcosqi-or-dct_t23i))
294+
interface dct_t23i
295+
procedure :: dcosqi
296+
end interface dct_t23i
297+
298+
!> Version: experimental
299+
!>
300+
!> Perform DCT type-2
301+
!> ([Specification](../page/specs/fftpack.html#compute-dct-2-dcosqb-or-dct_t2))
302+
interface dct_t2
303+
procedure :: dcosqb
304+
end interface dct_t2
305+
306+
!> Version: experimental
307+
!>
308+
!> Perform DCT type-3
309+
!> ([Specification](../page/specs/fftpack.html#compute-dct-3-dcosqf-or-dct_t3))
310+
interface dct_t3
311+
procedure :: dcosqf
312+
end interface dct_t3
313+
271314
!> Version: experimental
272315
!>
273316
!> Shifts zero-frequency component to center of spectrum.

test/test_fftpack_dct.f90

Lines changed: 31 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
module test_fftpack_dct
22

3-
use fftpack, only: rk, dcosti, dcost, dct, idct, dcosqi, dcosqf, dcosqb
3+
use fftpack
44
use testdrive, only: new_unittest, unittest_type, error_type, check
55
implicit none
66
private
@@ -24,27 +24,53 @@ end subroutine collect_dct
2424

2525
subroutine test_classic_dct(error)
2626
type(error_type), allocatable, intent(out) :: error
27-
real(kind=rk) :: w(3*4 + 15)
27+
real(kind=rk) :: w(3*4 + 15), w2(3*4 + 15)
2828
real(kind=rk) :: x(4) = [1, 2, 3, 4]
29+
real(kind=rk) :: x2(4)
2930
real(kind=rk) :: eps = 1.0e-10_rk
30-
31+
32+
x2 = x
3133
call dcosti(4, w)
3234
call dcost(4, x, w)
33-
call check(error, all(x == [real(kind=rk) :: 15, -4, 0, -1.0000000000000009_rk]), "`dcosti` failed.")
35+
call check(error, sum(abs(x - [real(kind=rk) :: 15, -4, 0, -1.0000000000000009_rk])) < eps, &
36+
"`dcost` failed.")
37+
if (allocated(error)) return
38+
39+
call dct_t1i(4, w2)
40+
call dct_t1(4, x2, w2)
41+
call check(error, maxval(abs(x2-x)) < eps, "dct_t1 failed")
3442
if (allocated(error)) return
43+
3544
call dcost(4, x, w)
36-
call check(error, all(x/(2*3) == [real(kind=rk) :: 1, 2, 3, 4]), "`dcost` failed.")
45+
call check(error, sum(abs(x/(2*3) - [real(kind=rk) :: 1, 2, 3, 4])) < eps, &
46+
"2nd `dcost` failed.")
47+
if (allocated(error)) return
48+
49+
call dct_t1(4, x2, w2)
50+
call check(error, maxval(abs(x2-x)) < eps, "2nd dct_t1 failed")
51+
if (allocated(error)) return
3752

3853
x = [1, 2, 3, 4]
54+
x2 = x
3955
call dcosqi(4, w)
4056
call dcosqf(4, x, w)
4157
call check(error, sum(abs(x - [11.999626276085150_rk, -9.1029432177492193_rk, &
4258
2.6176618435106480_rk, -1.5143449018465791_rk])) < eps, &
4359
"`dcosqf` failed.")
4460
if (allocated(error)) return
61+
62+
call dct_t23i(4, w2)
63+
call dct_t3(4, x2, w2)
64+
call check(error, maxval(abs(x2-x)) < eps, "dct_t3 failed")
65+
if (allocated(error)) return
66+
4567
call dcosqb(4, x, w)
4668
call check(error, sum(abs(x/(4*4) - [real(kind=rk) :: 1, 2, 3, 4])) < eps, &
4769
"`dcosqb` failed.")
70+
if (allocated(error)) return
71+
72+
call dct_t2(4, x2, w2)
73+
call check(error, maxval(abs(x2-x)) < eps, "dct_t2 failed")
4874

4975
end subroutine test_classic_dct
5076

0 commit comments

Comments
 (0)