diff options
author | Matti Picus <matti.picus@gmail.com> | 2021-02-15 16:56:32 +0200 |
---|---|---|
committer | GitHub <noreply@github.com> | 2021-02-15 16:56:32 +0200 |
commit | 9fd5ca55cc8c76759927831f44a98e2c65ea73d1 (patch) | |
tree | 8bd900ed589dd13b9806cbba0a9ccfdb4208f1ac | |
parent | ba4eab8681a12863269e5159910cbaeb06517ef9 (diff) | |
parent | a3af13823b309917b35953ba2e2ac861d6f7be58 (diff) | |
download | numpy-9fd5ca55cc8c76759927831f44a98e2c65ea73d1.tar.gz |
Merge pull request #18416 from pearu/gh-2848
BUG: Fix f2py parsing continued lines that follow comment lines.
-rwxr-xr-x | numpy/f2py/crackfortran.py | 20 | ||||
-rw-r--r-- | numpy/f2py/tests/test_crackfortran.py | 24 |
2 files changed, 41 insertions, 3 deletions
diff --git a/numpy/f2py/crackfortran.py b/numpy/f2py/crackfortran.py index 1149633c0..1bf724d12 100755 --- a/numpy/f2py/crackfortran.py +++ b/numpy/f2py/crackfortran.py @@ -341,7 +341,9 @@ def readfortrancode(ffile, dowithline=show, istop=1): if ffile == []: return localdolowercase = dolowercase - cont = 0 + # cont: set to True when the content of the last line read + # indicates statement continuation + cont = False finalline = '' ll = '' includeline = re.compile( @@ -392,14 +394,26 @@ def readfortrancode(ffile, dowithline=show, istop=1): if rl[:5].lower() == '!f2py': # f2py directive l, _ = split_by_unquoted(l + 4 * ' ' + rl[5:], '!') if l.strip() == '': # Skip empty line - cont = 0 + if sourcecodeform == 'free': + # In free form, a statement continues in the next line + # that is not a comment line [3.3.2.4^1], lines with + # blanks are comment lines [3.3.2.3^1]. Hence, the + # line continuation flag must retain its state. + pass + else: + # In fixed form, statement continuation is determined + # by a non-blank character at the 6-th position. Empty + # line indicates a start of a new statement + # [3.3.3.3^1]. Hence, the line continuation flag must + # be reset. + cont = False continue if sourcecodeform == 'fix': if l[0] in ['*', 'c', '!', 'C', '#']: if l[1:5].lower() == 'f2py': # f2py directive l = ' ' + l[5:] else: # Skip comment line - cont = 0 + cont = False continue elif strictf77: if len(l) > 72: diff --git a/numpy/f2py/tests/test_crackfortran.py b/numpy/f2py/tests/test_crackfortran.py index 827c71ae9..d26917f0c 100644 --- a/numpy/f2py/tests/test_crackfortran.py +++ b/numpy/f2py/tests/test_crackfortran.py @@ -115,3 +115,27 @@ class TestExternal(util.F2PyTest): return x + 123 r = self.module.external_as_attribute(incr) assert r == 123 + +class TestCrackFortran(util.F2PyTest): + + suffix = '.f90' + + code = textwrap.dedent(""" + subroutine gh2848( & + ! first 2 parameters + par1, par2,& + ! last 2 parameters + par3, par4) + + integer, intent(in) :: par1, par2 + integer, intent(out) :: par3, par4 + + par3 = par1 + par4 = par2 + + end subroutine gh2848 + """) + + def test_gh2848(self): + r = self.module.gh2848(1, 2) + assert r == (1, 2) |