DAMASK_EICMD/src/parallelization.f90

173 lines
6.5 KiB
Fortran
Raw Normal View History

!--------------------------------------------------------------------------------------------------
!> @author Martin Diehl, Max-Planck-Institut für Eisenforschung GmbH
!> @brief Inquires variables related to parallelization (openMP, MPI)
!--------------------------------------------------------------------------------------------------
module parallelization
2020-09-19 14:20:32 +05:30
use, intrinsic :: ISO_fortran_env, only: &
OUTPUT_UNIT
2021-07-08 18:35:01 +05:30
#ifdef PETSC
#include <petsc/finclude/petscsys.h>
use PETScSys
#if (PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>14) && !defined(PETSC_HAVE_MPI_F90MODULE_VISIBILITY)
use MPI_f08
#endif
!$ use OMP_LIB
#endif
2020-09-19 14:20:32 +05:30
use prec
implicit none
private
#ifndef PETSC
integer, parameter, public :: &
MPI_INTEGER_KIND = pI64
integer(MPI_INTEGER_KIND), parameter, public :: &
worldrank = 0_MPI_INTEGER_KIND, & !< MPI dummy worldrank
worldsize = 1_MPI_INTEGER_KIND !< MPI dummy worldsize
#else
integer(MPI_INTEGER_KIND), protected, public :: &
worldrank = 0_MPI_INTEGER_KIND, & !< MPI worldrank (/=0 for MPI simulations only)
worldsize = 1_MPI_INTEGER_KIND !< MPI worldsize (/=1 for MPI simulations only)
#endif
2020-09-13 13:49:38 +05:30
2021-07-27 12:24:17 +05:30
#ifndef PETSC
2022-01-13 16:33:22 +05:30
public :: parallelization_bcast_str
2021-07-27 12:24:17 +05:30
contains
subroutine parallelization_bcast_str(string)
character(len=:), allocatable, intent(inout) :: string
2021-07-27 12:24:17 +05:30
end subroutine parallelization_bcast_str
#else
public :: &
2021-07-27 12:24:17 +05:30
parallelization_init, &
parallelization_bcast_str
contains
!--------------------------------------------------------------------------------------------------
2021-05-23 03:40:46 +05:30
!> @brief Initialize shared memory (openMP) and distributed memory (MPI) parallelization.
!--------------------------------------------------------------------------------------------------
subroutine parallelization_init
integer(MPI_INTEGER_KIND) :: err_MPI, typeSize, version, subversion, devNull
character(len=4) :: rank_str
character(len=MPI_MAX_LIBRARY_VERSION_STRING) :: MPI_library_version
!$ integer :: got_env, threadLevel
!$ integer(pI32) :: OMP_NUM_THREADS
!$ character(len=6) NumThreadsString
PetscErrorCode :: err_PETSc
2020-09-13 13:49:38 +05:30
#ifdef _OPENMP
! If openMP is enabled, check if the MPI libary supports it and initialize accordingly.
! Otherwise, the first call to PETSc will do the initialization.
call MPI_Init_Thread(MPI_THREAD_FUNNELED,threadLevel,err_MPI)
if (err_MPI /= 0_MPI_INTEGER_KIND) error stop 'MPI init failed'
2021-05-23 03:40:46 +05:30
if (threadLevel<MPI_THREAD_FUNNELED) error stop 'MPI library does not support OpenMP'
2020-09-13 13:49:38 +05:30
#endif
#if defined(DEBUG)
call PetscInitialize(PETSC_NULL_CHARACTER,err_PETSc)
#else
call PetscInitializeNoArguments(err_PETSc)
#endif
CHKERRQ(err_PETSc)
2020-09-13 13:49:38 +05:30
#if defined(DEBUG) && defined(__INTEL_COMPILER)
call PetscSetFPTrap(PETSC_FP_TRAP_ON,err_PETSc)
2020-11-11 16:16:12 +05:30
#else
call PetscSetFPTrap(PETSC_FP_TRAP_OFF,err_PETSc)
2020-11-11 16:16:12 +05:30
#endif
CHKERRQ(err_PETSc)
2020-11-11 16:16:12 +05:30
call MPI_Comm_rank(MPI_COMM_WORLD,worldrank,err_MPI)
if (err_MPI /= 0_MPI_INTEGER_KIND) &
error stop 'Could not determine worldrank'
#ifdef LOGFILE
write(rank_str,'(i4.4)') worldrank
open(OUTPUT_UNIT,file='log.'//rank_str,status='replace',encoding='UTF-8')
#else
if (worldrank /= 0) then
close(OUTPUT_UNIT) ! disable output
open(OUTPUT_UNIT,file='/dev/null',status='replace') ! close() alone will leave some temp files in cwd
else
open(OUTPUT_UNIT,encoding='UTF-8') ! for special characters in output
endif
#endif
print'(/,1x,a)', '<<<+- parallelization init -+>>>'
call MPI_Get_library_version(MPI_library_version,devNull,err_MPI)
print'(/,1x,a)', trim(MPI_library_version)
call MPI_Get_version(version,subversion,err_MPI)
print'(1x,a,i0,a,i0)', 'MPI standard: ',version,'.',subversion
#ifdef _OPENMP
print'(1x,a,i0)', 'OpenMP version: ',openmp_version
#endif
call MPI_Comm_size(MPI_COMM_WORLD,worldsize,err_MPI)
if (err_MPI /= 0_MPI_INTEGER_KIND) &
error stop 'Could not determine worldsize'
if (worldrank == 0) print'(/,1x,a,i0)', 'MPI processes: ',worldsize
2020-09-13 13:49:38 +05:30
call MPI_Type_size(MPI_INTEGER,typeSize,err_MPI)
if (err_MPI /= 0_MPI_INTEGER_KIND) &
error stop 'Could not determine size of MPI_INTEGER'
if (typeSize*8_MPI_INTEGER_KIND /= int(bit_size(0),MPI_INTEGER_KIND)) &
error stop 'Mismatch between MPI_INTEGER and DAMASK default integer'
2020-09-13 13:49:38 +05:30
call MPI_Type_size(MPI_INTEGER8,typeSize,err_MPI)
if (err_MPI /= 0) &
error stop 'Could not determine size of MPI_INTEGER8'
if (typeSize*8_MPI_INTEGER_KIND /= int(bit_size(0_pI64),MPI_INTEGER_KIND)) &
error stop 'Mismatch between MPI_INTEGER8 and DAMASK pI64'
2022-01-12 21:26:24 +05:30
call MPI_Type_size(MPI_DOUBLE,typeSize,err_MPI)
if (err_MPI /= 0_MPI_INTEGER_KIND) &
error stop 'Could not determine size of MPI_DOUBLE'
if (typeSize*8_MPI_INTEGER_KIND /= int(storage_size(0.0_pReal),MPI_INTEGER_KIND)) &
error stop 'Mismatch between MPI_DOUBLE and DAMASK pReal'
2020-09-13 13:49:38 +05:30
2021-01-13 17:18:20 +05:30
!$ call get_environment_variable(name='OMP_NUM_THREADS',value=NumThreadsString,STATUS=got_env)
2020-09-13 13:49:38 +05:30
!$ if(got_env /= 0) then
!$ print'(1x,a)', 'Could not get $OMP_NUM_THREADS, using default'
!$ OMP_NUM_THREADS = 4_pI32
!$ else
2021-01-13 17:18:20 +05:30
!$ read(NumThreadsString,'(i6)') OMP_NUM_THREADS
!$ if (OMP_NUM_THREADS < 1_pI32) then
!$ print'(1x,a)', 'Invalid OMP_NUM_THREADS: "'//trim(NumThreadsString)//'", using default'
!$ OMP_NUM_THREADS = 4_pI32
!$ endif
!$ endif
!$ print'(1x,a,i0)', 'OMP_NUM_THREADS: ',OMP_NUM_THREADS
2021-01-13 17:18:20 +05:30
!$ call omp_set_num_threads(OMP_NUM_THREADS)
2020-09-13 13:49:38 +05:30
end subroutine parallelization_init
2021-07-27 12:24:17 +05:30
!--------------------------------------------------------------------------------------------------
!> @brief Broadcast a string from process 0.
!--------------------------------------------------------------------------------------------------
subroutine parallelization_bcast_str(string)
character(len=:), allocatable, intent(inout) :: string
integer(MPI_INTEGER_KIND) :: strlen, err_MPI
2021-07-27 12:24:17 +05:30
if (worldrank == 0) strlen = len(string,MPI_INTEGER_KIND)
call MPI_Bcast(strlen,1_MPI_INTEGER_KIND,MPI_INTEGER,0_MPI_INTEGER_KIND,MPI_COMM_WORLD, err_MPI)
2021-07-27 12:24:17 +05:30
if (worldrank /= 0) allocate(character(len=strlen)::string)
call MPI_Bcast(string,strlen,MPI_CHARACTER,0_MPI_INTEGER_KIND,MPI_COMM_WORLD, err_MPI)
2021-07-27 12:24:17 +05:30
end subroutine parallelization_bcast_str
#endif
end module parallelization