summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorCharles Harris <charlesr.harris@gmail.com>2021-01-19 14:32:51 -0700
committerGitHub <noreply@github.com>2021-01-19 14:32:51 -0700
commitefaf210f7128f32e6d3d193f79ccc6b43ff8c006 (patch)
treee16f04c0d09a784f43ae944302cbbeb55dd35c93
parent09416c2934a5f592970e1d97c0298924eed009f6 (diff)
parentc4014c7fb4f1753f296bac842098decf061e39bb (diff)
downloadnumpy-efaf210f7128f32e6d3d193f79ccc6b43ff8c006.tar.gz
Merge pull request #18184 from pearu/17797
BUG: Fix f2py bugs when wrapping F90 subroutines.
-rw-r--r--numpy/distutils/fcompiler/__init__.py2
-rw-r--r--numpy/f2py/auxfuncs.py5
-rwxr-xr-xnumpy/f2py/crackfortran.py11
-rw-r--r--numpy/f2py/func2subr.py9
-rwxr-xr-xnumpy/f2py/rules.py17
-rw-r--r--numpy/f2py/tests/test_callback.py25
6 files changed, 58 insertions, 11 deletions
diff --git a/numpy/distutils/fcompiler/__init__.py b/numpy/distutils/fcompiler/__init__.py
index 4730a5a09..812461538 100644
--- a/numpy/distutils/fcompiler/__init__.py
+++ b/numpy/distutils/fcompiler/__init__.py
@@ -976,7 +976,7 @@ def is_free_format(file):
with open(file, encoding='latin1') as f:
line = f.readline()
n = 10000 # the number of non-comment lines to scan for hints
- if _has_f_header(line):
+ if _has_f_header(line) or _has_fix_header(line):
n = 0
elif _has_f90_header(line):
n = 0
diff --git a/numpy/f2py/auxfuncs.py b/numpy/f2py/auxfuncs.py
index 80b150655..5250fea84 100644
--- a/numpy/f2py/auxfuncs.py
+++ b/numpy/f2py/auxfuncs.py
@@ -257,6 +257,7 @@ def ismodule(rout):
def isfunction(rout):
return 'block' in rout and 'function' == rout['block']
+
def isfunction_wrap(rout):
if isintent_c(rout):
return 0
@@ -284,6 +285,10 @@ def hasassumedshape(rout):
return False
+def requiresf90wrapper(rout):
+ return ismoduleroutine(rout) or hasassumedshape(rout)
+
+
def isroutine(rout):
return isfunction(rout) or issubroutine(rout)
diff --git a/numpy/f2py/crackfortran.py b/numpy/f2py/crackfortran.py
index d27845796..1149633c0 100755
--- a/numpy/f2py/crackfortran.py
+++ b/numpy/f2py/crackfortran.py
@@ -3113,7 +3113,7 @@ def crack2fortrangen(block, tab='\n', as_interface=False):
result = ' result (%s)' % block['result']
if block['result'] not in argsl:
argsl.append(block['result'])
- body = crack2fortrangen(block['body'], tab + tabchar)
+ body = crack2fortrangen(block['body'], tab + tabchar, as_interface=as_interface)
vars = vars2fortran(
block, block['vars'], argsl, tab + tabchar, as_interface=as_interface)
mess = ''
@@ -3231,8 +3231,13 @@ def vars2fortran(block, vars, args, tab='', as_interface=False):
show(vars)
outmess('vars2fortran: No definition for argument "%s".\n' % a)
continue
- if a == block['name'] and not block['block'] == 'function':
- continue
+ if a == block['name']:
+ if block['block'] != 'function' or block.get('result'):
+ # 1) skip declaring a variable that name matches with
+ # subroutine name
+ # 2) skip declaring function when its type is
+ # declared via `result` construction
+ continue
if 'typespec' not in vars[a]:
if 'attrspec' in vars[a] and 'external' in vars[a]['attrspec']:
if a in args:
diff --git a/numpy/f2py/func2subr.py b/numpy/f2py/func2subr.py
index e9976f43c..21d4c009c 100644
--- a/numpy/f2py/func2subr.py
+++ b/numpy/f2py/func2subr.py
@@ -130,7 +130,7 @@ def createfuncwrapper(rout, signature=0):
l = l + ', ' + fortranname
if need_interface:
for line in rout['saved_interface'].split('\n'):
- if line.lstrip().startswith('use '):
+ if line.lstrip().startswith('use ') and '__user__' not in line:
add(line)
args = args[1:]
@@ -222,7 +222,7 @@ def createsubrwrapper(rout, signature=0):
if need_interface:
for line in rout['saved_interface'].split('\n'):
- if line.lstrip().startswith('use '):
+ if line.lstrip().startswith('use ') and '__user__' not in line:
add(line)
dumped_args = []
@@ -247,7 +247,10 @@ def createsubrwrapper(rout, signature=0):
pass
else:
add('interface')
- add(rout['saved_interface'].lstrip())
+ for line in rout['saved_interface'].split('\n'):
+ if line.lstrip().startswith('use ') and '__user__' in line:
+ continue
+ add(line)
add('end interface')
sargs = ', '.join([a for a in args if a not in extra_args])
diff --git a/numpy/f2py/rules.py b/numpy/f2py/rules.py
index f1490527e..4e1cf0c7d 100755
--- a/numpy/f2py/rules.py
+++ b/numpy/f2py/rules.py
@@ -73,7 +73,7 @@ from .auxfuncs import (
issubroutine, issubroutine_wrap, isthreadsafe, isunsigned,
isunsigned_char, isunsigned_chararray, isunsigned_long_long,
isunsigned_long_longarray, isunsigned_short, isunsigned_shortarray,
- l_and, l_not, l_or, outmess, replace, stripcomma,
+ l_and, l_not, l_or, outmess, replace, stripcomma, requiresf90wrapper
)
from . import capi_maps
@@ -1184,9 +1184,12 @@ def buildmodule(m, um):
nb1['args'] = a
nb_list.append(nb1)
for nb in nb_list:
+ # requiresf90wrapper must be called before buildapi as it
+ # rewrites assumed shape arrays as automatic arrays.
+ isf90 = requiresf90wrapper(nb)
api, wrap = buildapi(nb)
if wrap:
- if ismoduleroutine(nb):
+ if isf90:
funcwrappers2.append(wrap)
else:
funcwrappers.append(wrap)
@@ -1288,7 +1291,10 @@ def buildmodule(m, um):
'C It contains Fortran 77 wrappers to fortran functions.\n')
lines = []
for l in ('\n\n'.join(funcwrappers) + '\n').split('\n'):
- if l and l[0] == ' ':
+ if 0 <= l.find('!') < 66:
+ # don't split comment lines
+ lines.append(l + '\n')
+ elif l and l[0] == ' ':
while len(l) >= 66:
lines.append(l[:66] + '\n &')
l = l[66:]
@@ -1310,7 +1316,10 @@ def buildmodule(m, um):
'! It contains Fortran 90 wrappers to fortran functions.\n')
lines = []
for l in ('\n\n'.join(funcwrappers2) + '\n').split('\n'):
- if len(l) > 72 and l[0] == ' ':
+ if 0 <= l.find('!') < 72:
+ # don't split comment lines
+ lines.append(l + '\n')
+ elif len(l) > 72 and l[0] == ' ':
lines.append(l[:72] + '&\n &')
l = l[72:]
while len(l) > 66:
diff --git a/numpy/f2py/tests/test_callback.py b/numpy/f2py/tests/test_callback.py
index 81650a819..37736af21 100644
--- a/numpy/f2py/tests/test_callback.py
+++ b/numpy/f2py/tests/test_callback.py
@@ -211,3 +211,28 @@ class TestF77CallbackPythonTLS(TestF77Callback):
compiler-provided
"""
options = ["-DF2PY_USE_PYTHON_TLS"]
+
+
+class TestF90Callback(util.F2PyTest):
+
+ suffix = '.f90'
+
+ code = textwrap.dedent(
+ """
+ function gh17797(f, y) result(r)
+ external f
+ integer(8) :: r, f
+ integer(8), dimension(:) :: y
+ r = f(0)
+ r = r + sum(y)
+ end function gh17797
+ """)
+
+ def test_gh17797(self):
+
+ def incr(x):
+ return x + 123
+
+ y = np.array([1, 2, 3], dtype=np.int64)
+ r = self.module.gh17797(incr, y)
+ assert r == 123 + 1 + 2 + 3