Public Instance methods
MPIWrapperChkTrue( lmax, a_LocalLogical, a_GlobalLogical )
Subroutine : |
|
lmax : | integer, intent(in )
|
a_LocalLogical(1:lmax) : | logical, intent(in )
|
a_GlobalLogical(1:lmax) : | logical, intent(out)
|
全球の最大値を探す
Find the maximum of the globe
Alias for MPIWrapperChkTrue_1d
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, [itag] )
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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 1D array
Alias for MPIWrapperIRecv_int_1d
MPIWrapperIRecv( idep, im, buf, ireq, [itag] )
Subroutine : |
|
idep : | integer , intent(in )
: | Process number of departure
|
|
im : | integer , intent(in )
: | Size of 1st dimension of received data
|
|
buf( im ) : | logical , intent(out)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 1D array
Alias for MPIWrapperIRecv_logical_1d
MPIWrapperIRecv( idep, im, buf, ireq, [itag] )
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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 1D array
Alias for MPIWrapperIRecv_dble_1d
MPIWrapperIRecv( idep, im, jm, buf, ireq, [itag] )
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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
2D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 2D array
Alias for MPIWrapperIRecv_dble_2d
MPIWrapperIRecv( idep, im, jm, km, buf, ireq, [itag] )
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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
3D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 3D array
Alias for MPIWrapperIRecv_dble_3d
MPIWrapperIRecv( idep, im, jm, km, lm, buf, ireq, [itag] )
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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
4D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 4D array
Alias for MPIWrapperIRecv_dble_4d
MPIWrapperISend( idest, im, buf, ireq, [itag] )
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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 1D array
Alias for MPIWrapperISend_int_1d
MPIWrapperISend( idest, im, buf, ireq, [itag] )
Subroutine : |
|
idest : | integer , intent(in )
: | Process number of destination
|
|
im : | integer , intent(in )
: | Size of 1st dimension of sent data
|
|
buf( im ) : | logical , intent(in )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 1D array
Alias for MPIWrapperISend_logical_1d
MPIWrapperISend( idest, im, buf, ireq, [itag] )
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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 1D array
Alias for MPIWrapperISend_dble_1d
MPIWrapperISend( idest, im, jm, buf, ireq, [itag] )
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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
2D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 2D array
Alias for MPIWrapperISend_dble_2d
MPIWrapperISend( idest, im, jm, km, buf, ireq, [itag] )
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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
3D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 3D array
Alias for MPIWrapperISend_dble_3d
MPIWrapperISend( idest, im, jm, km, lm, buf, ireq, [itag] )
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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
4D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 4D array
Alias for MPIWrapperISend_dble_4d
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)
|
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
|
nprocs
Variable : |
|
nprocs : | integer, save, public
: | Number of MPI processes
|
|
Private Instance methods
Subroutine : |
|
lmax : | integer, intent(in )
|
a_LocalLogical(1:lmax) : | logical, intent(in )
|
a_GlobalLogical(1:lmax) : | logical, intent(out)
|
全球の最大値を探す
Find the maximum of the globe
[Source]
subroutine MPIWrapperChkTrue_1d( lmax, a_LocalLogical, a_GlobalLogical )
!
! 全球の最大値を探す
!
! Find the maximum of the globe
!
! モジュール引用 ; USE statements
!
integer, intent(in ) :: lmax
logical, intent(in ) :: a_LocalLogical (1:lmax)
logical, intent(out) :: a_GlobalLogical(1:lmax)
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: idep ! Process number of departure
integer :: idest ! Process number of destination
logical :: a_Buf(1:lmax) ! Array to be received
integer :: ireq ! Request number
logical, allocatable :: aa_LocalLogical(:,:)
integer :: l
integer :: n
if ( myrank == 0 ) then
allocate( aa_LocalLogical( 1:lmax, 0:nprocs-1 ) )
aa_LocalLogical(:,0) = a_LocalLogical
do n = 1, nprocs-1
idep = n
call MPIWrapperIRecv( idep, lmax, aa_LocalLogical(:,n), ireq )
call MPIWrapperWait( ireq )
end do
do l = 1, lmax
n = 0
a_GlobalLogical(l) = aa_LocalLogical(l,n)
do n = 1, nprocs-1
if ( aa_LocalLogical(l,n) ) then
a_GlobalLogical(l) = aa_LocalLogical(l,n)
end if
end do
end do
a_Buf = a_GlobalLogical
do n = 1, nprocs-1
idest = n
call MPIWrapperISend( idest, lmax, a_Buf, ireq )
call MPIWrapperWait( ireq )
end do
deallocate( aa_LocalLogical )
else
idest = 0
a_Buf = a_LocalLogical
call MPIWrapperISend( idest, lmax, a_Buf, ireq )
call MPIWrapperWait( ireq )
idep = 0
call MPIWrapperIRecv( idep, lmax, a_Buf, ireq )
call MPIWrapperWait( ireq )
a_GlobalLogical = a_Buf
end if
#else
a_GlobalLogical = a_LocalLogical
#endif
end subroutine MPIWrapperChkTrue_1d
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
real(DP), allocatable :: aa_LocalMax(:,:)
integer :: l
integer :: n
if ( myrank == 0 ) then
allocate( aa_LocalMax( 1:lmax, 0:nprocs-1 ) )
aa_LocalMax(:,0) = a_LocalMax
do n = 1, nprocs-1
idep = n
call MPIWrapperIRecv( idep, lmax, aa_LocalMax(:,n), ireq )
call MPIWrapperWait( ireq )
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, ireq )
call MPIWrapperWait( ireq )
end do
deallocate( aa_LocalMax )
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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 1D array
[Source]
subroutine MPIWrapperIRecv_dble_1d( idep, im, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_irecv( buf, isize, mpi_double_precision, idep, tag, 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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
2D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 2D array
[Source]
subroutine MPIWrapperIRecv_dble_2d( idep, im, jm, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_irecv( buf, isize, mpi_double_precision, idep, tag, 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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
3D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 3D array
[Source]
subroutine MPIWrapperIRecv_dble_3d( idep, im, jm, km, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_irecv( buf, isize, mpi_double_precision, idep, tag, 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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
4D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 4D array
[Source]
subroutine MPIWrapperIRecv_dble_4d( idep, im, jm, km, lm, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_irecv( buf, isize, mpi_double_precision, idep, tag, 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)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 1D array
[Source]
subroutine MPIWrapperIRecv_int_1d( idep, im, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_irecv( buf, isize, mpi_integer, idep, tag, mpi_comm_world, ireq, ierr )
#endif
end subroutine MPIWrapperIRecv_int_1d
Subroutine : |
|
idep : | integer , intent(in )
: | Process number of departure
|
|
im : | integer , intent(in )
: | Size of 1st dimension of received data
|
|
buf( im ) : | logical , intent(out)
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(受信)
Non-blocking transfer (receive) of real(8) 1D array
[Source]
subroutine MPIWrapperIRecv_logical_1d( idep, im, buf, ireq, itag )
!
! 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
logical , intent(out) :: buf( im )
! Array to be received
integer , intent(out) :: ireq
! Request number
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_irecv( buf, isize, mpi_logical, idep, tag, mpi_comm_world, ireq, ierr )
#endif
end subroutine MPIWrapperIRecv_logical_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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 1D array
[Source]
subroutine MPIWrapperISend_dble_1d( idest, im, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_isend( buf, isize, mpi_double_precision, idest, tag, 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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
2D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 2D array
[Source]
subroutine MPIWrapperISend_dble_2d( idest, im, jm, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_isend( buf, isize, mpi_double_precision, idest, tag, 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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
3D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 3D array
[Source]
subroutine MPIWrapperISend_dble_3d( idest, im, jm, km, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_isend( buf, isize, mpi_double_precision, idest, tag, 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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
4D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 4D array
[Source]
subroutine MPIWrapperISend_dble_4d( idest, im, jm, km, lm, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_isend( buf, isize, mpi_double_precision, idest, tag, 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 )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 1D array
[Source]
subroutine MPIWrapperISend_int_1d( idest, im, buf, ireq, itag )
!
! 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
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_isend( buf, isize, mpi_integer, idest, tag, mpi_comm_world, ireq, ierr )
#endif
end subroutine MPIWrapperISend_int_1d
Subroutine : |
|
idest : | integer , intent(in )
: | Process number of destination
|
|
im : | integer , intent(in )
: | Size of 1st dimension of sent data
|
|
buf( im ) : | logical , intent(in )
|
ireq : | integer , intent(out)
|
itag : | integer , intent(in ), optional
: | Size of 1st dimension of sent data
|
|
1D 倍精度配列の非ブロッキング通信(送信)
Non-blocking transfer (send) of real(8) 1D array
[Source]
subroutine MPIWrapperISend_logical_1d( idest, im, buf, ireq, itag )
!
! 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
logical , intent(in ) :: buf( im )
! Array to be sent
integer , intent(out) :: ireq
! Request number
integer , intent(in ), optional :: itag
! Size of 1st dimension of sent data
#ifdef LIB_MPI
! 作業変数
! Work variables
!
integer :: ierr
integer :: isize
integer :: tag
isize = size( buf )
if ( present( itag ) ) then
tag = itag
else
tag = 1
end if
call mpi_isend( buf, isize, mpi_logical, idest, tag, mpi_comm_world, ireq, ierr )
#endif
end subroutine MPIWrapperISend_logical_1d
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