| Class | sosi_dynamics_extarr | 
| In: | sosi/sosi_dynamics_extarr.F90 | 
| Subroutine : | |||
| x_SinLonS(0:imax-1) : | real(DP), intent(in ) | ||
| x_CosLonS(0:imax-1) : | real(DP), intent(in ) | ||
| x_SinLonN(0:imax-1) : | real(DP), intent(in ) | ||
| x_CosLonN(0:imax-1) : | real(DP), intent(in ) | ||
| xyzf_QMix(0:imax-1, 1:jmax, 1:kmax, 1:ncmax) : | real(DP), intent(in ) | ||
| PM : | real(DP), intent(in ) 
 | ||
| xyzf_ExtQMixS(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax) : | real(DP), intent(out) | ||
| xyzf_ExtQMixN(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax) : | real(DP), intent(out) | 
  subroutine SLTTExtArrExt2( x_SinLonS, x_CosLonS, x_SinLonN, x_CosLonN, xyzf_QMix, PM, xyzf_ExtQMixS, xyzf_ExtQMixN )
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    use mpi_wrapper, only : nprocs, myrank, MPIWrapperISend, MPIWrapperIRecv, MPIWrapperWait
    use sltt_const , only : iexmin, iexmax, jexmin, jexmax
    real(DP), intent(in ) :: x_SinLonS(0:imax-1)
    real(DP), intent(in ) :: x_CosLonS(0:imax-1)
    real(DP), intent(in ) :: x_SinLonN(0:imax-1)
    real(DP), intent(in ) :: x_CosLonN(0:imax-1)
    real(DP), intent(in ) :: xyzf_QMix(0:imax-1, 1:jmax, 1:kmax, 1:ncmax)
    real(DP), intent(in ) :: PM   ! 配列拡張する際、極ごえ後に符号が変わる場合は -1.0を与える。そうでない場合は1.0を与える。
    real(DP), intent(out) :: xyzf_ExtQMixS(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax)
    real(DP), intent(out) :: xyzf_ExtQMixN(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax)
    !
    ! local variables
    !
    real(DP) :: xzfy_SdRecvBuf(0:imax-1,1:kmax,1:ncmax,jexmin_min:jexmax_max)
    real(DP) :: xzfy_NdRecvBuf(0:imax-1,1:kmax,1:ncmax,jexmin_min:jexmax_max)
    real(DP) :: xzfy_dSendBuf (0:imax-1,1:kmax,1:ncmax,1         :jmax      )
    integer  :: ya_ExtSiReq   (jexmin_min:jexmax_max, 0:nprocs-1)
    integer  :: ya_ExtNiReq   (jexmin_min:jexmax_max, 0:nprocs-1)
    logical  :: ya_ExtSMPIFlag(jexmin_min:jexmax_max, 0:nprocs-1)
    logical  :: ya_ExtNMPiFlag(jexmin_min:jexmax_max, 0:nprocs-1)
    integer :: idep
    integer :: idest
    integer :: irank
    integer :: jlocal
    integer :: itag
    integer  :: i
    integer  :: j
    integer  :: k
    integer  :: n
    integer  :: ii
    ! 初期化確認
    ! Initialization check
    !
    if ( .not. sosi_dynamics_extarr_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if
    do j = 1, jmax
      xzfy_dSendBuf(:,:,:,j) = xyzf_QMix(:,j,:,:)
    end do
    do n = 0, nprocs-1
      if ( n == myrank ) then
        !
        ! Receive
        !
        do j = a_jexmin(n), a_jexmax(n)
          ! South array, mpi receive request
          irank  = ya_ExtRankOrgDataS       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataS(j,n)
          if ( irank == n ) then
            xzfy_SdRecvBuf(:,:,:,j) = xyzf_QMix(:,jlocal,:,:)
            ya_ExtSMPIFlag(j,n) = .false.
          else
            itag = MkSendRecvDestTag( n, j, 'S' )
            idep = irank
            call MPIWrapperIRecv( idep, imax, kmax, ncmax, xzfy_SdRecvBuf(:,:,:,j), ya_ExtSiReq(j,n), itag )
            ya_ExtSMPIFlag(j,n) = .true.
          end if
          ! North array, mpi receive request
          irank  = ya_ExtRankOrgDataN       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataN(j,n)
          if ( irank == n ) then
            xzfy_NdRecvBuf(:,:,:,j) = xyzf_QMix(:,jlocal,:,:)
            ya_ExtNMPIFlag(j,n) = .false.
          else
            itag = MkSendRecvDestTag( n, j, 'N' )
            idep = irank
            call MPIWrapperIRecv( idep, imax, kmax, ncmax, xzfy_NdRecvBuf(:,:,:,j), ya_ExtNiReq(j,n), itag )
            ya_ExtNMPIFlag(j,n) = .true.
          end if
        end do
      else
        !
        ! Send
        !
        do j = a_jexmin(n), a_jexmax(n)
          ! South array, mpi send request
          irank  = ya_ExtRankOrgDataS       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataS(j,n)
          if ( irank == myrank ) then
            itag = MkSendRecvDestTag( n, j, 'S' )
            idest = n
            call MPIWrapperISend( idest, imax, kmax, ncmax, xzfy_dSendBuf(:,:,:,jlocal), ya_ExtSiReq(j,n), itag )
            ya_ExtSMPIFlag(j,n) = .true.
          else
            ya_ExtSMPIFlag(j,n) = .false.
          end if
          ! North array, mpi send request
          irank  = ya_ExtRankOrgDataN       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataN(j,n)
          if ( irank == myrank ) then
            itag = MkSendRecvDestTag( n, j, 'N' )
            idest = n
            call MPIWrapperISend( idest, imax, kmax, ncmax, xzfy_dSendBuf(:,:,:,jlocal), ya_ExtNiReq(j,n), itag )
            ya_ExtNMPIFlag(j,n) = .true.
          else
            ya_ExtNMPIFlag(j,n) = .false.
          end if
        end do
      end if
    end do
    do n = 0, nprocs-1
      do j = a_jexmin(n), a_jexmax(n)
        ! South array
        if ( ya_ExtSMPIFlag(j,n) ) call MPIWrapperWait( ya_ExtSiReq(j,n) )
        ! North array
        if ( ya_ExtNMPIFlag(j,n) ) call MPIWrapperWait( ya_ExtNiReq(j,n) )
      end do
    end do
    do j = jexmin, jexmax
      if ( ya_ExtFlagCrossPoleOrgDataS(j,myrank) ) then
        do i = 0, imax/2-1
          ii = i + imax/2
          xyzf_ExtQMixS(ii,j,:,:) = PM * xzfy_SdRecvBuf(i,:,:,j)
        end do
        do i = imax/2, imax-1
          ii = i - imax/2
          xyzf_ExtQMixS(ii,j,:,:) = PM * xzfy_SdRecvBuf(i,:,:,j)
        end do
      else
        xyzf_ExtQMixS(0:imax-1,j,:,:) = xzfy_SdRecvBuf(:,:,:,j)
      end if
      if ( ya_ExtFlagCrossPoleOrgDataN(j,myrank) ) then
        do i = 0, imax/2-1
          ii = i + imax/2
          xyzf_ExtQMixN(ii,j,:,:) = PM * xzfy_NdRecvBuf(i,:,:,j)
        end do
        do i = imax/2, imax-1
          ii = i - imax/2
          xyzf_ExtQMixN(ii,j,:,:) = PM * xzfy_NdRecvBuf(i,:,:,j)
        end do
      else
        xyzf_ExtQMixN(0:imax-1,j,:,:) = xzfy_NdRecvBuf(:,:,:,j)
      end if
    end do
    !===========================================
    ! set values at longitudinal edge
    !-------------------------------------------
#if defined(AXISYMMETRY) || defined(AXISYMMETRY_SJPACK)
    do n = 1, ncmax
      do k = 1, kmax
        do j = jexmin, jexmax
          do i = iexmin, 0-1
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(0,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(0,j,k,n)
          end do
          do i = imax-1+1, iexmax
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(0,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(0,j,k,n)
          end do
        end do
      end do
    end do
#else
    do n = 1, ncmax
      do k = 1, kmax
        do j = jexmin, jexmax
          do i = iexmin, 0-1
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(imax+i,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(imax+i,j,k,n)
          end do
          do i = imax-1+1, iexmax
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(i-imax,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(i-imax,j,k,n)
          end do
        end do
      end do
    end do
#endif
  end subroutine SLTTExtArrExt2
          | Subroutine : | |
| x_LonS( 0:imax-1 ) : | real(DP), intent(in ) | 
| y_LatS( 1:jmax/2 ) : | real(DP), intent(in ) | 
| x_LonN( 0:imax-1 ) : | real(DP), intent(in ) | 
| y_LatN( 1:jmax/2 ) : | real(DP), intent(in ) | 
| x_ExtLonS(iexmin:iexmax) : | real(DP), intent(out) | 
| y_ExtLatS(jexmin:jexmax) : | real(DP), intent(out) | 
| x_ExtLonN(iexmin:iexmax) : | real(DP), intent(out) | 
| y_ExtLatN(jexmin:jexmax) : | real(DP), intent(out) | 
  subroutine SLTTExtArrInit( x_LonS, y_LatS, x_LonN, y_LatN, x_ExtLonS, y_ExtLatS, x_ExtLonN, y_ExtLatN )
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    !
    ! MPI
    !
    use mpi_wrapper   , only : myrank, nprocs, MPIWrapperISend, MPIWrapperIRecv, MPIWrapperWait
    use constants0, only : PI
    use axesset   , only : y_Lat
    use sltt_const, only : PIx2, PIH, iexmin, iexmax, jexmin, jexmax, jmaxh, jexglobalmin, jexglobalmax
    real(DP), intent(in ) :: x_LonS   (   0:imax-1    )
    real(DP), intent(in ) :: y_LatS   (   1:jmax/2    )
    real(DP), intent(in ) :: x_LonN   (   0:imax-1    )
    real(DP), intent(in ) :: y_LatN   (   1:jmax/2    )
    real(DP), intent(out) :: x_ExtLonS(iexmin:iexmax)
    real(DP), intent(out) :: y_ExtLatS(jexmin:jexmax)
    real(DP), intent(out) :: x_ExtLonN(iexmin:iexmax)
    real(DP), intent(out) :: y_ExtLatN(jexmin:jexmax)
    !
    ! local variables
    !
    integer  :: idest
    integer  :: idep
    ! y_LatGlobal(jj)         : latitude of jj-th component
    !                         : (latitudinal global array)
    ! y_RankGlobal(jj)        : rank of jj-th component
    !                         : (latitudinal global array)
    ! y_JLocalIndexGlobal(jj) : j index for each rank array of jj-th component
    !                         : (latitudinal global array)
    real(DP) :: y_LatGlobal        (1:jmax_global)
    integer  :: y_RankGlobal       (1:jmax_global)
    integer  :: y_JLocalIndexGlobal(1:jmax_global)
    real(DP) :: y_ExtLatGlobal          (jexglobalmin:jexglobalmax)
    integer  :: y_ExtRankGlobal         (jexglobalmin:jexglobalmax)
    integer  :: y_ExtJLocalIndexGlobal  (jexglobalmin:jexglobalmax)
    logical  :: y_ExtFlagCrossPoleGlobal(jexglobalmin:jexglobalmax)
    integer, allocatable :: ya_ExtSiReq   (:,:)
    integer, allocatable :: ya_ExtNiReq   (:,:)
    logical, allocatable :: ya_ExtSMPIFlag(:,:)
    logical, allocatable :: ya_ExtNMPiFlag(:,:)
    real(DP), allocatable :: aa_SdRecvBuf(:,:)
    real(DP), allocatable :: aa_NdRecvBuf(:,:)
    real(DP)              :: aa_dSendBuf (1,1:jmax)
    integer  :: irank
    integer  :: jlocal
    integer  :: itag
    integer  :: j1
    integer  :: i
    integer  :: j
    integer  :: n
    if ( sosi_dynamics_extarr_inited ) return
#if defined(AXISYMMETRY) || defined(AXISYMMETRY_SJPACK)
    do i = iexmin, 0-1
      x_ExtLonS(i) = x_LonS(0) - ( PIx2 * dble(-i) - x_LonS(0) )
    end do
    do i = 0, imax-1
      x_ExtLonS(i) = x_LonS(i)
    end do
    x_ExtLonS(imax-1+1) = PIx2
    do i = imax-1+1+1, iexmax
      x_ExtLonS(i) = PIx2 + ( PIx2 * dble(i-(imax-1+1)) - x_LonS(0) )
    end do
    !
    do i = iexmin, 0-1
      x_ExtLonN(i) = x_LonN(0) - ( PIx2 * dble(-i) - x_LonN(0) )
    end do
    do i = 0, imax-1
      x_ExtLonN(i) = x_LonN(i)
    end do
    x_ExtLonN(imax-1+1) = PIx2
    do i = imax-1+1+1, iexmax
      x_ExtLonN(i) = PIx2 + ( PIx2 * dble(i-(imax-1+1)) - x_LonN(0) )
    end do
#else
    do i = iexmin, 0-1
      x_ExtLonS(i) = x_LonS(0) - ( x_LonS(-i) - x_LonS(0) )
    end do
    do i = 0, imax-1
      x_ExtLonS(i) = x_LonS(i)
    end do
    x_ExtLonS(imax-1+1) = PIx2
    do i = imax-1+1+1, iexmax
      x_ExtLonS(i) = PIx2 + ( x_LonS(i-(imax-1+1)) - x_LonS(0) )
    end do
    !
    do i = iexmin, 0-1
      x_ExtLonN(i) = x_LonN(0) - ( x_LonN(-i) - x_LonN(0) )
    end do
    do i = 0, imax-1
      x_ExtLonN(i) = x_LonN(i)
    end do
    x_ExtLonN(imax-1+1) = PIx2
    do i = imax-1+1+1, iexmax
      x_ExtLonN(i) = PIx2 + ( x_LonN(i-(imax-1+1)) - x_LonN(0) )
    end do
#endif
    !====================================================================
    ! Set y_LatGlobal, y_RankGlobal, y_JLocalIndexGlobal.
    call SLTTExtArrPrepGlobalArray( y_LatGlobal, y_RankGlobal, y_JLocalIndexGlobal )
    do j = jexglobalmin, 1-1
      y_ExtLatGlobal          (j) = -PIH - ( y_LatGlobal(1-j) - ( -PIH ) )
      y_ExtRankGlobal         (j) = y_RankGlobal       (1-j)
      y_ExtJLocalIndexGlobal  (j) = y_JLocalIndexGlobal(1-j)
      y_ExtFlagCrossPoleGlobal(j) = .true.
    end do
    do j = 1, jmax_global
      y_ExtLatGlobal          (j) = y_LatGlobal        (j)
      y_ExtRankGlobal         (j) = y_RankGlobal       (j)
      y_ExtJLocalIndexGlobal  (j) = y_JLocalIndexGlobal(j)
      y_ExtFlagCrossPoleGlobal(j) = .false.
    end do
    do j = jmax_global+1, jexglobalmax
      y_ExtLatGlobal          (j) =  PIH + ( PIH - y_LatGlobal(jmax_global-(j-(jmax_global+1))) )
      y_ExtRankGlobal         (j) = y_RankGlobal       (jmax_global-(j-(jmax_global+1)))
      y_ExtJLocalIndexGlobal  (j) = y_JLocalIndexGlobal(jmax_global-(j-(jmax_global+1)))
      y_ExtFlagCrossPoleGlobal(j) = .true.
    end do
    allocate( ya_ExtLatS                 (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtLatN                 (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtRankOrgDataS         (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtRankOrgDataN         (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtJLocalIndexOrgDataS  (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtJLocalIndexOrgDataN  (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtFlagCrossPoleOrgDataS(jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtFlagCrossPoleOrgDataN(jexmin_min:jexmax_max, 0:nprocs-1) )
    ! rank, (local) j index, and flag for cross pole for original data transfered to rank n
    do n = 0, nprocs-1
      ! j1: global j index for (local) j = 1 at rank n (first index of south array)
      j1 = FindGlobalJIndex( y_RankGlobal, y_JLocalIndexGlobal, n, 1 )
      do j = a_jexmin(n), a_jexmax(n)
        ya_ExtLatS                 (j,n) = y_ExtLatGlobal          (j1+j-1)
        ya_ExtRankOrgDataS         (j,n) = y_ExtRankGlobal         (j1+j-1)
        ya_ExtJLocalIndexOrgDataS  (j,n) = y_ExtJLocalIndexGlobal  (j1+j-1)
        ya_ExtFlagCrossPoleOrgDataS(j,n) = y_ExtFlagCrossPoleGlobal(j1+j-1)
      end do
      ! j1: global j index for (local) j = 1 at rank n (first index of north array)
      j1 = FindGlobalJIndex( y_RankGlobal, y_JLocalIndexGlobal, n, a_jmaxs(n)+1 )
      do j = a_jexmin(n), a_jexmax(n)
        ya_ExtLatN                 (j,n) = y_ExtLatGlobal          (j1+j-1)
        ya_ExtRankOrgDataN         (j,n) = y_ExtRankGlobal         (j1+j-1)
        ya_ExtJLocalIndexOrgDataN  (j,n) = y_ExtJLocalIndexGlobal  (j1+j-1)
        ya_ExtFlagCrossPoleOrgDataN(j,n) = y_ExtFlagCrossPoleGlobal(j1+j-1)
      end do
    end do
    allocate( aa_SdRecvBuf(1,jexmin_min:jexmax_max) )
    allocate( aa_NdRecvBuf(1,jexmin_min:jexmax_max) )
    allocate( ya_ExtSiReq   (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtNiReq   (jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtSMPIFlag(jexmin_min:jexmax_max, 0:nprocs-1) )
    allocate( ya_ExtNMPiFlag(jexmin_min:jexmax_max, 0:nprocs-1) )
    !
    ! Test transfer of data (latitude)
    ! Transfered data for test are y_ExtLatS and y_ExtLatN.
    ! But, those data will be overwritten after the test transfer.
    !
    aa_dSendBuf(1,:) = y_Lat
    do n = 0, nprocs-1
      if ( n == myrank ) then
        !
        ! Receive
        !
        do j = a_jexmin(n), a_jexmax(n)
          ! South array, mpi receive request
          irank  = ya_ExtRankOrgDataS       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataS(j,n)
          if ( irank == n ) then
            aa_SdRecvBuf(:,j) = y_Lat(jlocal)
            ya_ExtSMPIFlag(j,n) = .false.
          else
            itag = MkSendRecvDestTag( n, j, 'S' )
            idep = irank
            call MPIWrapperIRecv( idep, 1, aa_SdRecvBuf(:,j), ya_ExtSiReq(j,n), itag )
            ya_ExtSMPIFlag(j,n) = .true.
          end if
          ! North array, mpi receive request
          irank  = ya_ExtRankOrgDataN       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataN(j,n)
          if ( irank == n ) then
            aa_NdRecvBuf(:,j) = y_Lat(jlocal)
            ya_ExtNMPIFlag(j,n) = .false.
          else
            itag = MkSendRecvDestTag( n, j, 'N' )
            idep = irank
            call MPIWrapperIRecv( idep, 1, aa_NdRecvBuf(:,j), ya_ExtNiReq(j,n), itag )
            ya_ExtNMPIFlag(j,n) = .true.
          end if
        end do
      else
        !
        ! Send
        !
        do j = a_jexmin(n), a_jexmax(n)
          ! South array, mpi send request
          irank  = ya_ExtRankOrgDataS       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataS(j,n)
          if ( irank == myrank ) then
            itag = MkSendRecvDestTag( n, j, 'S' )
            idest = n
            call MPIWrapperISend( idest, 1, aa_dSendBuf(:,jlocal), ya_ExtSiReq(j,n), itag )
            ya_ExtSMPIFlag(j,n) = .true.
          else
            ya_ExtSMPIFlag(j,n) = .false.
          end if
          ! North array, mpi send request
          irank  = ya_ExtRankOrgDataN       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataN(j,n)
          if ( irank == myrank ) then
            itag = MkSendRecvDestTag( n, j, 'N' )
            idest = n
            call MPIWrapperISend( idest, 1, aa_dSendBuf(:,jlocal), ya_ExtNiReq(j,n), itag )
            ya_ExtNMPIFlag(j,n) = .true.
          else
            ya_ExtNMPIFlag(j,n) = .false.
          end if
        end do
      end if
    end do
    do n = 0, nprocs-1
      do j = a_jexmin(n), a_jexmax(n)
        ! South array
        if ( ya_ExtSMPIFlag(j,n) ) call MPIWrapperWait( ya_ExtSiReq(j,n) )
        ! North array
        if ( ya_ExtNMPIFlag(j,n) ) call MPIWrapperWait( ya_ExtNiReq(j,n) )
      end do
    end do
    y_ExtLatS = aa_SdRecvBuf(1,jexmin:jexmax)
    y_ExtLatN = aa_NdRecvBuf(1,jexmin:jexmax)
    deallocate( aa_SdRecvBuf )
    deallocate( aa_NdRecvBuf )
    deallocate( ya_ExtSiReq    )
    deallocate( ya_ExtNiReq    )
    deallocate( ya_ExtSMPIFlag )
    deallocate( ya_ExtNMPiFlag )
    ! Checking transfer
    !
    ! global j index for j = 1 at local (southern edge of south array)
    j1 = FindGlobalJIndex( y_RankGlobal, y_JLocalIndexGlobal, myrank, 1 )
    do j = jexmin, jexmax
      if ( .not. ya_ExtFlagCrossPoleOrgDataS(j,myrank) ) then
        if ( y_ExtLatS(j) /= y_ExtLatGlobal(j1+j-1) ) then
          call MessageNotify( 'E', module_name, 'Rank = %d, y_ExtLatS(%d) = %f is not same as y_ExtLatGlobal(%d) = %f.', i = (/myrank, j, j1+j-1/), d = (/y_ExtLatS(j), y_ExtLatGlobal(j1+j-1)/) )
        end if
      end if
    end do
    ! global j index for j = 1 at local (southern edge of north array)
    j1 = FindGlobalJIndex( y_RankGlobal, y_JLocalIndexGlobal, myrank, jmaxh+1 )
    do j = jexmin, jexmax
      if ( .not. ya_ExtFlagCrossPoleOrgDataN(j,myrank) ) then
        if ( y_ExtLatN(j) /= y_ExtLatGlobal(j1+j-1) ) then
          call MessageNotify( 'E', module_name, 'Rank = %d, y_ExtLatN(%d) = %f is not same as y_ExtLatGlobal(%d) = %f.', i = (/myrank, j, j1+j-1/), d = (/y_ExtLatN(j), y_ExtLatGlobal(j1+j-1)/) )
        end if
      end if
    end do
    ! global j index for j = 1 at local (southern edge of south array)
    j1 = FindGlobalJIndex( y_RankGlobal, y_JLocalIndexGlobal, myrank, 1 )
    do j = jexmin, jexmax
      y_ExtLatS(j) = y_ExtLatGlobal(j1+j-1)
    end do
    ! global j index for j = 1 at local (southern edge of north array)
    j1 = FindGlobalJIndex( y_RankGlobal, y_JLocalIndexGlobal, myrank, jmaxh+1 )
    do j = jexmin, jexmax
      y_ExtLatN(j) = y_ExtLatGlobal(j1+j-1)
    end do
    sosi_dynamics_extarr_inited = .true.
  end subroutine SLTTExtArrInit
          | Function : | |
| j : | integer | 
| y_RankGlobal(1:jmax_global) : | integer, intent(in) | 
| y_JLocalIndexGlobal(1:jmax_global) : | integer, intent(in) | 
| irank : | integer, intent(in) | 
| jindex : | integer, intent(in) | 
  function FindGlobalJIndex( y_RankGlobal, y_JLocalIndexGlobal, irank, jindex ) result( j )
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    integer, intent(in) :: y_RankGlobal       (1:jmax_global)
    integer, intent(in) :: y_JLocalIndexGlobal(1:jmax_global)
    integer, intent(in) :: irank
    integer, intent(in) :: jindex
    integer :: j
    search_j: do j = 1, jmax_global
      if ( ( y_RankGlobal(j) == irank ) .and. ( y_JLocalIndexGlobal(j) == jindex ) ) exit search_j
    end do search_j
    if ( j > jmax_global ) then
      call MessageNotify( 'E', module_name, 'Cannot find proper j.' )
    end if
  end function FindGlobalJIndex
          | Function : | |
| itag : | integer | 
| irank : | integer , intent(in) | 
| jlocalindex : | integer , intent(in) | 
| hemisphere : | character(*), intent(in) | 
  function MkSendRecvDestTag( irank, jlocalindex, hemisphere ) result( itag )
    ! irank       : rank at destination
    ! jlocalindex : j index at destination
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    !
    ! MPI
    !
    use mpi_wrapper, only : nprocs
    integer     , intent(in) :: irank
    integer     , intent(in) :: jlocalindex
    character(*), intent(in) :: hemisphere
    integer :: itag
    integer :: ioffset
    ioffset = jmax_max * nprocs
    select case( hemisphere )
    case ( 'S' )
      itag = jmax_max * (        + irank ) + jlocalindex
    case ( 'N' )
      itag = jmax_max * ( nprocs + irank ) + jlocalindex
    case default
      call MessageNotify( 'E', module_name, 'Unexpected value of hemisphere.' )
    end select
    itag = itag + ioffset
    if ( itag < 0 ) then
      write( 6, * ) 'tag is negative:', itag, irank, jlocalindex, trim( hemisphere )
    end if
  end function MkSendRecvDestTag
          | Subroutine : | |
| y_ExtLatS(jexmin:jexmax) : | real(DP), intent(in ) | 
| y_ExtLatN(jexmin:jexmax) : | real(DP), intent(in ) | 
| x_SinLonS(0:imax-1) : | real(DP), intent(in ) | 
| x_CosLonS(0:imax-1) : | real(DP), intent(in ) | 
| x_SinLonN(0:imax-1) : | real(DP), intent(in ) | 
| x_CosLonN(0:imax-1) : | real(DP), intent(in ) | 
| xyzf_QMix(0:imax-1, 1:jmax, 1:kmax, 1:ncmax) : | real(DP), intent(in ) | 
| xyz_U(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in ) | 
| xyz_V(0:imax-1, 1:jmax, 1:kmax) : | real(DP), intent(in ) | 
| xyzf_ExtDQMixDLatS(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax) : | real(DP), intent(in ) | 
| xyzf_ExtDQMixDLatN(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax) : | real(DP), intent(in ) | 
| xyzf_ExtQMixS(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax) : | real(DP), intent(out) | 
| xyzf_ExtQMixN(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax) : | real(DP), intent(out) | 
| xyz_ExtUS(iexmin:iexmax, jexmin:jexmax, 1:kmax) : | real(DP), intent(out) | 
| xyz_ExtUN(iexmin:iexmax, jexmin:jexmax, 1:kmax) : | real(DP), intent(out) | 
| xyz_ExtVS(iexmin:iexmax, jexmin:jexmax, 1:kmax) : | real(DP), intent(out) | 
| xyz_ExtVN(iexmin:iexmax, jexmin:jexmax, 1:kmax) : | real(DP), intent(out) | 
  subroutine SLTTExtArrExt( y_ExtLatS, y_ExtLatN, x_SinLonS, x_CosLonS, x_SinLonN, x_CosLonN, xyzf_QMix, xyz_U, xyz_V, xyzf_ExtDQMixDLatS, xyzf_ExtDQMixDLatN, xyzf_ExtQMixS, xyzf_ExtQMixN, xyz_ExtUS, xyz_ExtUN, xyz_ExtVS, xyz_ExtVN )
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    use mpi_wrapper, only : nprocs, myrank, MPIWrapperISend, MPIWrapperIRecv, MPIWrapperWait
    use sltt_const , only : iexmin, iexmax, jexmin, jexmax
!!$    use sltt_lagint, only : SLTTIrrHerIntQui1DNonUni
    real(DP), intent(in ) :: y_ExtLatS(jexmin:jexmax)
    real(DP), intent(in ) :: y_ExtLatN(jexmin:jexmax)
    real(DP), intent(in ) :: x_SinLonS(0:imax-1)
    real(DP), intent(in ) :: x_CosLonS(0:imax-1)
    real(DP), intent(in ) :: x_SinLonN(0:imax-1)
    real(DP), intent(in ) :: x_CosLonN(0:imax-1)
    real(DP), intent(in ) :: xyzf_QMix(0:imax-1, 1:jmax, 1:kmax, 1:ncmax)
    real(DP), intent(in ) :: xyz_U    (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyz_V    (0:imax-1, 1:jmax, 1:kmax)
    real(DP), intent(in ) :: xyzf_ExtDQMixDLatS(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax)
    real(DP), intent(in ) :: xyzf_ExtDQMixDLatN(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax)
    real(DP), intent(out) :: xyzf_ExtQMixS(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax)
    real(DP), intent(out) :: xyzf_ExtQMixN(iexmin:iexmax, jexmin:jexmax, 1:kmax, 1:ncmax)
    real(DP), intent(out) :: xyz_ExtUS    (iexmin:iexmax, jexmin:jexmax, 1:kmax)
    real(DP), intent(out) :: xyz_ExtUN    (iexmin:iexmax, jexmin:jexmax, 1:kmax)
    real(DP), intent(out) :: xyz_ExtVS    (iexmin:iexmax, jexmin:jexmax, 1:kmax)
    real(DP), intent(out) :: xyz_ExtVN    (iexmin:iexmax, jexmin:jexmax, 1:kmax)
    !
    ! local variables
    !
    real(DP) :: xzfy_SdRecvBuf(0:imax-1,1:kmax,1:ncmax+2,jexmin_min:jexmax_max)
    real(DP) :: xzfy_NdRecvBuf(0:imax-1,1:kmax,1:ncmax+2,jexmin_min:jexmax_max)
    real(DP) :: xzfy_dSendBuf (0:imax-1,1:kmax,1:ncmax+2,1         :jmax      )
    integer  :: ya_ExtSiReq   (jexmin_min:jexmax_max, 0:nprocs-1)
    integer  :: ya_ExtNiReq   (jexmin_min:jexmax_max, 0:nprocs-1)
    logical  :: ya_ExtSMPIFlag(jexmin_min:jexmax_max, 0:nprocs-1)
    logical  :: ya_ExtNMPiFlag(jexmin_min:jexmax_max, 0:nprocs-1)
    integer :: idep
    integer :: idest
    integer :: irank
    integer :: jlocal
    integer :: itag
    integer  :: i
    integer  :: j
    integer  :: k
    integer  :: n
    integer  :: ii
    ! 初期化確認
    ! Initialization check
    !
    if ( .not. sosi_dynamics_extarr_inited ) then
      call MessageNotify( 'E', module_name, 'This module has not been initialized.' )
    end if
    do j = 1, jmax
      xzfy_dSendBuf(:,:,1:ncmax  ,j) = xyzf_QMix(:,j,:,:)
      xzfy_dSendBuf(:,:,  ncmax+1,j) = xyz_U(:,j,:)
      xzfy_dSendBuf(:,:,  ncmax+2,j) = xyz_V(:,j,:)
    end do
    do n = 0, nprocs-1
      if ( n == myrank ) then
        !
        ! Receive
        !
        do j = a_jexmin(n), a_jexmax(n)
          ! South array, mpi receive request
          irank  = ya_ExtRankOrgDataS       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataS(j,n)
          if ( irank == n ) then
            xzfy_SdRecvBuf(:,:,1:ncmax  ,j) = xyzf_QMix(:,jlocal,:,:)
            xzfy_SdRecvBuf(:,:,  ncmax+1,j) = xyz_U(:,jlocal,:)
            xzfy_SdRecvBuf(:,:,  ncmax+2,j) = xyz_V(:,jlocal,:)
            ya_ExtSMPIFlag(j,n) = .false.
          else
            itag = MkSendRecvDestTag( n, j, 'S' )
            idep = irank
            call MPIWrapperIRecv( idep, imax, kmax, ncmax+2, xzfy_SdRecvBuf(:,:,:,j), ya_ExtSiReq(j,n), itag )
            ya_ExtSMPIFlag(j,n) = .true.
          end if
          ! North array, mpi receive request
          irank  = ya_ExtRankOrgDataN       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataN(j,n)
          if ( irank == n ) then
            xzfy_NdRecvBuf(:,:,1:ncmax  ,j) = xyzf_QMix(:,jlocal,:,:)
            xzfy_NdRecvBuf(:,:,  ncmax+1,j) = xyz_U(:,jlocal,:)
            xzfy_NdRecvBuf(:,:,  ncmax+2,j) = xyz_V(:,jlocal,:)
            ya_ExtNMPIFlag(j,n) = .false.
          else
            itag = MkSendRecvDestTag( n, j, 'N' )
            idep = irank
            call MPIWrapperIRecv( idep, imax, kmax, ncmax+2, xzfy_NdRecvBuf(:,:,:,j), ya_ExtNiReq(j,n), itag )
            ya_ExtNMPIFlag(j,n) = .true.
          end if
        end do
      else
        !
        ! Send
        !
        do j = a_jexmin(n), a_jexmax(n)
          ! South array, mpi send request
          irank  = ya_ExtRankOrgDataS       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataS(j,n)
          if ( irank == myrank ) then
            itag = MkSendRecvDestTag( n, j, 'S' )
            idest = n
            call MPIWrapperISend( idest, imax, kmax, ncmax+2, xzfy_dSendBuf(:,:,:,jlocal), ya_ExtSiReq(j,n), itag )
            ya_ExtSMPIFlag(j,n) = .true.
          else
            ya_ExtSMPIFlag(j,n) = .false.
          end if
          ! North array, mpi send request
          irank  = ya_ExtRankOrgDataN       (j,n)
          jlocal = ya_ExtJLocalIndexOrgDataN(j,n)
          if ( irank == myrank ) then
            itag = MkSendRecvDestTag( n, j, 'N' )
            idest = n
            call MPIWrapperISend( idest, imax, kmax, ncmax+2, xzfy_dSendBuf(:,:,:,jlocal), ya_ExtNiReq(j,n), itag )
            ya_ExtNMPIFlag(j,n) = .true.
          else
            ya_ExtNMPIFlag(j,n) = .false.
          end if
        end do
      end if
    end do
    do n = 0, nprocs-1
      do j = a_jexmin(n), a_jexmax(n)
        ! South array
        if ( ya_ExtSMPIFlag(j,n) ) call MPIWrapperWait( ya_ExtSiReq(j,n) )
        ! North array
        if ( ya_ExtNMPIFlag(j,n) ) call MPIWrapperWait( ya_ExtNiReq(j,n) )
      end do
    end do
    do j = jexmin, jexmax
      if ( ya_ExtFlagCrossPoleOrgDataS(j,myrank) ) then
        do i = 0, imax/2-1
          ii = i + imax/2
          xyzf_ExtQMixS(ii,j,:,:) =   xzfy_SdRecvBuf(i,:,1:ncmax  ,j)
          xyz_ExtUS    (ii,j,:)   = - xzfy_SdRecvBuf(i,:,  ncmax+1,j)
          xyz_ExtVS    (ii,j,:)   = - xzfy_SdRecvBuf(i,:,  ncmax+2,j)
        end do
        do i = imax/2, imax-1
          ii = i - imax/2
          xyzf_ExtQMixS(ii,j,:,:) =   xzfy_SdRecvBuf(i,:,1:ncmax  ,j)
          xyz_ExtUS    (ii,j,:)   = - xzfy_SdRecvBuf(i,:,  ncmax+1,j)
          xyz_ExtVS    (ii,j,:)   = - xzfy_SdRecvBuf(i,:,  ncmax+2,j)
        end do
      else
        xyzf_ExtQMixS(0:imax-1,j,:,:) = xzfy_SdRecvBuf(:,:,1:ncmax  ,j)
        xyz_ExtUS    (0:imax-1,j,:)   = xzfy_SdRecvBuf(:,:,  ncmax+1,j)
        xyz_ExtVS    (0:imax-1,j,:)   = xzfy_SdRecvBuf(:,:,  ncmax+2,j)
      end if
      if ( ya_ExtFlagCrossPoleOrgDataN(j,myrank) ) then
        do i = 0, imax/2-1
          ii = i + imax/2
          xyzf_ExtQMixN(ii,j,:,:) =   xzfy_NdRecvBuf(i,:,1:ncmax  ,j)
          xyz_ExtUN    (ii,j,:)   = - xzfy_NdRecvBuf(i,:,  ncmax+1,j)
          xyz_ExtVN    (ii,j,:)   = - xzfy_NdRecvBuf(i,:,  ncmax+2,j)
        end do
        do i = imax/2, imax-1
          ii = i - imax/2
          xyzf_ExtQMixN(ii,j,:,:) =   xzfy_NdRecvBuf(i,:,1:ncmax  ,j)
          xyz_ExtUN    (ii,j,:)   = - xzfy_NdRecvBuf(i,:,  ncmax+1,j)
          xyz_ExtVN    (ii,j,:)   = - xzfy_NdRecvBuf(i,:,  ncmax+2,j)
        end do
      else
        xyzf_ExtQMixN(0:imax-1,j,:,:) = xzfy_NdRecvBuf(:,:,1:ncmax  ,j)
        xyz_ExtUN    (0:imax-1,j,:)   = xzfy_NdRecvBuf(:,:,  ncmax+1,j)
        xyz_ExtVN    (0:imax-1,j,:)   = xzfy_NdRecvBuf(:,:,  ncmax+2,j)
      end if
    end do
    !===========================================
    ! set values at longitudinal edge
    !-------------------------------------------
#if defined(AXISYMMETRY) || defined(AXISYMMETRY_SJPACK)
    do n = 1, ncmax
      do k = 1, kmax
        do j = jexmin, jexmax
          do i = iexmin, 0-1
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(0,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(0,j,k,n)
          end do
          do i = imax-1+1, iexmax
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(0,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(0,j,k,n)
          end do
        end do
      end do
    end do
    do k = 1, kmax
      do j = jexmin, jexmax
        do i = iexmin, 0-1
          xyz_ExtUS(i,j,k) = xyz_ExtUS(0,j,k)
          xyz_ExtUN(i,j,k) = xyz_ExtUN(0,j,k)
          xyz_ExtVS(i,j,k) = xyz_ExtVS(0,j,k)
          xyz_ExtVN(i,j,k) = xyz_ExtVN(0,j,k)
        end do
        do i = imax-1+1, iexmax
          xyz_ExtUS(i,j,k) = xyz_ExtUS(0,j,k)
          xyz_ExtUN(i,j,k) = xyz_ExtUN(0,j,k)
          xyz_ExtVS(i,j,k) = xyz_ExtVS(0,j,k)
          xyz_ExtVN(i,j,k) = xyz_ExtVN(0,j,k)
        end do
      end do
    end do
#else
    do n = 1, ncmax
      do k = 1, kmax
        do j = jexmin, jexmax
          do i = iexmin, 0-1
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(imax+i,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(imax+i,j,k,n)
          end do
          do i = imax-1+1, iexmax
            xyzf_ExtQMixS(i,j,k,n) = xyzf_ExtQMixS(i-imax,j,k,n)
            xyzf_ExtQMixN(i,j,k,n) = xyzf_ExtQMixN(i-imax,j,k,n)
          end do
        end do
      end do
    end do
    do k = 1, kmax
      do j = jexmin, jexmax
        do i = iexmin, 0-1
          xyz_ExtUS(i,j,k) = xyz_ExtUS(imax+i,j,k)
          xyz_ExtUN(i,j,k) = xyz_ExtUN(imax+i,j,k)
          xyz_ExtVS(i,j,k) = xyz_ExtVS(imax+i,j,k)
          xyz_ExtVN(i,j,k) = xyz_ExtVN(imax+i,j,k)
        end do
        do i = imax-1+1, iexmax
          xyz_ExtUS(i,j,k) = xyz_ExtUS(i-imax,j,k)
          xyz_ExtUN(i,j,k) = xyz_ExtUN(i-imax,j,k)
          xyz_ExtVS(i,j,k) = xyz_ExtVS(i-imax,j,k)
          xyz_ExtVN(i,j,k) = xyz_ExtVN(i-imax,j,k)
        end do
      end do
    end do
#endif
  end subroutine SLTTExtArrExt
          | Subroutine : | |
| y_LatGlobal(1:jmax_global) : | real(DP), intent(out) | 
| y_RankGlobal(1:jmax_global) : | integer , intent(out) | 
| y_JLocalIndexGlobal(1:jmax_global) : | integer , intent(out) | 
  subroutine SLTTExtArrPrepGlobalArray( y_LatGlobal, y_RankGlobal, y_JLocalIndexGlobal )
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    !
    ! MPI
    !
    use mpi_wrapper, only : nprocs, myrank, MPIWrapperISend, MPIWrapperIRecv, MPIWrapperWait
    use axesset   , only : y_Lat
    use sltt_const, only : jexmin, jexmax
    real(DP), intent(out) :: y_LatGlobal        (1:jmax_global)
    integer , intent(out) :: y_RankGlobal       (1:jmax_global)
    integer , intent(out) :: y_JLocalIndexGlobal(1:jmax_global)
    !
    ! local variables
    !
    real(DP)              :: a_dSendBuf (1:jmax)
    real(DP), allocatable :: aa_dRecvBuf(:,:)
    integer  :: a_iReqSend (0:nprocs-1)
    integer  :: a_iReqRecv (0:nprocs-1)
    integer              :: idep
    integer              :: idest
    integer              :: j
    integer              :: jglobal
    integer              :: n
    allocate( a_jmax (0:nprocs-1) )
    allocate( a_jmaxs(0:nprocs-1) )
    allocate( a_jmaxn(0:nprocs-1) )
    allocate( a_jexmin(0:nprocs-1) )
    allocate( a_jexmax(0:nprocs-1) )
    ! share jmax
    call SLTTExtShareiScalar( jmax, a_jmax )
    ! share jexmin
    call SLTTExtShareiScalar( jexmin, a_jexmin )
    ! share jexmax
    call SLTTExtShareiScalar( jexmax, a_jexmax )
    jmax_max = 1
    do n = 0, nprocs-1
      jmax_max = max( jmax_max, a_jmax(n) )
    end do
    a_jmaxs = a_jmax / 2
    a_jmaxn = a_jmax / 2
    jexmin_min = 1
    do n = 0, nprocs-1
      jexmin_min = min( jexmin_min, a_jexmin(n) )
    end do
    jexmax_max = 1
    do n = 0, nprocs-1
      jexmax_max = max( jexmax_max, a_jexmax(n) )
    end do
    !----------
    allocate( aa_dRecvBuf(1:jmax_max, 0:nprocs-1) )
    a_dSendBuf = y_Lat
    do n = 0, nprocs-1
      if ( n == myrank ) cycle
      idest = n
      call MPIWrapperISend( idest, jmax, a_dSendBuf, a_IReqSend(n) )
    end do
    do n = 0, nprocs-1
      if ( n == myrank ) cycle
      idep = n
      call MPIWrapperIRecv( idep, a_jmax(n), aa_dRecvBuf(1:a_jmax(n), n), a_IReqRecv(n) )
    end do
    n = myrank
    aa_dRecvBuf(1:jmax, n) = a_dSendBuf
    do n = 0, nprocs-1
      if ( n == myrank ) cycle
      call MPIWrapperWait( a_IReqSend(n) )
      call MPIWrapperWait( a_IReqRecv(n) )
    end do
    !
    jglobal = 0
    do n = nprocs-1, 0, -1
      do j = 1, a_jmaxs(n)
        jglobal = jglobal + 1
        y_LatGlobal        (jglobal) = aa_dRecvBuf(j,n)
        y_RankGlobal       (jglobal) = n
        y_JLocalIndexGlobal(jglobal) = j
      end do
    end do
    do n = 0, nprocs-1
      do j = a_jmaxs(n)+1, a_jmax(n)
        jglobal = jglobal + 1
        y_LatGlobal        (jglobal) = aa_dRecvBuf(j,n)
        y_RankGlobal       (jglobal) = n
        y_JLocalIndexGlobal(jglobal) = j
      end do
    end do
    deallocate( aa_dRecvBuf )
  end subroutine SLTTExtArrPrepGlobalArray
          | Subroutine : | |
| iVal : | integer , intent(in ) | 
| a_iVal(0:nprocs-1) : | integer , intent(out) | 
  subroutine SLTTExtShareiScalar( iVal, a_iVal )
    ! メッセージ出力
    ! Message output
    !
    use dc_message, only: MessageNotify
    !
    ! MPI
    !
    use mpi_wrapper, only : nprocs, myrank, MPIWrapperISend, MPIWrapperIRecv, MPIWrapperWait
    integer , intent(in ) :: iVal
    integer , intent(out) :: a_iVal(0:nprocs-1)
    !
    ! local variables
    !
    integer :: a_iSendBuf (1)
    integer :: aa_iRecvBuf(1,0:nprocs-1)
    integer :: a_iReqSend (0:nprocs-1)
    integer :: a_iReqRecv (0:nprocs-1)
    integer :: idep
    integer :: idest
    integer :: n
    a_iSendBuf = iVal
    do n = 0, nprocs-1
      if ( n == myrank ) then
        aa_iRecvBuf(:,n) = a_iSendBuf
      else
        idest = n
        call MPIWrapperISend( idest, 1, a_iSendBuf      , a_iReqSend(n) )
        idep = n
        call MPIWrapperIRecv( idep, 1, aa_iRecvBuf(:,n), a_iReqRecv(n) )
      end if
    end do
    do n = 0, nprocs-1
      if ( n == myrank ) cycle
      call MPIWrapperWait( a_iReqSend(n) )
      call MPIWrapperWait( a_iReqRecv(n) )
    end do
    !
    a_iVal = aa_iRecvBuf(1,:)
  end subroutine SLTTExtShareiScalar
          | Constant : | |||
| module_name = ‘sosi_dynamics_extarr‘ : | character(*), parameter 
 | 
| Variable : | |||
| sosi_dynamics_extarr_inited = .false. : | logical, save 
 | 
| Variable : | |
| ya_ExtFlagCrossPoleOrgDataN(:,:) : | logical , save, allocatable | 
| Variable : | |
| ya_ExtFlagCrossPoleOrgDataS(:,:) : | logical , save, allocatable |