!= MPI 関連ルーチン ! != MPI related routines ! ! Authors:: Yoshiyuki O. Takahashi ! Version:: $Id: mpi_wrapper.f90,v 1.1 2011-03-28 01:50:37 sugiyama Exp $ ! Tag Name:: $Name: arare4-20120911 $ ! Copyright:: Copyright (C) GFD Dennou Club, 2008. All rights reserved. ! License:: See COPYRIGHT[link:../../../COPYRIGHT] ! module mpi_wrapper ! != 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 ! ! モジュール引用 ; USE statements ! ! 種別型パラメタ ! Kind type parameter ! use dc_types, only: DP, & ! 倍精度実数型. Double precision. & STRING, & ! 文字列. Strings. & TOKEN ! キーワード. Keywords. ! MPI ! use mpi ! 宣言文 ; Declaration statements ! implicit none private ! 公開手続き ! Public procedure ! public :: MPIWrapperInit public :: MPIWrapperFinalize public :: MPIWrapperISend public :: MPIWrapperIRecv public :: MPIWrapperWait ! 公開変数 ! Public variables ! integer, save, public :: nprocs = 1 ! Number of MPI processes integer, save, public :: myrank = 0 ! My rank ! 非公開変数 ! Private variables ! interface MPIWrapperISend module procedure & MPIWrapperISend_dble_1d, & MPIWrapperISend_dble_2d, & MPIWrapperISend_dble_3d, & MPIWrapperISend_dble_4d end interface interface MPIWrapperIRecv module procedure & MPIWrapperIRecv_dble_1d, & MPIWrapperIRecv_dble_2d, & MPIWrapperIRecv_dble_3d, & MPIWrapperIRecv_dble_4d end interface interface MPIWrapperAbort module procedure & MPIWrapperStop end interface contains !-------------------------------------------------------------------------------------- subroutine MPIWrapperInit ! ! MPI の初期化 ! ! Initialization of MPI ! ! モジュール引用 ; USE statements ! ! 作業変数 ! Work variables ! integer :: ierr ! nprocs = 1 ! myrank = 0 call mpi_init( ierr ) call mpi_comm_size( mpi_comm_world, nprocs, ierr ) call mpi_comm_rank( mpi_comm_world, myrank, ierr ) end subroutine MPIWrapperInit !-------------------------------------------------------------------------------------- subroutine MPIWrapperFinalize ! ! MPI の終了処理 ! ! Finalization of MPI ! ! モジュール引用 ; USE statements ! ! 作業変数 ! Work variables ! integer :: ierr call mpi_finalize( ierr ) end subroutine MPIWrapperFinalize !-------------------------------------------------------------------------------------- subroutine MPIWrapperStop ! ! MPI の異常終了処理 ! ! Abort of MPI ! ! モジュール引用 ; USE statements ! ! 作業変数 ! Work variables ! integer :: errorcode = 9 integer :: ierr call mpi_abort( mpi_comm_world, errorcode, ierr ) call MPIWrapperFinalize stop end subroutine MPIWrapperstop !-------------------------------------------------------------------------------------- subroutine MPIWrapperWait( ireq ) ! ! MPI 通信終了まで待機 ! ! Wait finishing MPI transfer ! ! モジュール引用 ; USE statements ! integer, intent(inout) :: ireq ! request number ! 作業変数 ! Work variables ! integer :: ierr integer :: istatus( MPI_STATUS_SIZE ) call mpi_wait( ireq, istatus, ierr ) end subroutine MPIWrapperWait !-------------------------------------------------------------------------------------- subroutine MPIWrapperISend_dble_1d( & & idest, im, & ! (in) & buf, & ! (in) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_isend( buf, isize, & mpi_double_precision, idest, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperISend_dble_1d !-------------------------------------------------------------------------------------- subroutine MPIWrapperIRecv_dble_1d( & & idep, im, & ! (in) & buf, & ! (out) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_irecv( buf, isize, & mpi_double_precision, idep, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperIRecv_dble_1d !-------------------------------------------------------------------------------------- subroutine MPIWrapperISend_dble_2d( & & idest, im, jm, & ! (in) & buf, & ! (in) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_isend( buf, isize, & mpi_double_precision, idest, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperISend_dble_2d !-------------------------------------------------------------------------------------- subroutine MPIWrapperIRecv_dble_2d( & & idep, im, jm, & ! (in) & buf, & ! (out) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_irecv( buf, isize, & mpi_double_precision, idep, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperIRecv_dble_2d !-------------------------------------------------------------------------------------- subroutine MPIWrapperISend_dble_3d( & & idest, im, jm, km, & ! (in) & buf, & ! (in) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_isend( buf, isize, & mpi_double_precision, idest, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperISend_dble_3d !-------------------------------------------------------------------------------------- subroutine MPIWrapperIRecv_dble_3d( & & idep, im, jm, km, & ! (in) & buf, & ! (out) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_irecv( buf, isize, & mpi_double_precision, idep, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperIRecv_dble_3d !-------------------------------------------------------------------------------------- subroutine MPIWrapperISend_dble_4d( & & idest, im, jm, km, lm, & ! (in) & buf, & ! (in) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_isend( buf, isize, & mpi_double_precision, idest, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperISend_dble_4d !-------------------------------------------------------------------------------------- subroutine MPIWrapperIRecv_dble_4d( & & idep, im, jm, km, lm, & ! (in) & buf, & ! (out) & ireq & ! (out) & ) ! ! 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 ! 作業変数 ! Work variables ! integer :: ierr integer :: isize isize = size( buf ) call mpi_irecv( buf, isize, & mpi_double_precision, idep, 1, mpi_comm_world, & ireq, ierr ) end subroutine MPIWrapperIRecv_dble_4d !-------------------------------------------------------------------------------------- end module mpi_wrapper