Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

bug: F08 makes array 0-based instead of the original 1-based #4170

Closed
hzhou opened this issue Nov 13, 2019 · 11 comments
Closed

bug: F08 makes array 0-based instead of the original 1-based #4170

hzhou opened this issue Nov 13, 2019 · 11 comments

Comments

@hzhou
Copy link
Contributor

hzhou commented Nov 13, 2019

It is reported by a user on discuss@mpich.org:

Hi there, My name is Pedro Ricardo and I'm currently a researcher at a Fluid Dynamics Laboratory in Brazil.
We have a CFD software with 12 continuous years of development mostly in Fortran, and we've been using Mpich and Openmpi since the very beginning.

Recently we started to refactor the code and write it using Fortran2008 standard and this includes the usage of the the newest MPI module for Fortran "mpi_f08".

So far no errors were present in Openmpi implementations and everything worked like a charm, but when we started compiling with Mpich lots of errors happened during execution.

I was able so isolate and reproduce the error in a simpler code. Whenever "use mpi_f08" is present, the arrays passed to a MPI subroutine return with C-like indexes. If I allocate an standard array in Fortran (that starts with index 1), after the call of some MPI routine this array come back starting with index 0.
And you can imagine that all access to this array will be wrong after that.
It is a funny behavior because if you put "use mpi" instead of "use mpi_f08", the indexes do not change.

The test is simple, allocate and array, print it's bound indexes, calls a MPI_Reduce routine and print the array bounds again. I've also tested other routines like MPI_ISend and the error also happens.

program check_bounds
    !Check Works
    ! use mpi
    !Check Fails
    use mpi_f08
   
    implicit none
    integer:: mpi_ierr, proc_id, n_proc
    real, allocatable, dimension(:):: t1, t1_0 

    call MPI_Init(mpi_ierr)  
    call MPI_Comm_rank(MPI_COMM_WORLD, proc_id, mpi_ierr)
    call MPI_Comm_size(MPI_COMM_WORLD, n_proc, mpi_ierr) 

    allocate(t1(5), t1_0(5))
    t1 = 1.0; t1_0 = 1.0 

    if (proc_id==0) then
        write(*,*)'Array Configuration Before:'
        write(*,'(a6,2(1x,i2))')'Low:',lbound(t1),lbound(t1_0)
        write(*,'(a6,2(1x,i2))')'Upper:',ubound(t1),ubound(t1_0)
    end if 

    CALL MPI_REDUCE(t1, t1_0, 5, MPI_REAL, MPI_SUM, 0, MPI_COMM_WORLD, mpi_ierr) 

    if (proc_id==0) then
        write(*,*)'Array Configuration After:'
        write(*,'(a6,2(1x,i2))')'Low:',lbound(t1),lbound(t1_0)
        write(*,'(a6,2(1x,i2))')'Upper:',ubound(t1),ubound(t1_0)
    end if 

    deallocate(t1, t1_0)
    call MPI_Finalize(mpi_ierr) 

end program check_bounds
@hzhou
Copy link
Contributor Author

hzhou commented Nov 21, 2019

Follow up:
Slightly modified the test program to also print "1st" element:

        ! explicitly initialize arrays
        t1 = (/ 1, 2, 3, 4, 5 /);
        t1_0 = (/ 0, 0, 0, 0, 0 /);
        ...
        write(*,'(a6,2(1x,i2),f4.0)')'Low:',lbound(t1),lbound(t1_0), t1(1)
        write(*,'(a6,2(1x,i2),f4.0)')'Upper:',ubound(t1),ubound(t1_0), t1_0(1)

Test with:

mpifort test.f90 && ./a.out
  • Tried with gcc-9.2, it confirms the issue:
Array Configuration Before:
  Low:  1  1 1.
Upper:  5  5 0.
 Array Configuration After:
  Low:  0  0 2.
Upper:  4  4 2.
  • Tried with icc (ver. 19.0.3):
 Array Configuration Before:
  Low:  1  1 1.
Upper:  5  5 0.
 Array Configuration After:
  Low:  1  1 1.
Upper:  5  5 1.

@hzhou
Copy link
Contributor Author

hzhou commented Nov 21, 2019

I probably will look into it further and identify the exact issue, but it appears to be an issue with current gfortran implementation.

Cross referencing in issue #3820, there reveals also a discrepancy related to the Fortran iso-c-binding between gfortran and ifort. In that case, our code works for gfortran but fails ifort (for legacy reason), even though I believe ifort is more correct in its implementation.

@pedro-ricardo
Copy link

Hello there, any updates on this?
Can I help with something?

@hzhou
Copy link
Contributor Author

hzhou commented Mar 3, 2020

@pedro-ricardo We believe at this point this is an issue of gfortran. Intel compilers (ifort) seems ok.

@pedro-ricardo
Copy link

pedro-ricardo commented Mar 3, 2020

Can it be isolated from MPICH and submitted as a issue to gcc?

Do you intend on making a workaround to enable mpi_f08 with gfortran while the bug isn't solved there?

@hzhou
Copy link
Contributor Author

hzhou commented Mar 3, 2020

Can it be isolated from MPICH and submitted as a issue to gcc?

The answer is probably yes. It will take some time (and unfortunately, some priority). Any help will be appreciated. If you would like to help, feel free to feed back to us on your effort. We'll be happy to discuss and provide pointers.

Do you intend on making a workaround to enable mpi_f08 with gfortran while the bug isn't solved there?

The work-around is to skip iso-c-binding and do it in a similar way as the f90 modules. Of course, doing that defeats the original purpose of the f08 interface -- for its stronger type checking.

@gcorbin
Copy link

gcorbin commented Jan 14, 2022

@pedro-ricardo: Have you submitted the bug to GCC yet?

@pedro-ricardo
Copy link

Hello @gcorbin, No I did not.

@hzhou
Copy link
Contributor Author

hzhou commented Jan 18, 2022

Minimum reproducer:
t.f90:

program t
    INTERFACE
        SUBROUTINE F_cdesc(buf) bind(C, name="F_cdesc")
            IMPLICIT NONE
            TYPE(*), DIMENSION(..), INTENT(in) :: buf
        END SUBROUTINE
    END INTERFACE

    ! Only reproducible with allocatable array. "INTEGER:: A(5)" would be fine.
    INTEGER, allocatable, dimension(:):: A

    allocate(A(5))

    A = (/1, 2, 3, 4, 5/);
    print*, A
    write(*,'(a6,2(i2))') 'Before:', lbound(A), ubound(A)

    ! "CALL F_cdesc(A)" would be fine
    CALL f08ts(A)

    print*, A
    write(*,'(a6,2(i2))') 'After:', lbound(A), ubound(A)

    deallocate(A)

    contains
    SUBROUTINE f08ts(buf)
        IMPLICIT NONE
        TYPE(*), DIMENSION(..), INTENT(in) :: buf
        call F_cdesc(buf)
    END SUBROUTINE f08ts
end

cdesc.c:

#include <ISO_Fortran_binding.h>

void F_cdesc(CFI_cdesc_t *buf) { }

Run with:

gcc-9 -c -o cdesc.o cdesc.c && gfortran-9 t.f90 cdesc.o -o t && ./t

Results:

           1           2           3           4           5
Before 1 5
           1           2           3           4           5
After: 0 4

EDIT: also reproduced with gcc 11

@hzhou
Copy link
Contributor Author

hzhou commented Jan 18, 2022

Reported to GCC bugzilla: https://gcc.gnu.org/bugzilla/show_bug.cgi?id=104100

@hzhou
Copy link
Contributor Author

hzhou commented Jan 18, 2022

It will be fixed in the upcoming gcc-12.0 according to the reply in the above bugzilla report.

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

No branches or pull requests

3 participants