Class mpi_wrapper
In: 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.

Procedures List

!$ ! RadiationFluxDennouAGCM :放射フラックスの計算
!$ ! RadiationFinalize :終了処理 (モジュール内部の変数の割り付け解除)
!$ ! ———— :————
!$ ! RadiationFluxDennouAGCM :Calculate radiation flux
!$ ! RadiationFinalize :Termination (deallocate variables in this module)

NAMELIST

!$ ! NAMELIST#radiation_DennouAGCM_nml

Methods

Included Modules

dc_types mpi

Public Instance methods

Subroutine :

MPI の終了処理

Finalization of MPI

[Source]

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

    ! モジュール引用 ; USE statements
    !


#ifdef LIB_MPI

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


    call mpi_finalize( ierr )

#endif

  end subroutine MPIWrapperFinalize
MPIWrapperFindMaxVal( lmax, a_LocalMax, a_GlobalMax )
Subroutine :
lmax :integer , intent(in )
a_LocalMax(1:lmax) :real(DP), intent(in )
a_GlobalMax(1:lmax) :real(DP), intent(out)

全球の最大値を探す

Find the maximum of the globe

Alias for MPIWrapperFindMaxVal_dble_1d

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

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

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

Alias for MPIWrapperIRecv_int_1d

MPIWrapperIRecv( idep, im, buf, ireq )
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

Alias for MPIWrapperIRecv_dble_1d

MPIWrapperIRecv( idep, im, jm, buf, ireq )
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

Alias for MPIWrapperIRecv_dble_2d

MPIWrapperIRecv( idep, im, jm, km, buf, ireq )
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(out)
: Array to be received
ireq :integer , intent(out)
: Request number

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

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

Alias for MPIWrapperIRecv_dble_3d

MPIWrapperIRecv( idep, im, jm, km, lm, buf, ireq )
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

Alias for MPIWrapperIRecv_dble_4d

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

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

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

Alias for MPIWrapperISend_int_1d

MPIWrapperISend( idest, im, buf, ireq )
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

Alias for MPIWrapperISend_dble_1d

MPIWrapperISend( idest, im, jm, buf, ireq )
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

Alias for MPIWrapperISend_dble_2d

MPIWrapperISend( idest, im, jm, km, buf, ireq )
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

Alias for MPIWrapperISend_dble_3d

MPIWrapperISend( idest, im, jm, km, lm, buf, ireq )
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

Alias for MPIWrapperISend_dble_4d

Subroutine :

MPI の初期化

Initialization of MPI

[Source]

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

    ! モジュール引用 ; USE statements
    !


#ifdef LIB_MPI

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

#endif

    nprocs = 1
    myrank = 0

#ifdef LIB_MPI

    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
    !

    ! モジュール引用 ; USE statements
    !


    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

Private Instance methods

MPIWrapperAbort( )
Subroutine :

MPI の異常終了処理

Abort of MPI

Alias for MPIWrapperStop

Subroutine :
lmax :integer , intent(in )
a_LocalMax(1:lmax) :real(DP), intent(in )
a_GlobalMax(1:lmax) :real(DP), intent(out)

全球の最大値を探す

Find the maximum of the globe

[Source]

  subroutine MPIWrapperFindMaxVal_dble_1d( lmax, a_LocalMax, a_GlobalMax )
    !
    ! 全球の最大値を探す
    !
    ! Find the maximum of the globe
    !

    ! モジュール引用 ; USE statements
    !

    integer , intent(in ) :: lmax
    real(DP), intent(in ) :: a_LocalMax (1:lmax)
    real(DP), intent(out) :: a_GlobalMax(1:lmax)


#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer  :: idep          ! Process number of departure
    integer  :: idest         ! Process number of destination
    real(DP) :: a_Buf(1:lmax) ! Array to be received
    integer  :: ireq          ! Request number
    integer , allocatable :: a_iReqSend(:)
    integer , allocatable :: a_iReqRecv(:)
    real(DP), allocatable :: aa_LocalMax(:,:)
    integer               :: l
    integer               :: n


    if ( myrank == 0 ) then

      allocate( aa_LocalMax( 1:lmax, 0:nprocs-1 ) )
      allocate( a_iReqRecv ( 0:nprocs-1 ) )
      allocate( a_iReqSend ( 0:nprocs-1 ) )

      aa_LocalMax(:,0) = a_LocalMax

      do n = 1, nprocs-1
        idep = n
        call MPIWrapperIRecv( idep, lmax, aa_LocalMax(:,n), a_iReqRecv(n) )
      end do
      do n = 1, nprocs-1
        call MPIWrapperWait( a_iReqRecv(n) )
      end do

      do l = 1, lmax
        n = 0
        a_GlobalMax(l) = aa_LocalMax(l,n)
        do n = 1, nprocs-1
          if ( a_GlobalMax(l) < aa_LocalMax(l,n) ) then
            a_GlobalMax(l) = aa_LocalMax(l,n)
          end if
        end do
      end do

      a_Buf = a_GlobalMax
      do n = 1, nprocs-1
        idest = n
        call MPIWrapperISend( idest, lmax, a_Buf, a_iReqSend(n) )
      end do
      do n = 1, nprocs-1
        call MPIWrapperWait( a_iReqSend(n) )
      end do

      deallocate( aa_LocalMax )
      deallocate( a_iReqRecv )
      deallocate( a_iReqSend )

    else

      idest = 0
      a_Buf = a_LocalMax
      call MPIWrapperISend( idest, lmax, a_Buf, ireq )
      call MPIWrapperWait( ireq )

      idep = 0
      call MPIWrapperIRecv( idep, lmax, a_Buf, ireq )
      call MPIWrapperWait( ireq )

      a_GlobalMax = a_Buf

    end if

#else

    a_GlobalMax = a_LocalMax

#endif


  end subroutine MPIWrapperFindMaxVal_dble_1d
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
    !

    ! モジュール引用 ; USE statements
    !


    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
    !

    ! モジュール引用 ; USE statements
    !


    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(out)
: 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
    !

    ! モジュール引用 ; USE statements
    !

    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(out) :: 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 )

    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
    !

    ! モジュール引用 ; USE statements
    !


    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 :
idep :integer , intent(in )
: Process number of departure
im :integer , intent(in )
: Size of 1st dimension of received data
buf( im ) :integer , 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_int_1d( idep, im, buf, ireq )
    !
    ! 1D 倍精度配列の非ブロッキング通信(受信)
    !
    ! Non-blocking transfer (receive) of real(8) 1D array
    !

    ! モジュール引用 ; USE statements
    !


    integer , intent(in ) :: idep
                              ! Process number of departure
    integer , intent(in ) :: im
                              ! Size of 1st dimension of received data
    integer , 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_int_1d
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
    !

    ! モジュール引用 ; USE statements
    !


    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
    !

    ! モジュール引用 ; USE statements
    !


    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
    !

    ! モジュール引用 ; USE statements
    !


    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


    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
    !

    ! モジュール引用 ; USE statements
    !


    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 :
idest :integer , intent(in )
: Process number of destination
im :integer , intent(in )
: Size of 1st dimension of sent data
buf( im ) :integer , 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_int_1d( idest, im, buf, ireq )
    !
    ! 1D 倍精度配列の非ブロッキング通信(送信)
    !
    ! Non-blocking transfer (send) of real(8) 1D array
    !

    ! モジュール引用 ; USE statements
    !


    integer , intent(in ) :: idest
                              ! Process number of destination
    integer , intent(in ) :: im
                              ! Size of 1st dimension of sent data
    integer , 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_int_1d
Subroutine :

MPI の異常終了処理

Abort of MPI

[Source]

  subroutine MPIWrapperStop
    !
    ! MPI の異常終了処理
    !
    ! Abort of MPI
    !

    ! モジュール引用 ; USE statements
    !

#ifdef LIB_MPI

    ! 作業変数
    ! Work variables
    !
    integer :: errorcode = 9
    integer :: ierr


    call mpi_abort( mpi_comm_world, errorcode, ierr )
    call MPIWrapperFinalize
    stop

#endif

  end subroutine MPIWrapperstop