Skip to content

Navigation Menu

Sign in
Appearance settings

Search code, repositories, users, issues, pull requests...

Provide feedback

We read every piece of feedback, and take your input very seriously.

Saved searches

Use saved searches to filter your results more quickly

Appearance settings

BUG: Modules from Fortran interface not accessible in Numpy 2.x (current 2.1.2) #27622

Copy link
Copy link
Closed
@2sn

Description

@2sn
Issue body actions

Describe the issue:

Using Numpy 2.x (tried 2.02, 2.1.2) I can no longer access module data.

After f2py I get the pyf file

!    -*- f90 -*-
! Note: the context of this file is case sensitive.

python module _solver ! in 
    interface  ! in :_solver
        module types ! in :_solver:solver.f90
            integer, parameter,optional :: int32=selected_int_kind(8)
            integer, parameter,optional :: int64=selected_int_kind(16)
            integer, parameter,optional :: real64=selected_real_kind(15)
        end module types
        module lanedata ! in :_solver:solver.f90
            use types, only: real64,int32
            integer(kind=4), parameter,optional :: maxdata=1048575
            real(kind=8), allocatable,dimension(:,:) :: theta
            integer(kind=4) :: ndata
        end module lanedata
        subroutine freelanedata ! in :_solver:solver.f90
            use lanedata, only: ndata,theta
        end subroutine freelanedata
        subroutine lane(dx,n,w) ! in :_solver:solver.f90
            use types, only: int32,real64
            use lanedata, only: maxdata,theta,ndata
            real(kind=8) intent(in) :: dx
            real(kind=8) intent(in) :: n
            real(kind=8) intent(in) :: w
        end subroutine lane
        subroutine rk4(x0,y0,y1,dx,n,w,z0,z1) ! in :_solver:solver.f90
            use types, only: real64
            real(kind=8) intent(in) :: x0
            real(kind=8) intent(in) :: y0
            real(kind=8) intent(in) :: y1
            real(kind=8) intent(in) :: dx
            real(kind=8) intent(in) :: n
            real(kind=8) intent(in) :: w
            real(kind=8) intent(out) :: z0
            real(kind=8) intent(out) :: z1
        end subroutine rk4
    end interface 
end python module _solver

! This file was auto-generated with f2py (version:2.1.2).
! See:
! https://web.archive.org/web/20140822061353/http://cens.ioc.ee/projects/f2py2e

but when I try to access the data, as suggested in the online doc (just checking that this is still as it used to be, https://numpy.org/doc/2.1/f2py/python-usage.html#fortran-90-module-data), I get

 In [2]: print(laneemden._solver.__doc__)
This module '_solver' is auto-generated with f2py (version:2.1.2).
Functions:
    freelanedata()
    lane(dx,n,w)
    z0,z1 = rk4(x0,y0,y1,dx,n,w)
.

or

In [3]: dir(laneemden._solver)
Out[3]: 
['__doc__',
 '__f2py_numpy_version__',
 '__file__',
 '__loader__',
 '__name__',
 '__package__',
 '__solver_error',
 '__spec__',
 '__version__',
 'freelanedata',
 'lane',
 'rk4']

so, not trace of module lanedata. This used to work (and still does) in Numpy 1.x.

Reproduce the code example:

from laneemden._solver import lanedata

Here my code, solver.f90

! 20111212 Alexander Heger

! 1/z**2 d/dz (z**2 d/dz theta(z)) + theta(z)**n = 0
! for small z one can approximate
! theta(z) = 1. + (-1/6.)*z**2 + (n/120.)*z**4 + O(z**6)
! Therefore lim(z)-->0 d**2 theta(z)/d z**2 = -1/3

! 20120423 Alexander Heger

! if we include a constant rotation rate Omega the equation becomes
! 1/z**2 d/dz (z**2 d/dz theta(z)) + theta(z)**n - w = 0
! where
! w = W/rho_c
! W = 2 Omega**2 / 4 pi G
! for small z one can approximate
! theta(z) = 1. + (w - 1)/6. *z**2 + (1.-w)*(n/120.)*z**4 + O(z**6)
! Therefore lim(z)-->0 d**2 theta(z)/d z**2 = (w-1.)/3

module types

  implicit none

  INTEGER, PARAMETER :: int32 = SELECTED_INT_KIND(8)
  INTEGER, PARAMETER :: int64 = SELECTED_INT_KIND(16)
  INTEGER, PARAMETER :: real64 = SELECTED_REAL_KIND(15)

end module types


module lanedata

  use types, only: &
       real64, int32

  implicit none

  save

  integer(kind=int32), parameter :: &
       maxdata = 2**20-1

  real(kind=real64), dimension(:,:), allocatable :: &
       theta

  integer(kind=int32) :: &
       ndata

end module lanedata


subroutine freelanedata()

  use lanedata, only: &
       ndata, theta

  implicit none

  if (allocated(theta)) deallocate(theta)
  ndata = -1

end subroutine freelanedata


subroutine lane(dx, n, w)

  ! lane emden integration, return data,
  ! including first invalid (rho < 0) point for interpolation.

  use types, only: &
       int32, real64

  use lanedata, only: &
       maxdata, theta, ndata

  implicit none

!f2py real(kind=real64), intent(in) :: dx, n, w

  real(kind=real64), intent(in) :: &
       dx, n, w

  real(kind=real64) :: &
       x
  integer(kind=int32) :: &
       i, j

  real(kind=real64), dimension(0:1) :: &
       y, z

  real(kind=real64), dimension(:,:), allocatable :: &
       theta_

  if (allocated(theta)) deallocate(theta)
  j = maxdata
  allocate(theta_(0:j, 0:1))

  x = 0.d0
  y(0:1) = (/1.d0, 0.d0/)

  i = 0
  theta_(i,:) = y(:)
  do while (.True.)
     call rk4(x,y(0),y(1),dx,n,w,z(0),z(1))
     if (i == j) then
        j = j + maxdata
        allocate(theta(0:j,0:1))
        theta(0:i,:) = theta_(0:i,:)
        deallocate(theta_)
        call move_alloc(theta, theta_)
     endif

     i = i+1
     theta_(i, :) = z(:)
     x = x + dx
     y(:) = z(:)

     if (y(0) < 0.d0) exit
  enddo
  ndata = i

  allocate(theta(0:ndata, 0:1))
  theta(0:ndata,0:1) = theta_(0:ndata,0:1)
  deallocate(theta_)

end subroutine lane


subroutine rk4(x0, y0, y1, dx, n, w, z0, z1)

  use types, only: &
       real64

  implicit none

  real(kind=real64), parameter :: &
       p13 = 1.d0 / 3.d0, &
       p16 = 1.d0 / 6.d0

  real(kind=real64), intent(in)  :: &
       x0, y0, y1
  real(kind=real64), intent(in)  :: &
       dx, n, w
  real(kind=real64), intent(out) :: &
       z0, z1

  real(kind=real64) :: &
       xh, dh
  real(kind=real64) :: &
       k10, k11, k20, k21, k30, k31, k40, k41

!f2py real(8), intent(in) :: x0, dx, n, w
!f2py real(8), intent(in) :: y0, y1
!f2py real(8), intent(out) :: z0, z1

  xh = x0 + 0.5d0 * dx
  dh = 0.5d0 * dx

  k10 = y1
  if (x0 == 0) then
     k11 = (w - 1.d0) * p13
  else
     k11 = -2.d0 / x0 * y1 - (max(y0, 0.d0))**n + w
  endif

  k20 = y1 + dh*k11
  k21 = -2.d0 / xh * k20 - (max(y0 + dh*k10, 0.d0))**n

  k30 = y1 + dh*k21
  k31 = -2.d0 / xh * k30 - (max(y0 + dh*k20, 0.d0))**n

  k40 = y1 + dx*k31
  k41 = -2.d0 / (x0+dx) * k40 - (max(y0 + dx*k30,0.d0))**n

  z0 = y0 + dx*(k10 + 2.d0 * (k20 + k30) + k40) * p16
  z1 = y1 + dx*(k11 + 2.d0 * (k21 + k31) + k41) * p16

end subroutine rk4

I use a custom build package that does not work with meson (and sadly, the usually very helpful Numpy developers out of principle refused to include fixes to obvious bugs in the build script to this day and a small patch to make my script continue to work, for both of which I had posted patches) and that are lengthy, so I won't include.

The issue seems to be that there is on interface in module lanedata. If I include an interface into the definition of module lanedata so it becomes, e.g.,

module lanedata

  use types, only: &
       real64, int32

  implicit none

  save

  integer(kind=int32), parameter :: &
       maxdata = 2**20-1

  real(kind=real64), dimension(:,:), allocatable :: &
       theta

  integer(kind=int32) :: &
       ndata

contains

  subroutine phoney
  end subroutine phoney
  
end module lanedata

then the module is included and an interface generated. This behaviour to skip modules w/o interface contents (just adding the contains line but an empty section does not suffice).

Maybe the Python 1.x behaviour can be restored, as it is common to have plain data modules w/o definition of functions or subroutines.

It is also possible that I missed a flag / behaviour change, my apologies in this case and please advise.

Error message:

---------------------------------------------------------------------------
ImportError                               Traceback (most recent call last)
Cell In[4], line 1
----> 1 from laneemden._solver import lanedata

ImportError: cannot import name 'lanedata' from 'laneemden._solver' (/home/alex/python/source/laneemden/_solver.cpython-311-x86_64-linux-gnu.so)

Python and NumPy Versions:

2.1.2
3.11.10 (main, Sep 8 2024, 14:25:06) [GCC 14.2.1 20240801 (Red Hat 14.2.1-1)]

Runtime Environment:

[{'numpy_version': '2.1.2',
'python': '3.11.10 (main, Sep 8 2024, 14:25:06) [GCC 14.2.1 20240801 (Red '
'Hat 14.2.1-1)]',
'uname': uname_result(system='Linux', node='w.2sn.net', release='6.11.3-200.fc40.x86_64', version='#1 SMP PREEMPT_DYNAMIC Thu Oct 10 22:31:19 UTC 2024', machine='x86_64')},
{'simd_extensions': {'baseline': ['SSE', 'SSE2', 'SSE3'],
'found': ['SSSE3',
'SSE41',
'POPCNT',
'SSE42',
'AVX',
'F16C',
'FMA3',
'AVX2',
'AVX512F',
'AVX512CD',
'AVX512_SKX'],
'not_found': ['AVX512_KNL',
'AVX512_KNM',
'AVX512_CLX',
'AVX512_CNL',
'AVX512_ICL']}},
{'architecture': 'SkylakeX',
'filepath': '/home/alex/Python_3.11.10/lib/python3.11/site-packages/numpy.libs/libscipy_openblas64_-ff651d7f.so',
'internal_api': 'openblas',
'num_threads': 16,
'prefix': 'libscipy_openblas',
'threading_layer': 'pthreads',
'user_api': 'blas',
'version': '0.3.27'}]

Context for the issue:

breaks code, can no longer access module data that worked in Numpy < 2.0.

Metadata

Metadata

Assignees

No one assigned

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions

      Morty Proxy This is a proxified and sanitized view of the page, visit original site.