diff options
author | Pearu Peterson <pearu.peterson@gmail.com> | 2021-01-18 20:23:18 +0200 |
---|---|---|
committer | Pearu Peterson <pearu.peterson@gmail.com> | 2021-01-18 20:23:18 +0200 |
commit | 1ef11c93a52e18ee2c45d4173ec8e685ceb53e9a (patch) | |
tree | 3eeefb5b3bd9f4907d247b23ce8134d43f5dd1ac | |
parent | 36439a7d0e0896703823c68d116b9900e49df998 (diff) | |
download | numpy-1ef11c93a52e18ee2c45d4173ec8e685ceb53e9a.tar.gz |
Add test for gh17797.
-rw-r--r-- | numpy/f2py/auxfuncs.py | 5 | ||||
-rwxr-xr-x | numpy/f2py/crackfortran.py | 11 | ||||
-rwxr-xr-x | numpy/f2py/rules.py | 13 | ||||
-rw-r--r-- | numpy/f2py/tests/test_callback.py | 24 |
4 files changed, 47 insertions, 6 deletions
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..2aa8dc420 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': + if block.get('result'): + # skip declaring function if its result is already declared + continue + else: + 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/rules.py b/numpy/f2py/rules.py index 12ed1f1ca..05fba9c4f 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) or issubroutine_wrap(nb): + if isf90: funcwrappers2.append(wrap) else: funcwrappers.append(wrap) @@ -1288,7 +1291,11 @@ 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] == ' ': + i = l.find('!') + if i >= 0 and i < 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:] diff --git a/numpy/f2py/tests/test_callback.py b/numpy/f2py/tests/test_callback.py index 81650a819..acf8c2392 100644 --- a/numpy/f2py/tests/test_callback.py +++ b/numpy/f2py/tests/test_callback.py @@ -211,3 +211,27 @@ class TestF77CallbackPythonTLS(TestF77Callback): compiler-provided """ options = ["-DF2PY_USE_PYTHON_TLS"] + + +class TestF90Callback(util.F2PyTest): + + suffix = '.f90' + + code = """ +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 |