Class mpi_wrapper
In: ../src/setup/mpi_wrapper.F90

MPI 関連ルーチン

MPI related routines

Note that Japanese and English are described in parallel.

MPI 関係の変数の管理と MPI 関係ラッパールーチンのモジュール.

This is a module containing MPI-related variables and wrapper routines.

Methods

Included Modules

dc_types mpi

Public Instance methods

FLAG_LIB_MPI
Variable :
FLAG_LIB_MPI :logical, save, public
: FLAG
Subroutine :
xsub :integer, intent(in)
: CPU nums of X and Y direction
ysub :integer, intent(in)
: CPU nums of X and Y direction
flag_periodic :logical, intent(in)
mpi_comm_cart :integer, intent(inout)

MPI cart の初期化

Initialization of MPI cart

[Source]

  subroutine MPIWrapperCartCreate( xsub, ysub, flag_periodic, mpi_comm_cart )
    !
    ! MPI cart の初期化
    !
    ! Initialization of MPI cart
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    !
    integer, intent(in)    :: xsub, ysub  ! CPU nums of X and Y direction
    logical, intent(in)    :: flag_periodic
    integer, intent(inout) :: mpi_comm_cart

    ! 作業変数
    !
    integer, parameter :: disp = 1
    integer, parameter :: ndim = 2
    logical, parameter :: reorder = .false.
    logical            :: periodic(ndim)
    integer            :: idivid(ndim)
    integer            :: direction
    integer            :: ierr

#ifdef LIB_MPI
      
    idivid(1) = ysub
    idivid(2) = xsub
    periodic(1) = flag_periodic
    periodic(2) = flag_periodic
    
    call mpi_cart_create( mpi_comm_world, ndim, idivid, periodic, reorder, mpi_comm_cart, ierr )    

!    if (.not. Flag_MPI_Cart) then 
!      Flag_MPI_Cart = .true.
!      mpi_comm_cart_save = mpi_comm_cart
!    end if

#endif
  end subroutine MPIWrapperCartCreate
Subroutine :
mpi_comm_cart :integer, intent(in)
direction :integer, intent(in)
disp :integer, intent(in)
rank1 :integer, intent(out)
rank2 :integer, intent(out)

[Source]

  subroutine MPIWrapperCartShift(mpi_comm_cart, direction, disp, rank1, rank2)

    integer, intent(in)   :: mpi_comm_cart
    integer, intent(in)   :: direction 
    integer, intent(in)   :: disp
    integer, intent(out)  :: rank1
    integer, intent(out)  :: rank2
    integer               :: ierr

    ! デフォルト値. 
    !
    rank1 = -1
    rank2 = -1

#ifdef LIB_MPI

    call mpi_cart_shift( mpi_comm_cart, direction, disp, rank1, rank2, ierr )

#endif
    
  end subroutine MPIWrapperCartShift
Subroutine :
mpi_comm_cart :integer, intent(in)

MPI の終了処理

Finalization of MPI

[Source]

  subroutine MPIWrapperCommFree(mpi_comm_cart)
    !
    ! MPI の終了処理
    !
    ! Finalization of MPI
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数 
    !
    integer, intent(in) :: mpi_comm_cart

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr

    call mpi_comm_free( mpi_comm_cart, ierr )

#endif

  end subroutine MPIWrapperCommFree
Subroutine :

MPI の終了処理

Finalization of MPI

[Source]

  subroutine MPIWrapperFinalize
    !
    ! MPI の終了処理
    !
    ! Finalization of MPI
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr

!    if (FLAG_MPI_Cart) then 
!      call mpi_comm_free( mpi_comm_cart_save, ierr )
!    end if
    call mpi_finalize( ierr )

#endif

  end subroutine MPIWrapperFinalize
Subroutine :
idep :integer , intent(in )
: Process number of departure
im :integer , intent(in )
: Size of 1st dimension of received data
buf( im ) :real(DP), intent(out)
: Array to be received
ireq :integer , intent(out)
: Request number

1D 倍精度配列の非ブロッキング通信(受信)

Non-blocking transfer (receive) of real(8) 1D array

[Source]

  subroutine MPIWrapperIRecv_dble_1d( idep, im, buf, ireq )
    !
    ! 1D 倍精度配列の非ブロッキング通信(受信)
    !
    ! Non-blocking transfer (receive) of real(8) 1D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idep
                              ! Process number of departure
    integer , intent(in ) :: im
                              ! Size of 1st dimension of received data
    real(DP), intent(out) :: buf( im )
                              ! Array to be received
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

    isize = size( buf )

    call mpi_irecv( buf, isize, mpi_double_precision, idep, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperIRecv_dble_1d
Subroutine :
idep :integer , intent(in )
: Process number of destination
im :integer , intent(in )
: Size of 1st dimension of received data
jm :integer , intent(in )
: Size of 2nd dimension of received data
buf( im, jm ) :real(DP), intent(out)
: Array to be received
ireq :integer , intent(out)
: Request number

2D 倍精度配列の非ブロッキング通信(受信)

Non-blocking transfer (receive) of real(8) 2D array

[Source]

  subroutine MPIWrapperIRecv_dble_2d( idep, im, jm, buf, ireq )
    !
    ! 2D 倍精度配列の非ブロッキング通信(受信)
    !
    ! Non-blocking transfer (receive) of real(8) 2D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idep
                              ! Process number of destination
    integer , intent(in ) :: im
                              ! Size of 1st dimension of received data
    integer , intent(in ) :: jm
                              ! Size of 2nd dimension of received data
    real(DP), intent(out) :: buf( im, jm )
                              ! Array to be received
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

    isize = size( buf )

    call mpi_irecv( buf, isize, mpi_double_precision, idep, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperIRecv_dble_2d
Subroutine :
idep :integer , intent(in )
: Process number of departure
im :integer , intent(in )
: Size of 1st dimension of received data
jm :integer , intent(in )
: Size of 2nd dimension of received data
km :integer , intent(in )
: Size of 3rd dimension of received data
buf( im, jm, km ) :real(DP), intent(inout)
: Array to be received
ireq :integer , intent(out)
: Request number

3D 倍精度配列の非ブロッキング通信(受信)

Non-blocking transfer (receive) of real(8) 3D array

[Source]

  subroutine MPIWrapperIRecv_dble_3d( idep, im, jm, km, buf, ireq )
    !
    ! 3D 倍精度配列の非ブロッキング通信(受信)
    !
    ! Non-blocking transfer (receive) of real(8) 3D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idep
                              ! Process number of departure
    integer , intent(in ) :: im
                              ! Size of 1st dimension of received data
    integer , intent(in ) :: jm
                              ! Size of 2nd dimension of received data
    integer , intent(in ) :: km
                              ! Size of 3rd dimension of received data
    real(DP), intent(inout) :: buf( im, jm, km )
                              ! Array to be received
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

    isize = size( buf )

!    write(*,*) "*** iRECV *** ", myrank, " <= ", idep

    call mpi_irecv( buf, isize, mpi_double_precision, idep, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperIRecv_dble_3d
Subroutine :
idep :integer , intent(in )
: Process number of departure
im :integer , intent(in )
: Size of 1st dimension of received data
jm :integer , intent(in )
: Size of 2nd dimension of received data
km :integer , intent(in )
: Size of 3rd dimension of received data
lm :integer , intent(in )
: Size of 4th dimension of received data
buf( im, jm, km, lm ) :real(DP), intent(out)
: Array to be received
ireq :integer , intent(out)
: Request number

4D 倍精度配列の非ブロッキング通信(受信)

Non-blocking transfer (receive) of real(8) 4D array

[Source]

  subroutine MPIWrapperIRecv_dble_4d( idep, im, jm, km, lm, buf, ireq )
    !
    ! 4D 倍精度配列の非ブロッキング通信(受信)
    !
    ! Non-blocking transfer (receive) of real(8) 4D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idep
                              ! Process number of departure
    integer , intent(in ) :: im
                              ! Size of 1st dimension of received data
    integer , intent(in ) :: jm
                              ! Size of 2nd dimension of received data
    integer , intent(in ) :: km
                              ! Size of 3rd dimension of received data
    integer , intent(in ) :: lm
                              ! Size of 4th dimension of received data
    real(DP), intent(out) :: buf( im, jm, km, lm )
                              ! Array to be received
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

    isize = size( buf )

    call mpi_irecv( buf, isize, mpi_double_precision, idep, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperIRecv_dble_4d
Subroutine :
idest :integer , intent(in )
: Process number of destination
im :integer , intent(in )
: Size of 1st dimension of sent data
buf( im ) :real(DP), intent(in )
: Array to be sent
ireq :integer , intent(out)
: Request number

1D 倍精度配列の非ブロッキング通信(送信)

Non-blocking transfer (send) of real(8) 1D array

[Source]

  subroutine MPIWrapperISend_dble_1d( idest, im, buf, ireq )
    !
    ! 1D 倍精度配列の非ブロッキング通信(送信)
    !
    ! Non-blocking transfer (send) of real(8) 1D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idest
                              ! Process number of destination
    integer , intent(in ) :: im
                              ! Size of 1st dimension of sent data
    real(DP), intent(in ) :: buf( im )
                              ! Array to be sent
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI
    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

    isize = size( buf )

    call mpi_isend( buf, isize, mpi_double_precision, idest, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperISend_dble_1d
Subroutine :
idest :integer , intent(in )
: Process number of destination
im :integer , intent(in )
: Size of 1st dimension of sent data
jm :integer , intent(in )
: Size of 2nd dimension of sent data
buf( im, jm ) :real(DP), intent(in )
: Array to be sent
ireq :integer , intent(out)
: Request number

2D 倍精度配列の非ブロッキング通信(送信)

Non-blocking transfer (send) of real(8) 2D array

[Source]

  subroutine MPIWrapperISend_dble_2d( idest, im, jm, buf, ireq )
    !
    ! 2D 倍精度配列の非ブロッキング通信(送信)
    !
    ! Non-blocking transfer (send) of real(8) 2D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idest
                              ! Process number of destination
    integer , intent(in ) :: im
                              ! Size of 1st dimension of sent data
    integer , intent(in ) :: jm
                              ! Size of 2nd dimension of sent data
    real(DP), intent(in ) :: buf( im, jm )
                              ! Array to be sent
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

    isize = size( buf )

    call mpi_isend( buf, isize, mpi_double_precision, idest, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperISend_dble_2d
Subroutine :
idest :integer , intent(in )
: Process number of destination
im :integer , intent(in )
: Size of 1st dimension of sent data
jm :integer , intent(in )
: Size of 2nd dimension of sent data
km :integer , intent(in )
: Size of 3rd dimension of sent data
buf( im, jm, km ) :real(DP), intent(in )
: Array to be sent
ireq :integer , intent(out)
: Request number

3D 倍精度配列の非ブロッキング通信(送信)

Non-blocking transfer (send) of real(8) 3D array

[Source]

  subroutine MPIWrapperISend_dble_3d( idest, im, jm, km, buf, ireq )
    !
    ! 3D 倍精度配列の非ブロッキング通信(送信)
    !
    ! Non-blocking transfer (send) of real(8) 3D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idest
                              ! Process number of destination
    integer , intent(in ) :: im
                              ! Size of 1st dimension of sent data
    integer , intent(in ) :: jm
                              ! Size of 2nd dimension of sent data
    integer , intent(in ) :: km
                              ! Size of 3rd dimension of sent data
    real(DP), intent(in ) :: buf( im, jm, km )
                              ! Array to be sent
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

!    write(*,*) "*** iSEND *** ", myrank, " => ", idest

    isize = size( buf )

    call mpi_isend( buf, isize, mpi_double_precision, idest, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperISend_dble_3d
Subroutine :
idest :integer , intent(in )
: Process number of destination
im :integer , intent(in )
: Size of 1st dimension of sent data
jm :integer , intent(in )
: Size of 2nd dimension of sent data
km :integer , intent(in )
: Size of 3rd dimension of sent data
lm :integer , intent(in )
: Size of 4th dimension of sent data
buf( im, jm, km, lm ) :real(DP), intent(in )
: Array to be sent
ireq :integer , intent(out)
: Request number

4D 倍精度配列の非ブロッキング通信(送信)

Non-blocking transfer (send) of real(8) 4D array

[Source]

  subroutine MPIWrapperISend_dble_4d( idest, im, jm, km, lm, buf, ireq )
    !
    ! 4D 倍精度配列の非ブロッキング通信(送信)
    !
    ! Non-blocking transfer (send) of real(8) 4D array
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 入出力変数
    ! input/output variables
    !
    integer , intent(in ) :: idest
                              ! Process number of destination
    integer , intent(in ) :: im
                              ! Size of 1st dimension of sent data
    integer , intent(in ) :: jm
                              ! Size of 2nd dimension of sent data
    integer , intent(in ) :: km
                              ! Size of 3rd dimension of sent data
    integer , intent(in ) :: lm
                              ! Size of 4th dimension of sent data
    real(DP), intent(in ) :: buf( im, jm, km, lm )
                              ! Array to be sent
    integer , intent(out) :: ireq
                              ! Request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: isize

    isize = size( buf )

    call mpi_isend( buf, isize, mpi_double_precision, idest, 1, mpi_comm_world, ireq, ierr )

#endif 

  end subroutine MPIWrapperISend_dble_4d
Subroutine :

MPI の初期化

Initialization of MPI

[Source]

  subroutine MPIWrapperInit
    !
    ! MPI の初期化
    !
    ! Initialization of MPI
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr

#endif

    ! set default values
    ! 
    FLAG_LIB_MPI = .false.
    nprocs = 1
    myrank = 0

#ifdef LIB_MPI

    ! call MPI initial subroutine
    ! 
    FLAG_LIB_MPI = .true.
    call mpi_init( ierr )
    call mpi_comm_size( mpi_comm_world, nprocs, ierr )
    call mpi_comm_rank( mpi_comm_world, myrank, ierr )

#endif

  end subroutine MPIWrapperInit
Subroutine :
ireq :integer, intent(inout)
: request number

MPI 通信終了まで待機

Wait finishing MPI transfer

[Source]

  subroutine MPIWrapperWait( ireq )
    !
    ! MPI 通信終了まで待機
    !
    ! Wait finishing MPI transfer
    !

    ! 暗黙の型宣言禁止
    !
    implicit none

    ! 作業変数
    ! Work variables
    !
    integer, intent(inout) :: ireq ! request number

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: ierr
    integer :: istatus( MPI_STATUS_SIZE )

    call mpi_wait( ireq, istatus, ierr )

#endif 

  end subroutine MPIWrapperWait
myrank
Variable :
myrank :integer, save, public
: My rank
nprocs
Variable :
nprocs :integer, save, public
: Number of MPI processes