!!****************************************************************************** !! !! This file is part of the AMUN source code, a program to perform !! Newtonian or relativistic magnetohydrodynamical simulations on uniform or !! adaptive mesh. !! !! Copyright (C) 2008-2020 Grzegorz Kowal !! !! 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. !! !! This program is distributed in the hope that it will be useful, !! 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 !! along with this program. If not, see . !! !!****************************************************************************** !! !! module: MPITOOLS !! !! This module provides wrapper subroutines handling the parallel execution !! with the Message Passing Interface protocol. !! !! !!****************************************************************************** ! module mpitools ! include external subroutines ! #ifdef MPI use mpi_f08 #endif /* MPI */ use timers, only : set_timer, start_timer, stop_timer ! module variables are not implicit by default ! implicit none ! subroutine interfaces ! #ifdef MPI interface reduce_minimum module procedure reduce_minimum_double_array end interface interface reduce_maximum module procedure reduce_maximum_integer module procedure reduce_maximum_double module procedure reduce_maximum_double_array end interface interface reduce_sum module procedure reduce_sum_integer_array module procedure reduce_sum_double_array module procedure reduce_sum_complex_array end interface #endif /* MPI */ ! timer indices ! integer, save :: imi, imc #ifdef PROFILE integer, save :: imb, imm, ims, imr, ime #endif /* PROFILE */ ! MPI global variables ! integer(kind=4), save :: nproc, nprocs, npmax, npairs logical , save :: master = .true. ! allocatable array for processor pairs ! integer(kind=4), dimension(:,:), allocatable, save :: pairs ! by default everything is private ! private ! declare public subroutines ! public :: initialize_mpitools, finalize_mpitools public :: check_status #ifdef MPI public :: reduce_minimum, reduce_maximum, reduce_sum public :: send_array, receive_array public :: exchange_arrays #endif /* MPI */ ! declare public variables ! public :: master, nproc, nprocs, npmax, npairs, pairs !- - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - ! contains ! !=============================================================================== !! !!*** PUBLIC SUBROUTINES ***************************************************** !! !=============================================================================== ! ! subroutine INITIALIZE_MPITOOLS: ! ------------------------------ ! ! Subroutine initializes the MPITOOLS modules. ! ! Arguments: ! ! status - the return value; if it is 0 everything went successfully, ! otherwise there was a problem; ! !=============================================================================== ! subroutine initialize_mpitools(status) ! include external procedures and variables ! use iso_fortran_env, only : error_unit ! local variables are not implicit by default ! implicit none ! subroutine arguments ! integer, intent(out) :: status #ifdef MPI ! local variables ! integer :: mprocs, i, j, l, n integer :: ierror ! allocatable array for processors order ! integer(kind=4), dimension(:), allocatable :: procs #endif /* MPI */ ! local parameters ! character(len=*), parameter :: loc = 'MPITOOLS::initialize_mpitools()' ! !------------------------------------------------------------------------------- ! #ifdef MPI call set_timer('MPI initialization' , imi) call set_timer('MPI communication' , imc) #ifdef PROFILE call set_timer('mpitools:: broadcast', imb) call set_timer('mpitools:: reduce' , imm) call set_timer('mpitools:: send' , ims) call set_timer('mpitools:: receive' , imr) call set_timer('mpitools:: exchange' , ime) #endif /* PROFILE */ call start_timer(imi) #endif /* MPI */ status = 0 nproc = 0 nprocs = 1 npmax = 0 npairs = 0 #ifdef MPI ! initialize the MPI parallelization ! call MPI_Init(ierror) if (ierror == MPI_SUCCESS) then call MPI_Comm_size(MPI_COMM_WORLD, nprocs, ierror) if (ierror == MPI_SUCCESS) then call MPI_Comm_rank(MPI_COMM_WORLD, nproc, ierror) if (ierror == MPI_SUCCESS) then ! set the master flag ! master = nproc == 0 ! allocate and fill up the array of processes and process pairs ! npmax = nprocs - 1 mprocs = nprocs + mod(nprocs, 2) npairs = nprocs * npmax / 2 allocate(procs(mprocs), pairs(2 * npairs, 2), stat = status) if (status == 0) then procs(:) = (/(l, l = 0, mprocs - 1)/) n = 0 do l = 1, mprocs - 1 do i = 1, mprocs / 2 j = mprocs - i + 1 if (procs(i) < nprocs .and. procs(j) < nprocs) then n = n + 1 pairs(n,1:2) = (/ procs(i), procs(j) /) end if ! max(procs(i), procs(j)) < nprocs end do ! i = 1, mprocs / 2 procs(2:mprocs) = cshift(procs(2:mprocs), -1) end do ! l = 1, mprocs - 1 pairs(npairs+1:2*npairs,1:2) = pairs(1:npairs,2:1:-1) deallocate(procs, stat = status) end if ! allocate else write(error_unit,"('[', a, ']: ', a)") trim(loc), & "Could not get the MPI process ID!" status = 1 end if else write(error_unit,"('[', a, ']: ', a)") trim(loc), & "Could not get the number of MPI processes!" status = 1 end if else write(error_unit,"('[', a, ']: ', a)") trim(loc), & "Could not initialize the MPI interface!" status = 1 end if call stop_timer(imi) #endif /* MPI */ !------------------------------------------------------------------------------- ! end subroutine initialize_mpitools ! !=============================================================================== ! ! subroutine FINALIZE_MPITOOLS: ! ---------------------------- ! ! Subroutine finalizes the MPITOOLS modules. ! ! Arguments: ! ! status - the return value; if it is 0 everything went successfully, ! otherwise there was a problem; ! !=============================================================================== ! subroutine finalize_mpitools(status) ! include external procedures and variables ! use iso_fortran_env, only : error_unit ! local variables are not implicit by default ! implicit none ! subroutine arguments ! integer, intent(out) :: status #ifdef MPI ! local variables ! integer :: ierror #endif /* MPI */ ! local parameters ! character(len=*), parameter :: loc = 'MPITOOLS::finalize_mpitools()' ! !------------------------------------------------------------------------------- ! status = 0 #ifdef MPI call start_timer(imi) ! deallocate space used for processor pairs ! if (allocated(pairs)) deallocate(pairs, stat = status) ! initialize the MPI interface ! call MPI_Finalize(ierror) if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "Could not finalize the MPI interface!" status = 1 end if call stop_timer(imi) #endif /* MPI */ !------------------------------------------------------------------------------- ! end subroutine finalize_mpitools ! !=============================================================================== ! ! subroutine CHECK_STATUS: ! ----------------------- ! ! Subroutine calculates the logical OR for input values from all MPI ! processes, if MPI is used, otherwise, just returns the input value. ! ! Arguments: ! ! ibuf - the logical buffer; ! !=============================================================================== ! logical function check_status(ibuf) result(obuf) ! include external procedures and variables ! use iso_fortran_env, only : error_unit ! local variables are not implicit by default ! implicit none ! subroutine arguments ! logical, intent(in) :: ibuf #ifdef MPI ! local variables ! integer :: ierror ! local parameters ! character(len=*), parameter :: loc = 'MPITOOLS::check_status()' #endif /* MPI */ ! !------------------------------------------------------------------------------- ! #ifdef MPI call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(ibuf, obuf, 1, & MPI_LOGICAL, MPI_LOR, MPI_COMM_WORLD, ierror) if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of logical buffer failed!" end if #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ call stop_timer(imc) #else /* MPI */ ! no MPI, so just copy the input to output ! obuf = ibuf #endif /* MPI */ !------------------------------------------------------------------------------- ! end function check_status ! #ifdef MPI ! !=============================================================================== ! ! subroutine SEND_ARRAY: ! --------------------- ! ! Subroutine sends an arrays of real values to another process. ! ! Arguments: ! ! dst - the ID of the destination process; ! tag - the tag identifying this operation; ! buf - the buffer of real values to send; ! !=============================================================================== ! subroutine send_array(dst, tag, buf) use iso_fortran_env, only : error_unit implicit none integer , intent(in) :: dst, tag real(kind=8), dimension(..), intent(in) :: buf integer :: ierror character(len=*), parameter :: loc = 'MPITOOLS::send_array()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(ims) #endif /* PROFILE */ call MPI_Send(buf, size(buf), MPI_REAL8, dst, tag, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(ims) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', 2(a, i9))") trim(loc) & , "Could not send real array from ", nproc, " to ", dst end if call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine send_array ! !=============================================================================== ! ! subroutine RECEIVE_ARRAY: ! ------------------------ ! ! Subroutine receives an arrays of real values from another process. ! ! Arguments: ! ! src - the ID of the source process; ! tag - the tag identifying this operation; ! buf - the received real array; ! !=============================================================================== ! subroutine receive_array(src, tag, buf) use iso_fortran_env, only : error_unit implicit none integer , intent(in) :: src, tag real(kind=8), dimension(..), intent(out) :: buf integer :: ierror character(len=*), parameter :: loc = 'MPITOOLS::receive_array()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imr) #endif /* PROFILE */ call MPI_Recv(buf, size(buf), MPI_REAL8, src, tag, & MPI_COMM_WORLD, MPI_STATUS_IGNORE, ierror) #ifdef PROFILE call stop_timer(imr) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', 2(a, i9))") trim(loc) & , "Could not receive real array from ", src, " to ", nproc end if call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine receive_array ! !=============================================================================== ! ! subroutine EXCHANGE_ARRAYS: ! -------------------------- ! ! Subroutine exchanges real data buffers between two processes. ! ! Arguments: ! ! proc - the remote process number to which send the buffer sbuf, ! and from which receive the buffer rbuf; ! tag - the tag identifying the send operation; ! sbuf - the real array buffer to send; ! rbuf - the real array buffer to receive; ! !=============================================================================== ! subroutine exchange_arrays(proc, tag, sbuf, rbuf) use iso_fortran_env, only : error_unit implicit none integer , intent(in) :: proc, tag real(kind=8), dimension(..), intent(in) :: sbuf real(kind=8), dimension(..), intent(in) :: rbuf integer :: ierror character(len=*), parameter :: loc = 'MPITOOLS::exchange_arrays()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(ime) #endif /* PROFILE */ call MPI_Sendrecv(sbuf, size(sbuf), MPI_REAL8, proc, tag, & rbuf, size(rbuf), MPI_REAL8, proc, tag, & MPI_COMM_WORLD, MPI_STATUS_IGNORE, ierror) #ifdef PROFILE call stop_timer(ime) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', 2(a, i9),'.')") trim(loc) & , "Could not exchange real data buffers between " & , proc, "and", nproc end if call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine exchange_arrays ! !=============================================================================== !! !!*** PRIVATE SUBROUTINES **************************************************** !! !=============================================================================== ! !=============================================================================== ! ! subroutine REDUCE_MINIMUM_DOUBLE_ARRAY: ! -------------------------------------- ! ! Subroutine find the minimum value for each double precision array element ! among the corresponding values from all processes. ! ! Argument: ! ! buf - a buffer to be reduced; ! !=============================================================================== ! subroutine reduce_minimum_double_array(buf) use iso_fortran_env, only : error_unit implicit none real(kind=8), dimension(:), intent(inout) :: buf integer :: ierror real(kind=8), dimension(size(buf)) :: tmp character(len=*), parameter :: loc = 'MPITOOLS::reduce_minimum_double_array()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(buf, tmp, size(buf), & MPI_REAL8, MPI_MIN, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of a real array failed!" end if buf(:) = tmp(:) call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine reduce_minimum_double_array ! !=============================================================================== ! ! subroutine REDUCE_MAXIMUM_INTEGER: ! --------------------------------- ! ! Subroutine find the maximum value among the integer values from all ! processes. ! ! Argument: ! ! buf - a buffer to be reduced; ! !=============================================================================== ! subroutine reduce_maximum_integer(buf) use iso_fortran_env, only : error_unit implicit none integer, intent(inout) :: buf integer :: ierror integer :: tmp character(len=*), parameter :: loc = 'MPITOOLS::reduce_maximum_integer()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(buf, tmp, 1, & MPI_INTEGER, MPI_MAX, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of an integer value failed!" end if buf = tmp call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine reduce_maximum_integer ! !=============================================================================== ! ! subroutine REDUCE_MAXIMUM_DOUBLE: ! -------------------------------- ! ! Subroutine find the maximum value among the double precision values ! from all processes. ! ! Argument: ! ! buf - a buffer to be reduced; ! !=============================================================================== ! subroutine reduce_maximum_double(buf) use iso_fortran_env, only : error_unit implicit none real(kind=8), intent(inout) :: buf integer :: ierror real(kind=8) :: tmp character(len=*), parameter :: loc = 'MPITOOLS::reduce_maximum_double()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(buf, tmp, 1, MPI_REAL8, MPI_MAX, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of a real value failed!" end if buf = tmp call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine reduce_maximum_double ! !=============================================================================== ! ! subroutine REDUCE_MAXIMUM_DOUBLE_ARRAY: ! -------------------------------------- ! ! Subroutine find the maximum value for each double plrecision array element ! among the corresponding values from all processes. ! ! Argument: ! ! buf - a buffer to be reduced; ! !=============================================================================== ! subroutine reduce_maximum_double_array(buf) use iso_fortran_env, only : error_unit implicit none real(kind=8), dimension(:), intent(inout) :: buf integer :: ierror real(kind=8), dimension(size(buf)) :: tmp character(len=*), parameter :: loc = 'MPITOOLS::reduce_maximum_double_array()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(buf, tmp, size(buf), & MPI_REAL8, MPI_MAX, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of a real array failed!" end if buf(:) = tmp(:) call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine reduce_maximum_double_array ! !=============================================================================== ! ! subroutine REDUCE_SUM_INTEGER_ARRAY: ! ----------------------------------- ! ! Subroutine sums the values for each array element from the corresponding ! values from all processes. ! ! Argument: ! ! buf - a buffer to be reduced; ! !=============================================================================== ! subroutine reduce_sum_integer_array(buf) use iso_fortran_env, only : error_unit implicit none integer, dimension(:), intent(inout) :: buf integer :: ierror integer, dimension(size(buf)) :: tmp character(len=*), parameter :: loc = 'MPITOOLS::reduce_sum_integer_array()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(buf, tmp, size(buf), & MPI_INTEGER, MPI_SUM, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of an integer array failed!" end if buf(:) = tmp(:) call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine reduce_sum_integer_array ! !=============================================================================== ! ! subroutine REDUCE_SUM_DOUBLE_ARRAY: ! ---------------------------------- ! ! Subroutine sums the values for each double precision array element from ! the corresponding values from all processes. ! ! Argument: ! ! buf - a buffer to be reduced; ! !=============================================================================== ! subroutine reduce_sum_double_array(buf) use iso_fortran_env, only : error_unit implicit none real(kind=8), dimension(:), intent(inout) :: buf integer :: ierror real(kind=8), dimension(size(buf)) :: tmp character(len=*), parameter :: loc = 'MPITOOLS::reduce_sum_double_array()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(buf, tmp, size(buf), & MPI_REAL8, MPI_SUM, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of a real array failed!" end if buf(:) = tmp(:) call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine reduce_sum_double_array ! !=============================================================================== ! ! subroutine REDUCE_SUM_COMPLEX_ARRAY: ! ----------------------------------- ! ! Subroutine sums the values for each array element from the corresponding ! complex values from all processes. ! ! Argument: ! ! buf - a buffer to be reduced; ! !=============================================================================== ! subroutine reduce_sum_complex_array(buf) use iso_fortran_env, only : error_unit implicit none complex(kind=8), dimension(:,:), intent(inout) :: buf integer :: ierror complex(kind=8), dimension(size(buf,1),size(buf,2)) :: tmp character(len=*), parameter :: loc = 'MPITOOLS::reduce_sum_complex_array()' ! !------------------------------------------------------------------------------- ! call start_timer(imc) #ifdef PROFILE call start_timer(imm) #endif /* PROFILE */ call MPI_Allreduce(buf, tmp, size(buf), & MPI_DOUBLE_COMPLEX, MPI_SUM, MPI_COMM_WORLD, ierror) #ifdef PROFILE call stop_timer(imm) #endif /* PROFILE */ if (ierror /= MPI_SUCCESS) then write(error_unit,"('[', a, ']: ', a)") trim(loc), & "MPI_Allreduce of a complex array failed!" end if buf(:,:) = tmp(:,:) call stop_timer(imc) !------------------------------------------------------------------------------- ! end subroutine reduce_sum_complex_array #endif /* MPI */ !=============================================================================== ! end module mpitools