PSyclone icon indicating copy to clipboard operation
PSyclone copied to clipboard

HoistLocalArraysTrans ignores size/bound changes on different calls to the same procedure

Open nmnobre opened this issue 2 years ago • 3 comments

In the following excerpt, kjcnt changes on different calls to mpp_init_readbot_strip causing a runtime error with the current implementation. In addition to checking if the array is allocated, we should also confirm its size needs no change.

module mppini
  real(kind=wp), allocatable, dimension(:,:), private :: zbot
  
  contains
  subroutine mpp_init_readbot_strip(kjstr, kjcnt, ldoce)
    integer, intent(in) :: kjstr
    integer, intent(in) :: kjcnt
    logical, dimension(jpiglo,kjcnt), intent(out) :: ldoce
    integer :: inumsave
    integer :: idx
    integer :: idx_1
    integer :: idx_2
    integer :: idx_3
    integer :: idx_4
    integer :: idx_5
    integer :: th_idx
    integer :: nthreads
    ! REAL(wp), DIMENSION(jpiglo,kjcnt) ::   zbot
  
    if (.NOT.ALLOCATED(zbot)) then
      ALLOCATE(zbot(1 : jpiglo, 1 : kjcnt))
    end if
  
    ...
  end subroutine mpp_init_readbot_strip
end module mppini

Cheers, -Nuno

nmnobre avatar Oct 04 '22 09:10 nmnobre

So the proper check could be something like:

    if (.NOT.ALLOCATED(X) .or. LBOUND(X,1) != L1D .or.  UBOUND(X,1) != U1D .or. LBOUND(X,2) != L2D .or. UBOUND(X,2) != U2D ) then
      if (ALLOCATED(X)) DEALLOCATE(X)
      ALLOCATE(X(L1D:U1D, L2D:U2D))
    end if

in the previous case with X=zbot, L1D=1, U1D=jpiglo, L2D=1, U2D=kjcnt

Would this work? Is it worth using this manually in mpp_init_readbot_strip first to check this works as expected?

sergisiso avatar Oct 04 '22 10:10 sergisiso

I can confirm that the suggested fix,

if (.NOT.ALLOCATED(zbot) .or. LBOUND(zbot,1) /= 1 .or.  UBOUND(zbot,1) /= jpiglo .or. LBOUND(zbot,2) /= 1 .or. UBOUND(zbot,2) /= kjcnt ) then
  if (ALLOCATED(zbot)) DEALLOCATE(zbot)
  ALLOCATE(zbot(1:jpiglo, 1:kjcnt))
end if

works as expected.

nmnobre avatar Oct 04 '22 11:10 nmnobre

We could perhaps tweak this in two ways:

  1. Using an if...else... to avoid the second conditional in case the array is not allocated yet.
  2. Check all the ubounds before any lbound since it's more likely they have changed, i.e. since it's more likely the array indexing still starts at one.

nmnobre avatar Oct 04 '22 11:10 nmnobre