Skip to content

Commit

Permalink
Add support for different real kinds
Browse files Browse the repository at this point in the history
  • Loading branch information
zoziha committed May 7, 2022
1 parent 8f1f4e4 commit b86491d
Show file tree
Hide file tree
Showing 6 changed files with 44 additions and 17 deletions.
5 changes: 4 additions & 1 deletion src/Makefile
Original file line number Diff line number Diff line change
Expand Up @@ -60,7 +60,7 @@ SRCF90 = \
fftpack_qct.f90\
fftpack_iqct.f90\
fftpack_dct.f90\
rk.f90
rk.F90

OBJF := $(SRCF:.f90=.o)
OBJF90 := $(SRCF90:.f90=.o)
Expand All @@ -74,6 +74,9 @@ shared: $(OBJ)
clean:
rm -f -r *.o *.a *.so *.mod *.smod

%.o: %.F90
$(FC) $(FFLAGS) -c $<

%.o: %.f90
$(FC) $(FFLAGS) -c $<

Expand Down
16 changes: 16 additions & 0 deletions src/rk.F90
Original file line number Diff line number Diff line change
@@ -0,0 +1,16 @@
!> fftpack kind
module fftpack_kind
implicit none

!> fftpack real kind
#if defined(fftpack_sp)
integer, parameter :: rk = selected_real_kind(6)
#elif defined(fftpack_xdp)
integer, parameter :: rk = selected_real_kind(18)
#elif defined(fftpack_qp)
integer, parameter :: rk = selected_real_kind(33)
#else
integer, parameter :: rk = selected_real_kind(15)
#endif

end module fftpack_kind
4 changes: 0 additions & 4 deletions src/rk.f90

This file was deleted.

2 changes: 1 addition & 1 deletion test/Makefile
Original file line number Diff line number Diff line change
Expand Up @@ -29,7 +29,7 @@ testdrive.F90:
$(FETCH) https://github.com/fortran-lang/test-drive/raw/v0.4.0/src/testdrive.F90 > $@

%.o: %.F90
$(FC) $(FFLAGS) -c $<
$(FC) $(FFLAGS) -I../src -c $<

%.o: %.f90
$(FC) $(FFLAGS) -I../src -c $<
Expand Down
25 changes: 17 additions & 8 deletions test/test_fftpack_dct.f90 → test/test_fftpack_dct.F90
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,12 @@ module test_fftpack_dct

public :: collect_dct

#if defined(fftpack_sp)
real(kind=rk) :: eps = 1.0e-5_rk
#else
real(kind=rk) :: eps = 1.0e-10_rk
#endif

contains

!> Collect all exported unit tests
Expand All @@ -26,15 +32,16 @@ subroutine test_classic_dct(error)
type(error_type), allocatable, intent(out) :: error
real(kind=rk) :: w(3*4 + 15)
real(kind=rk) :: x(4) = [1, 2, 3, 4]
real(kind=rk) :: eps = 1.0e-10_rk

call dcosti(4, w)
call dcost(4, x, w)
call check(error, all(x == [real(kind=rk) :: 15, -4, 0, -1.0000000000000009_rk]), "`dcosti` failed.")
call check(error, sum(abs(x - [real(kind=rk) :: 15, -4, 0, -1.0000000000000009_rk])) < eps, &
"`dcosti` failed.")
if (allocated(error)) return

call dcost(4, x, w)
call check(error, all(x/(2.0_rk*(4.0_rk - 1.0_rk)) == [real(kind=rk) :: 1, 2, 3, 4]), "`dcost` failed.")
call check(error, sum(abs(x/(2.0_rk*(4.0_rk - 1.0_rk)) - &
[real(kind=rk) :: 1, 2, 3, 4])) < eps, "`dcost` failed.")

end subroutine test_classic_dct

Expand All @@ -46,23 +53,25 @@ subroutine test_modernized_dct(error)
if (allocated(error)) return
call check(error, all(dct(x, 3) == dct(x)), "`dct(x, 3)` failed.")
if (allocated(error)) return
call check(error, all(dct(x, 4) == [real(kind=rk) :: -3, -3.0000000000000036_rk, 15, 33]), "`dct(x, 4)` failed.")
call check(error, sum(abs(dct(x, 4) - [real(kind=rk) :: -3, -3.0000000000000036_rk, 15, 33])) &
< eps, "`dct(x, 4)` failed.")

end subroutine test_modernized_dct

subroutine test_modernized_idct(error)
type(error_type), allocatable, intent(out) :: error
real(kind=rk) :: eps = 1.0e-10_rk
real(kind=rk) :: x(4) = [1, 2, 3, 4]

call check(error, all(idct(dct(x))/(2.0_rk*(4.0_rk - 1.0_rk)) == [real(kind=rk) :: 1, 2, 3, 4]), &
call check(error, sum(abs(idct(dct(x))/(2.0_rk*(4.0_rk - 1.0_rk)) - &
[real(kind=rk) :: 1, 2, 3, 4])) < eps, &
"`idct(dct(x))/(2.0_rk*(4.0_rk-1.0_rk))` failed.")
if (allocated(error)) return
call check(error, all(idct(dct(x), 2)/(2.0_rk*(2.0_rk - 1.0_rk)) == [real(kind=rk) :: 5.5, 9.5]), &
"`idct(dct(x), 2)/(2.0_rk*(2.0_rk-1.0_rk))` failed.")
if (allocated(error)) return
call check(error, all(idct(dct(x, 2), 4)/(2.0_rk*(4.0_rk - 1.0_rk)) == &
[0.16666666666666666_rk, 0.33333333333333331_rk, 0.66666666666666663_rk, 0.83333333333333315_rk]), &
call check(error, sum(abs(idct(dct(x, 2), 4)/(2.0_rk*(4.0_rk - 1.0_rk)) - &
[0.16666666666666666_rk, 0.33333333333333331_rk, &
0.66666666666666663_rk, 0.83333333333333315_rk])) < eps, &
"`idct(dct(x, 2), 4)/(2.0_rk*(4.0_rk-1.0_rk))` failed.")

end subroutine test_modernized_idct
Expand Down
9 changes: 6 additions & 3 deletions test/test_fftpack_qct.f90 → test/test_fftpack_qct.F90
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,12 @@ module test_fftpack_qct

public :: collect_qct

#if defined(fftpack_sp)
real(kind=rk) :: eps = 1.0e-5_rk
#else
real(kind=rk) :: eps = 1.0e-10_rk
#endif

contains

!> Collect all exported unit tests
Expand All @@ -26,7 +32,6 @@ subroutine test_classic_qct(error)
type(error_type), allocatable, intent(out) :: error
real(kind=rk) :: w(3*4 + 15)
real(kind=rk) :: x(4) = [1, 2, 3, 4]
real(kind=rk) :: eps = 1.0e-10_rk

call dcosqi(4, w)
call dcosqf(4, x, w)
Expand All @@ -42,7 +47,6 @@ end subroutine test_classic_qct

subroutine test_modernized_qct(error)
type(error_type), allocatable, intent(out) :: error
real(kind=rk) :: eps = 1.0e-10_rk
real(kind=rk) :: x(3) = [9, -9, 3]

call check(error, sum(abs(qct(x, 2) - [-3.7279220613578570_rk, 21.727922061357859_rk])) < eps, &
Expand All @@ -59,7 +63,6 @@ end subroutine test_modernized_qct

subroutine test_modernized_iqct(error)
type(error_type), allocatable, intent(out) :: error
real(kind=rk) :: eps = 1.0e-10_rk
real(kind=rk) :: x(4) = [1, 2, 3, 4]

call check(error, sum(abs(iqct(qct(x))/(4.0_rk*4.0_rk) - [real(kind=rk) :: 1, 2, 3, 4])) < eps, &
Expand Down

0 comments on commit b86491d

Please sign in to comment.