2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-12-19 17:33:29 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This file is part of the AMUN source code, a program to perform
|
|
|
|
!! Newtonian or relativistic magnetohydrodynamical simulations on uniform or
|
|
|
|
!! adaptive mesh.
|
2008-12-19 17:33:29 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! Copyright (C) 2008-2012 Grzegorz Kowal <grzegorz@amuncode.org>
|
2008-12-19 17:33:29 -06:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!! This program is free software: you can redistribute it and/or modify
|
|
|
|
!! it under the terms of the GNU General Public License as published by
|
|
|
|
!! the Free Software Foundation, either version 3 of the License, or
|
|
|
|
!! (at your option) any later version.
|
2008-12-19 17:33:29 -06:00
|
|
|
!!
|
2011-04-29 11:21:30 -03:00
|
|
|
!! This program is distributed in the hope that it will be useful,
|
2008-12-19 17:33:29 -06:00
|
|
|
!! but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
!! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
!! GNU General Public License for more details.
|
|
|
|
!!
|
|
|
|
!! You should have received a copy of the GNU General Public License
|
2012-07-22 12:30:20 -03:00
|
|
|
!! along with this program. If not, see <http://www.gnu.org/licenses/>.
|
2008-12-19 17:33:29 -06:00
|
|
|
!!
|
2010-10-13 03:32:10 -03:00
|
|
|
!!******************************************************************************
|
2008-12-19 17:33:29 -06:00
|
|
|
!!
|
2012-07-22 19:01:27 -03:00
|
|
|
!! module: MPITOOLS
|
|
|
|
!!
|
|
|
|
!! This module provides wrapper subroutines handling the parallel execution
|
|
|
|
!! with the Message Passing Interface protocol.
|
2012-07-22 12:30:20 -03:00
|
|
|
!!
|
2012-07-28 12:29:44 -03:00
|
|
|
!!
|
2012-07-22 12:30:20 -03:00
|
|
|
!!******************************************************************************
|
2008-12-19 17:33:29 -06:00
|
|
|
!
|
|
|
|
module mpitools
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external subroutines
|
|
|
|
!
|
|
|
|
use timers, only : set_timer, start_timer, stop_timer
|
|
|
|
|
|
|
|
! module variables are not implicit by default
|
|
|
|
!
|
2008-12-19 17:33:29 -06:00
|
|
|
implicit none
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! timer indices
|
|
|
|
!
|
|
|
|
integer , save :: imi, imc
|
|
|
|
|
2008-12-19 17:33:29 -06:00
|
|
|
! MPI global variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer(kind=4), save :: comm3d
|
|
|
|
integer(kind=4), save :: nproc, nprocs
|
|
|
|
integer(kind=4), save, dimension(3) :: pdims, pcoords, pparity
|
|
|
|
integer(kind=4), save, dimension(3,2) :: pneighs
|
|
|
|
logical , save, dimension(3) :: periodic
|
|
|
|
logical , save :: master = .true.
|
2008-12-22 14:57:31 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! by default everything is public
|
|
|
|
!
|
|
|
|
public
|
|
|
|
|
|
|
|
!- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
contains
|
|
|
|
!
|
|
|
|
!===============================================================================
|
2012-07-28 12:29:44 -03:00
|
|
|
!!
|
|
|
|
!!*** PUBLIC SUBROUTINES *****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
! subroutine INITIALIZE_MPITOOLS:
|
|
|
|
! ------------------------------
|
2012-07-22 19:01:27 -03:00
|
|
|
!
|
|
|
|
! Subroutine initializes the MPITOOLS modules.
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
subroutine initialize_mpitools()
|
2008-12-22 14:57:31 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
#ifdef MPI
|
2012-07-22 19:01:27 -03:00
|
|
|
use mpi, only : mpi_comm_world, mpi_success
|
2008-12-22 14:57:31 -06:00
|
|
|
#endif /* MPI */
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
implicit none
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2008-12-22 14:57:31 -06:00
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
#ifdef MPI
|
|
|
|
integer :: iret
|
2008-12-22 14:57:31 -06:00
|
|
|
#endif /* MPI */
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
#ifdef MPI
|
|
|
|
! set timer descriptions
|
|
|
|
!
|
|
|
|
call set_timer('MPI initialization', imi)
|
|
|
|
call set_timer('MPI communication' , imc)
|
|
|
|
|
|
|
|
! start time accounting for the MPI initialization
|
|
|
|
!
|
|
|
|
call start_timer(imi)
|
|
|
|
#endif /* MPI */
|
|
|
|
|
|
|
|
! initialize parralel execution parameters
|
|
|
|
!
|
|
|
|
nproc = 0
|
|
|
|
nprocs = 1
|
|
|
|
pdims(:) = 1
|
|
|
|
pcoords(:) = 0
|
|
|
|
pparity(:) = 0
|
|
|
|
pneighs(:,:) = -1
|
|
|
|
periodic(:) = .false.
|
2008-12-22 14:57:31 -06:00
|
|
|
|
|
|
|
#ifdef MPI
|
2012-07-22 19:01:27 -03:00
|
|
|
! initialize the MPI interface
|
|
|
|
!
|
|
|
|
call mpi_init(iret)
|
|
|
|
|
|
|
|
! check if the MPI interface was initialized successfully
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success) then
|
|
|
|
write(*,*) 'The MPI interface could not be initializes! Exiting...'
|
|
|
|
write(*,*)
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
|
|
|
|
! obtain the total number of processes
|
|
|
|
!
|
|
|
|
call mpi_comm_size(mpi_comm_world, nprocs, iret)
|
|
|
|
|
|
|
|
! check if the total number of processes could be obtained
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success) then
|
|
|
|
write(*,*) 'The MPI process ID could not be obtained! Exiting...'
|
|
|
|
write(*,*)
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
|
|
|
|
! obtain the current process identificator
|
|
|
|
!
|
|
|
|
call mpi_comm_rank(mpi_comm_world, nproc , iret)
|
|
|
|
|
|
|
|
! check if the process ID was return successfully
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
if (iret .ne. mpi_success) then
|
|
|
|
write(*,*) 'The MPI process ID could not be obtained! Exiting...'
|
|
|
|
write(*,*)
|
|
|
|
stop
|
|
|
|
end if
|
2008-12-22 14:57:31 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! set the master flag
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
master = (nproc .eq. 0)
|
2008-12-22 14:57:31 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! store the MPI pool handles
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
comm3d = mpi_comm_world
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! stop time accounting for the MPI initialization
|
|
|
|
!
|
|
|
|
call stop_timer(imi)
|
2008-12-22 14:57:31 -06:00
|
|
|
#endif /* MPI */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
end subroutine initialize_mpitools
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
! subroutine FINALIZE_MPITOOLS:
|
|
|
|
! ----------------------------
|
2012-07-22 19:01:27 -03:00
|
|
|
!
|
|
|
|
! Subroutine finalizes the MPITOOLS modules.
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
subroutine finalize_mpitools()
|
2008-12-22 14:57:31 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
#ifdef MPI
|
2012-07-22 19:01:27 -03:00
|
|
|
use mpi, only : mpi_comm_world, mpi_success
|
|
|
|
#endif /* MPI */
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
2008-12-22 14:57:31 -06:00
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
#ifdef MPI
|
|
|
|
integer :: iret
|
2008-12-22 14:57:31 -06:00
|
|
|
#endif /* MPI */
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
#ifdef MPI
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI initialization
|
|
|
|
!
|
|
|
|
call start_timer(imi)
|
|
|
|
|
|
|
|
! initialize the MPI interface
|
|
|
|
!
|
|
|
|
call mpi_finalize(iret)
|
|
|
|
|
|
|
|
! check if the MPI interface was finalizes successfully
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
if (iret .ne. mpi_success) then
|
|
|
|
if (master) then
|
|
|
|
write(*,*) 'The MPI interface could not be finalized! Exiting...'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
stop
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI initialization
|
|
|
|
!
|
|
|
|
call stop_timer(imi)
|
2008-12-22 14:57:31 -06:00
|
|
|
#endif /* MPI */
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-28 12:29:44 -03:00
|
|
|
end subroutine finalize_mpitools
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine SETUP_MPI:
|
|
|
|
! --------------------
|
|
|
|
!
|
|
|
|
! Subroutine sets the MPI geometry.
|
2008-12-31 16:07:23 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
subroutine setup_mpi(div, per, set)
|
2008-12-31 16:07:23 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
2008-12-31 16:07:23 -06:00
|
|
|
#ifdef MPI
|
2012-07-22 19:01:27 -03:00
|
|
|
use mpi, only : mpi_comm_world, mpi_success
|
|
|
|
#endif /* MPI */
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! input arguments
|
|
|
|
!
|
|
|
|
integer, dimension(3), intent(in) :: div
|
|
|
|
logical, dimension(3), intent(in) :: per
|
2012-07-28 13:48:15 -03:00
|
|
|
logical , intent(in) :: set
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2008-12-31 16:07:23 -06:00
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer :: iret
|
2008-12-31 16:07:23 -06:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI initialization
|
|
|
|
!
|
|
|
|
call start_timer(imi)
|
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
! set the periodic flag
|
|
|
|
!
|
|
|
|
periodic(:) = per(:)
|
|
|
|
|
2012-08-02 15:09:44 -03:00
|
|
|
#ifdef MPI
|
2012-07-28 13:48:15 -03:00
|
|
|
! if set = .true. set the MPI domain
|
|
|
|
!
|
|
|
|
if (set) then
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! check if the total number of chunks in division corresponds to the number of
|
|
|
|
! processes, if not try to find the best division
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
if (nprocs .ne. product(div(:))) then
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
if (master) then
|
|
|
|
write(*,*) 'The number of MPI processes does not correspond to' &
|
|
|
|
// ' the number of domain chunks!'
|
|
|
|
write(*,*) 'Looking for the best division...'
|
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! try to find the best division
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
pdims(:) = 1
|
|
|
|
iret = 0
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
do while(product(pdims(:)) .lt. nprocs)
|
2012-07-22 19:01:27 -03:00
|
|
|
#ifdef R3D
|
2012-07-28 13:48:15 -03:00
|
|
|
iret = mod(iret, 3) + 1
|
2012-07-22 19:01:27 -03:00
|
|
|
#else /* R3D */
|
2012-07-28 13:48:15 -03:00
|
|
|
iret = mod(iret, 2) + 1
|
2012-07-22 19:01:27 -03:00
|
|
|
#endif /* R3D */
|
2012-07-28 13:48:15 -03:00
|
|
|
pdims(iret) = 2 * pdims(iret)
|
|
|
|
end do
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! check if the best division found
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
if (product(pdims(:)) .ne. nprocs) then
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
if (master) then
|
|
|
|
write(*,*) 'Improssible to find the best domain division! Exiting...'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
call finalize_mpitools()
|
|
|
|
stop
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
if (master) then
|
|
|
|
write(*,*) 'Found the best division:', pdims(:)
|
|
|
|
write(*,*)
|
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
else
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! substitute div(:) to pdims(:)
|
2008-12-31 16:07:23 -06:00
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
pdims(:) = div(:)
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! set up the Cartesian geometry
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
call mpi_cart_create(mpi_comm_world, 3, pdims(:), periodic(:) &
|
2012-07-22 19:01:27 -03:00
|
|
|
, .true., comm3d, iret)
|
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
if (iret .ne. mpi_success) then
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
if (master) then
|
|
|
|
write(*,*) 'The MPI could not create the Cartesian geometry! Exiting...'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
stop
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! assign process coordinate
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
call mpi_cart_coords(comm3d, nproc, 3, pcoords(:), iret)
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
if (iret .ne. mpi_success) then
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
if (master) then
|
|
|
|
write(*,*) 'The MPI could not assign process coordinates! Exiting...'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
stop
|
2012-07-22 19:01:27 -03:00
|
|
|
|
2012-07-28 13:48:15 -03:00
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! set the neighbors
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
if (pdims(1) .gt. 1) then
|
|
|
|
call mpi_cart_shift(comm3d, 0, 1, pneighs(1,1), pneighs(1,2), iret)
|
|
|
|
end if
|
|
|
|
if (pdims(2) .gt. 1) then
|
|
|
|
call mpi_cart_shift(comm3d, 1, 1, pneighs(2,1), pneighs(2,2), iret)
|
|
|
|
end if
|
|
|
|
if (pdims(3) .gt. 1) then
|
|
|
|
call mpi_cart_shift(comm3d, 2, 1, pneighs(3,1), pneighs(3,2), iret)
|
|
|
|
end if
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! set parity flag
|
|
|
|
!
|
2012-07-28 13:48:15 -03:00
|
|
|
pparity(1) = mod(pcoords(1), 2)
|
|
|
|
pparity(2) = mod(pcoords(2), 2)
|
|
|
|
pparity(3) = mod(pcoords(3), 2)
|
|
|
|
|
|
|
|
end if
|
2012-08-02 15:09:44 -03:00
|
|
|
#endif /* MPI */
|
2012-07-22 19:01:27 -03:00
|
|
|
|
|
|
|
! stop time accounting for the MPI initialization
|
|
|
|
!
|
|
|
|
call stop_timer(imi)
|
2008-12-31 16:07:23 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine setup_mpi
|
|
|
|
#ifdef MPI
|
2008-12-31 16:07:23 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine BCAST_INTEGER_VARIABLE:
|
|
|
|
! ---------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine broadcast an integer variable from the master process to all
|
|
|
|
! other processes.
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine bcast_integer_variable(ibuf, iret)
|
2008-12-22 14:57:31 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_integer, mpi_success
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2008-12-22 14:57:31 -06:00
|
|
|
implicit none
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine arguments
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer, intent(inout) :: ibuf
|
|
|
|
integer, intent(inout) :: iret
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-22 14:57:31 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_bcast(ibuf, 1, mpi_integer, 0, comm3d, iret)
|
|
|
|
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not broadcast an integer variable!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2008-12-22 14:57:31 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine bcast_integer_variable
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine BCAST_REAL_VARIABLE:
|
|
|
|
! ------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine broadcast a real variable from the master process to all
|
|
|
|
! other processes.
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine bcast_real_variable(rbuf, iret)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_success
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2008-12-31 12:02:36 -06:00
|
|
|
implicit none
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
real , intent(inout) :: rbuf
|
|
|
|
integer, intent(inout) :: iret
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
call mpi_bcast(rbuf, 1, mpi_real8, 0, comm3d, iret)
|
|
|
|
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not broadcast an integer variable!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine bcast_real_variable
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine BCAST_STRING_VARIABLE:
|
|
|
|
! --------------------------------
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! Subroutine broadcast a string variable from the master process to all
|
|
|
|
! other processes.
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine bcast_string_variable(sbuf, iret)
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_character, mpi_success
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
character(len=*), intent(inout) :: sbuf
|
|
|
|
integer , intent(out) :: iret
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_bcast(sbuf, len(sbuf), mpi_character, 0, comm3d, iret)
|
|
|
|
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not broadcast a string variable!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine bcast_string_variable
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_MINIMUM_INTEGER:
|
|
|
|
! ---------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine finds the minimum value among the integer values from all
|
|
|
|
! processes.
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_minimum_integer(ibuf, iret)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_integer, mpi_min, mpi_success
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer, intent(inout) :: ibuf
|
|
|
|
integer, intent(out) :: iret
|
2008-12-31 12:02:36 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer :: tbuf
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(ibuf, tbuf, 1, mpi_integer, mpi_min, comm3d, iret)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
ibuf = tbuf
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the minimum value!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine reduce_minimum_integer
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_MINIMUM_REAL:
|
|
|
|
! ------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine finds the minimum value among the real values from all processes.
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_minimum_real(rbuf, iret)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_min, mpi_success
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
2008-12-31 12:02:36 -06:00
|
|
|
implicit none
|
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine arguments
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real , intent(inout) :: rbuf
|
|
|
|
integer, intent(out) :: iret
|
2008-12-31 12:02:36 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real :: tbuf
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(rbuf, tbuf, 1, mpi_real8, mpi_min, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
rbuf = tbuf
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the minimum value!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_minimum_real
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_MAXIMUM_INTEGER:
|
|
|
|
! ---------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine find the maximum value among the integer values from all
|
|
|
|
! processes.
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_maximum_integer(ibuf, iret)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_integer, mpi_max, mpi_success
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine arguments
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer, intent(inout) :: ibuf
|
|
|
|
integer, intent(out) :: iret
|
2008-12-31 12:02:36 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer :: tbuf
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(ibuf, tbuf, 1, mpi_integer, mpi_max, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
ibuf = tbuf
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the maximum value!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine reduce_maximum_integer
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_MAXIMUM_REAL:
|
|
|
|
! ------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine find the maximum value among the values from all processes.
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_maximum_real(rbuf, iret)
|
2009-01-08 00:08:52 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_max, mpi_success
|
2009-01-08 00:08:52 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
real , intent(inout) :: rbuf
|
|
|
|
integer, intent(out) :: iret
|
2009-01-08 00:08:52 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real :: tbuf
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(rbuf, tbuf, 1, mpi_real8, mpi_max, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
rbuf = tbuf
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the maximum value!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2009-01-08 00:08:52 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_maximum_real
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_SUM_INTEGER:
|
|
|
|
! -----------------------------
|
|
|
|
!
|
|
|
|
! Subroutine finds the sum from all integer values from all processes.
|
2011-03-02 14:56:22 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_sum_integer(ibuf, iret)
|
2011-03-02 14:56:22 -03:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_integer, mpi_sum, mpi_success
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
2011-03-02 14:56:22 -03:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine arguments
|
2011-03-02 14:56:22 -03:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer, intent(inout) :: ibuf
|
|
|
|
integer, intent(out) :: iret
|
2011-03-02 14:56:22 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer :: tbuf
|
2011-03-02 14:56:22 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(ibuf, tbuf, 1, mpi_integer, mpi_sum, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
ibuf = tbuf
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the maximum value!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2011-03-02 14:56:22 -03:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_sum_integer
|
2011-03-02 14:56:22 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_SUM_REAL:
|
|
|
|
! --------------------------
|
|
|
|
!
|
|
|
|
! Subroutine sums the values from all processes.
|
2011-03-07 00:28:58 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_sum_real(rbuf, iret)
|
2011-03-07 00:28:58 -03:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_sum, mpi_success
|
2011-03-07 00:28:58 -03:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
2011-03-07 00:28:58 -03:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
real , intent(inout) :: rbuf
|
|
|
|
integer, intent(out) :: iret
|
2011-03-07 00:28:58 -03:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real :: tbuf
|
2011-03-07 00:28:58 -03:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(rbuf, tbuf, 1, mpi_real8, mpi_sum, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
rbuf = tbuf
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not sum the values from all processes!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2011-03-07 00:28:58 -03:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_sum_real
|
2011-03-07 00:28:58 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_MINIMUM_REAL_ARRAY:
|
|
|
|
! ------------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine find the minimum value for each array element among the
|
|
|
|
! corresponding values from all processes.
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_minimum_real_array(n, rbuf, iret)
|
2009-01-08 00:08:52 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_min, mpi_success
|
2009-01-08 00:08:52 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
|
|
|
integer , intent(in) :: n
|
2012-07-22 19:01:27 -03:00
|
|
|
real , dimension(n), intent(inout) :: rbuf
|
|
|
|
integer , intent(out) :: iret
|
2009-01-08 00:08:52 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real(kind=8), dimension(n) :: tbuf
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(rbuf, tbuf, n, mpi_real8, mpi_min, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
rbuf(:) = tbuf(:)
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the minima for all array elements!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2009-01-08 00:08:52 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_minimum_real_array
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_MAXIMUM_REAL_ARRAY:
|
|
|
|
! ------------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine find the maximum value for each array element among the
|
|
|
|
! corresponding values from all processes.
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_maximum_real_array(n, rbuf, iret)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_max, mpi_success
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
2008-12-31 12:02:36 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine arguments
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
|
|
|
integer , intent(in) :: n
|
2012-07-22 19:01:27 -03:00
|
|
|
real , dimension(n), intent(inout) :: rbuf
|
|
|
|
integer , intent(out) :: iret
|
2008-12-31 12:02:36 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real(kind=8), dimension(n) :: tbuf
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2008-12-31 12:02:36 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(rbuf, tbuf, n, mpi_real8, mpi_max, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
rbuf(:) = tbuf(:)
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the maxima for all array elements!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2008-12-31 12:02:36 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_maximum_real_array
|
2009-01-02 20:18:57 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_SUM_INTEGER_ARRAY:
|
|
|
|
! -----------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine sums the values for each array element from the corresponding
|
|
|
|
! values from all processes.
|
2009-01-02 20:18:57 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_sum_integer_array(n, ibuf, iret)
|
2009-01-02 20:18:57 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_integer, mpi_sum, mpi_success
|
2009-01-02 20:18:57 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
2009-01-02 20:18:57 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: n
|
|
|
|
integer, dimension(n), intent(inout) :: ibuf
|
|
|
|
integer , intent(out) :: iret
|
2009-01-02 20:18:57 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer, dimension(n) :: tbuf
|
2009-01-02 20:18:57 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2009-01-02 20:18:57 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(ibuf, tbuf, n, mpi_integer, mpi_sum, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
ibuf(:) = tbuf(:)
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the maxima for all array elements!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2009-01-02 20:18:57 -06:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_sum_integer_array
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_SUM_REAL_ARRAY:
|
|
|
|
! --------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine sums the values for each array element from the corresponding
|
|
|
|
! values from all processes.
|
2010-12-03 17:58:20 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_sum_real_array(n, rbuf, iret)
|
2010-12-03 17:58:20 -02:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_sum, mpi_success
|
2010-12-03 17:58:20 -02:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
2010-12-03 17:58:20 -02:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: n
|
|
|
|
real , dimension(n), intent(inout) :: rbuf
|
|
|
|
integer , intent(out) :: iret
|
2010-12-03 17:58:20 -02:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real(kind=8), dimension(n) :: tbuf
|
2010-12-03 17:58:20 -02:00
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_allreduce(rbuf, tbuf, n, mpi_real8, mpi_sum, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
rbuf(:) = tbuf(:)
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the maxima for all array elements!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2010-12-03 17:58:20 -02:00
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_sum_real_array
|
2010-12-03 17:58:20 -02:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! subroutine REDUCE_SUM_COMPLEX_ARRAY:
|
|
|
|
! -----------------------------------
|
|
|
|
!
|
|
|
|
! Subroutine sums the values for each array element from the corresponding
|
|
|
|
! complex values from all processes.
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
subroutine reduce_sum_complex_array(n, cbuf, iret)
|
2009-01-08 00:08:52 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_sum, mpi_success
|
2009-01-08 00:08:52 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: n
|
|
|
|
complex, dimension(n), intent(inout) :: cbuf
|
|
|
|
integer , intent(out) :: iret
|
2009-01-08 00:08:52 -06:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
real(kind=8), dimension(n) :: rbuf, ibuf, tbuf
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2010-10-13 03:32:10 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
2009-01-08 00:08:52 -06:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
tbuf(:) = real(cbuf(:))
|
|
|
|
call mpi_allreduce(tbuf, rbuf, n, mpi_real8, mpi_sum, comm3d, iret)
|
|
|
|
tbuf(:) = aimag(cbuf(:))
|
|
|
|
call mpi_allreduce(tbuf, ibuf, n, mpi_real8, mpi_sum, comm3d, iret)
|
|
|
|
|
|
|
|
! substitute the result
|
|
|
|
!
|
|
|
|
cbuf(1:n) = cmplx(rbuf(1:n), ibuf(1:n))
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not find the maxima for all array elements!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
|
|
|
|
2009-01-08 00:08:52 -06:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
end subroutine reduce_sum_complex_array
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine SEND_REAL_ARRAY:
|
|
|
|
! --------------------------
|
|
|
|
!
|
|
|
|
! Subroutine sends an arrays of real values to another process.
|
|
|
|
!
|
|
|
|
! Arguments:
|
2010-12-08 11:37:56 -02:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! n - the number of array elements;
|
|
|
|
! dst - the ID of the destination process;
|
|
|
|
! tag - the tag identifying this operation;
|
|
|
|
! rbuf - the real array to send;
|
|
|
|
! iret - the result flag identifying if the operation was successful;
|
2010-12-08 11:37:56 -02:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine send_real_array(n, dst, tag, rbuf, iret)
|
2010-12-08 11:37:56 -02:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_success
|
2010-12-08 11:37:56 -02:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
! local variables are not implicit by default
|
2010-12-08 11:37:56 -02:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: n, dst, tag
|
|
|
|
real , dimension(n), intent(in) :: rbuf
|
|
|
|
integer , intent(out) :: iret
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
! start time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_send(rbuf, n, mpi_real8, dst, tag, comm3d, iret)
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not send the real array to another process!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
|
|
|
|
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine send_real_array
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
! subroutine RECEIVE_REAL_ARRAY:
|
|
|
|
! -----------------------------
|
|
|
|
!
|
|
|
|
! Subroutine receives an arrays of real values from another process.
|
|
|
|
!
|
|
|
|
! Arguments:
|
|
|
|
!
|
|
|
|
! n - the number of array elements;
|
|
|
|
! src - the ID of the source process;
|
|
|
|
! tag - the tag identifying this operation;
|
|
|
|
! rbuf - the received real array;
|
|
|
|
! iret - the result flag identifying if the operation was successful;
|
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
subroutine receive_real_array(n, src, tag, rbuf, iret)
|
|
|
|
|
|
|
|
! include external procedures and variables
|
|
|
|
!
|
|
|
|
use mpi, only : mpi_real8, mpi_success, mpi_status_size
|
|
|
|
|
|
|
|
! local variables are not implicit by default
|
|
|
|
!
|
|
|
|
implicit none
|
|
|
|
|
|
|
|
! subroutine arguments
|
|
|
|
!
|
|
|
|
integer , intent(in) :: n, src, tag
|
|
|
|
real , dimension(n), intent(out) :: rbuf
|
|
|
|
integer , intent(out) :: iret
|
2010-12-08 11:37:56 -02:00
|
|
|
|
|
|
|
! local variables
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
integer :: status(mpi_status_size)
|
|
|
|
!
|
|
|
|
!-------------------------------------------------------------------------------
|
2010-12-08 11:37:56 -02:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
! start time accounting for the MPI communication
|
2010-12-08 11:37:56 -02:00
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
call start_timer(imc)
|
|
|
|
|
|
|
|
call mpi_recv(rbuf, n, mpi_real8, src, tag, comm3d, status, iret)
|
|
|
|
|
|
|
|
! check if the operation was successful
|
|
|
|
!
|
|
|
|
if (iret .ne. mpi_success .and. master) then
|
|
|
|
write(*,*) 'The MPI could not send the real array to another process!'
|
|
|
|
write(*,*)
|
|
|
|
end if
|
|
|
|
|
|
|
|
! stop time accounting for the MPI communication
|
|
|
|
!
|
|
|
|
call stop_timer(imc)
|
2010-12-08 11:37:56 -02:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
!-------------------------------------------------------------------------------
|
|
|
|
!
|
|
|
|
end subroutine receive_real_array
|
2012-07-28 12:29:44 -03:00
|
|
|
!
|
|
|
|
!===============================================================================
|
|
|
|
!!
|
|
|
|
!!*** PRIVATE SUBROUTINES ****************************************************
|
|
|
|
!!
|
|
|
|
!===============================================================================
|
|
|
|
!
|
2012-07-22 19:01:27 -03:00
|
|
|
#endif /* MPI */
|
2008-12-19 17:33:29 -06:00
|
|
|
|
2012-07-22 19:01:27 -03:00
|
|
|
!===============================================================================
|
|
|
|
!
|
|
|
|
end module mpitools
|