summaryrefslogtreecommitdiff
path: root/numpy/f2py/tests/src
diff options
context:
space:
mode:
authorNick Papior <nickpapior@gmail.com>2016-11-23 09:11:52 +0100
committerNick Papior <nickpapior@gmail.com>2016-11-27 21:10:04 +0100
commit8a7e7a44c194e3c748ccee0862f4ae01e371a208 (patch)
tree45eaa09569e3fd84f1d9899cd8de476344b561ba /numpy/f2py/tests/src
parentce6d3ffc991c63c1a91fd806b4eb3d0008caace2 (diff)
downloadnumpy-8a7e7a44c194e3c748ccee0862f4ae01e371a208.tar.gz
BUG: fixed kind specifications for parameters
Fortran sources with parameters having kind-specifiers where not correctly intercepted in the crackfortran.py source. The reason was a restrictive check for only integer specifiers which did not split real's into the correct number. Furthermore, several tests has been added which tests the different kind specifiers and their use in codes, also all of them together. Signed-off-by: Nick Papior <nickpapior@gmail.com>
Diffstat (limited to 'numpy/f2py/tests/src')
-rw-r--r--numpy/f2py/tests/src/parameter/constant_both.f9057
-rw-r--r--numpy/f2py/tests/src/parameter/constant_integer.f9022
-rw-r--r--numpy/f2py/tests/src/parameter/constant_real.f9023
3 files changed, 102 insertions, 0 deletions
diff --git a/numpy/f2py/tests/src/parameter/constant_both.f90 b/numpy/f2py/tests/src/parameter/constant_both.f90
new file mode 100644
index 000000000..ac90cedc5
--- /dev/null
+++ b/numpy/f2py/tests/src/parameter/constant_both.f90
@@ -0,0 +1,57 @@
+! Check that parameters are correct intercepted.
+! Constants with comma separations are commonly
+! used, for instance Pi = 3._dp
+subroutine foo(x)
+ implicit none
+ integer, parameter :: sp = selected_real_kind(6)
+ integer, parameter :: dp = selected_real_kind(15)
+ integer, parameter :: ii = selected_int_kind(9)
+ integer, parameter :: il = selected_int_kind(18)
+ real(dp), intent(inout) :: x
+ dimension x(3)
+ real(sp), parameter :: three_s = 3._sp
+ real(dp), parameter :: three_d = 3._dp
+ integer(ii), parameter :: three_i = 3_ii
+ integer(il), parameter :: three_l = 3_il
+ x(1) = x(1) + x(2) * three_s * three_i + x(3) * three_d * three_l
+ x(2) = x(2) * three_s
+ x(3) = x(3) * three_l
+ return
+end subroutine
+
+
+subroutine foo_no(x)
+ implicit none
+ integer, parameter :: sp = selected_real_kind(6)
+ integer, parameter :: dp = selected_real_kind(15)
+ integer, parameter :: ii = selected_int_kind(9)
+ integer, parameter :: il = selected_int_kind(18)
+ real(dp), intent(inout) :: x
+ dimension x(3)
+ real(sp), parameter :: three_s = 3.
+ real(dp), parameter :: three_d = 3.
+ integer(ii), parameter :: three_i = 3
+ integer(il), parameter :: three_l = 3
+ x(1) = x(1) + x(2) * three_s * three_i + x(3) * three_d * three_l
+ x(2) = x(2) * three_s
+ x(3) = x(3) * three_l
+ return
+end subroutine
+
+subroutine foo_sum(x)
+ implicit none
+ integer, parameter :: sp = selected_real_kind(6)
+ integer, parameter :: dp = selected_real_kind(15)
+ integer, parameter :: ii = selected_int_kind(9)
+ integer, parameter :: il = selected_int_kind(18)
+ real(dp), intent(inout) :: x
+ dimension x(3)
+ real(sp), parameter :: three_s = 2._sp + 1._sp
+ real(dp), parameter :: three_d = 1._dp + 2._dp
+ integer(ii), parameter :: three_i = 2_ii + 1_ii
+ integer(il), parameter :: three_l = 1_il + 2_il
+ x(1) = x(1) + x(2) * three_s * three_i + x(3) * three_d * three_l
+ x(2) = x(2) * three_s
+ x(3) = x(3) * three_l
+ return
+end subroutine
diff --git a/numpy/f2py/tests/src/parameter/constant_integer.f90 b/numpy/f2py/tests/src/parameter/constant_integer.f90
new file mode 100644
index 000000000..aaa83d2eb
--- /dev/null
+++ b/numpy/f2py/tests/src/parameter/constant_integer.f90
@@ -0,0 +1,22 @@
+! Check that parameters are correct intercepted.
+! Constants with comma separations are commonly
+! used, for instance Pi = 3._dp
+subroutine foo_int(x)
+ implicit none
+ integer, parameter :: ii = selected_int_kind(9)
+ integer(ii), intent(inout) :: x
+ dimension x(3)
+ integer(ii), parameter :: three = 3_ii
+ x(1) = x(1) + x(2) + x(3) * three
+ return
+end subroutine
+
+subroutine foo_long(x)
+ implicit none
+ integer, parameter :: ii = selected_int_kind(18)
+ integer(ii), intent(inout) :: x
+ dimension x(3)
+ integer(ii), parameter :: three = 3_ii
+ x(1) = x(1) + x(2) + x(3) * three
+ return
+end subroutine
diff --git a/numpy/f2py/tests/src/parameter/constant_real.f90 b/numpy/f2py/tests/src/parameter/constant_real.f90
new file mode 100644
index 000000000..02ac9dd99
--- /dev/null
+++ b/numpy/f2py/tests/src/parameter/constant_real.f90
@@ -0,0 +1,23 @@
+! Check that parameters are correct intercepted.
+! Constants with comma separations are commonly
+! used, for instance Pi = 3._dp
+subroutine foo_single(x)
+ implicit none
+ integer, parameter :: rp = selected_real_kind(6)
+ real(rp), intent(inout) :: x
+ dimension x(3)
+ real(rp), parameter :: three = 3._rp
+ x(1) = x(1) + x(2) + x(3) * three
+ return
+end subroutine
+
+subroutine foo_double(x)
+ implicit none
+ integer, parameter :: rp = selected_real_kind(15)
+ real(rp), intent(inout) :: x
+ dimension x(3)
+ real(rp), parameter :: three = 3._rp
+ x(1) = x(1) + x(2) + x(3) * three
+ return
+end subroutine
+