44#if CAF_ENABLED || MPI_ENABLED
79 logical(
LK) :: first
= .false._LK
80 logical(
LK) :: extra
= .false._LK
81 logical(
LK) :: leader
= .false._LK
84 logical(
LK) :: rooter
= .false._LK
134 integer(
IK) :: count
= -huge(
1_IK)
135 integer(
IK) ::
id = -huge(
1_IK)
137 character(:,
SK),
allocatable :: label
234 PURE module subroutine setForkJoinScaling_RK5(conProb, seqSecTime, parSecTime, comSecTime, scaling, numproc, scalingMaxVal, scalingMaxLoc, scalingMinLen)
235#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
236 !DEC$ ATTRIBUTES DLLEXPORT :: setForkJoinScaling_RK5
239 integer(IK) ,
intent(out) :: scalingMaxLoc
240 integer(IK) ,
intent(in) ,
optional :: scalingMinLen
241 integer(IK) ,
intent(out) ,
allocatable :: numproc(:)
242 real(RKG) ,
intent(out) ,
allocatable :: scaling(:)
243 real(RKG) ,
intent(in) :: conProb, seqSecTime, parSecTime, comSecTime
244 real(RKG) ,
intent(out) :: scalingMaxVal
249 PURE module subroutine setForkJoinScaling_RK4(conProb, seqSecTime, parSecTime, comSecTime, scaling, numproc, scalingMaxVal, scalingMaxLoc, scalingMinLen)
250#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
251 !DEC$ ATTRIBUTES DLLEXPORT :: setForkJoinScaling_RK4
254 integer(IK) ,
intent(out) :: scalingMaxLoc
255 integer(IK) ,
intent(in) ,
optional :: scalingMinLen
256 integer(IK) ,
intent(out) ,
allocatable :: numproc(:)
257 real(RKG) ,
intent(out) ,
allocatable :: scaling(:)
258 real(RKG) ,
intent(in) :: conProb, seqSecTime, parSecTime, comSecTime
259 real(RKG) ,
intent(out) :: scalingMaxVal
264 PURE module subroutine setForkJoinScaling_RK3(conProb, seqSecTime, parSecTime, comSecTime, scaling, numproc, scalingMaxVal, scalingMaxLoc, scalingMinLen)
265#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
266 !DEC$ ATTRIBUTES DLLEXPORT :: setForkJoinScaling_RK3
269 integer(IK) ,
intent(out) :: scalingMaxLoc
270 integer(IK) ,
intent(in) ,
optional :: scalingMinLen
271 integer(IK) ,
intent(out) ,
allocatable :: numproc(:)
272 real(RKG) ,
intent(out) ,
allocatable :: scaling(:)
273 real(RKG) ,
intent(in) :: conProb, seqSecTime, parSecTime, comSecTime
274 real(RKG) ,
intent(out) :: scalingMaxVal
279 PURE module subroutine setForkJoinScaling_RK2(conProb, seqSecTime, parSecTime, comSecTime, scaling, numproc, scalingMaxVal, scalingMaxLoc, scalingMinLen)
280#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
281 !DEC$ ATTRIBUTES DLLEXPORT :: setForkJoinScaling_RK2
284 integer(IK) ,
intent(out) :: scalingMaxLoc
285 integer(IK) ,
intent(in) ,
optional :: scalingMinLen
286 integer(IK) ,
intent(out) ,
allocatable :: numproc(:)
287 real(RKG) ,
intent(out) ,
allocatable :: scaling(:)
288 real(RKG) ,
intent(in) :: conProb, seqSecTime, parSecTime, comSecTime
289 real(RKG) ,
intent(out) :: scalingMaxVal
294 PURE module subroutine setForkJoinScaling_RK1(conProb, seqSecTime, parSecTime, comSecTime, scaling, numproc, scalingMaxVal, scalingMaxLoc, scalingMinLen)
295#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
296 !DEC$ ATTRIBUTES DLLEXPORT :: setForkJoinScaling_RK1
299 integer(IK) ,
intent(out) :: scalingMaxLoc
300 integer(IK) ,
intent(in) ,
optional :: scalingMinLen
301 integer(IK) ,
intent(out) ,
allocatable :: numproc(:)
302 real(RKG) ,
intent(out) ,
allocatable :: scaling(:)
303 real(RKG) ,
intent(in) :: conProb, seqSecTime, parSecTime, comSecTime
304 real(RKG) ,
intent(out) :: scalingMaxVal
316 logical(LK) ,
save :: mv_failed
= .false._LK
369 image
%label
= SK_
"@process("//getStr(image
%id)
//SK_
")"
370 image
%is
%first
= image
%id
== 1_IK
371 image
%is
%extra
= image
%id
/= 1_IK
413#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
422 call mpi_barrier(
mpi_comm_world, ierrMPI)
466#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
475 logical :: isFinalized
476 call mpi_finalized(isFinalized, ierrMPI)
477 if (
.not. isFinalized)
then
478 call mpi_barrier(
mpi_comm_world, ierrMPI)
479 call mpi_finalize(ierrMPI)
524#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
527 integer(IK) :: imageID
529 imageID
= this_image()
533 integer :: rank, ierrMPI
534 logical :: isinit, isfinit
535 call mpi_initialized(isinit, ierrMPI)
536 if (
.not. isinit)
then
537 call mpi_finalized(isfinit, ierrMPI)
538 if (isfinit)
error stop MODULE_NAME//"@getImageID(): Error occurred. A finalized MPI library cannot be reinitialized."
539 call mpi_init(ierrMPI)
541 call mpi_comm_rank(
mpi_comm_world, rank, ierrMPI)
542 if (ierrMPI
/= 0)
error stop "Failed to fetch the MPI process counts."
543 imageID
= int(rank, IK)
+ 1_IK
547 use omp_lib,
only: omp_get_thread_num
548 imageID
= int(omp_get_thread_num()
+ 1, IK)
600#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
603 integer(IK) :: imageCount
605 imageCount
= num_images()
655#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
663 integer(IK) :: imageCount
665 call mpi_initialized(isinit, ierrMPI)
666 if (ierrMPI
/= 0)
return
667 if (
.not. isinit)
then
668 call mpi_init(ierrMPI)
669 if (ierrMPI
/= 0)
return
671 call mpi_comm_size(
mpi_comm_world, nproc, ierrMPI)
672 if (ierrMPI
/= 0)
return
673 imageCount
= int(nproc, IK)
675 integer(IK) :: imageCount
713#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
717 use omp_lib,
only: omp_get_num_threads
718 integer(IK) :: imageCount
719 imageCount
= omp_get_num_threads()
721 integer(IK) :: imageCount
762#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
765 integer(IK),
intent(in),
optional :: count
768 use omp_lib,
only: omp_get_num_procs, omp_set_num_threads
770 count_def
= omp_get_num_procs()
771 if (
present(count))
then
772 if (
0_IK < count) count_def
= int(count)
774 call omp_set_num_threads(count_def)
829#if __INTEL_COMPILER && DLL_ENABLED && (_WIN32 || _WIN64)
832 logical(LK),
intent(in) :: failed
833 logical(LK) :: failedParallelism
835 logical(LK),
allocatable,
save :: failure(:)[:]
836 integer :: iid, imageCount
837 imageCount
= num_images()
838 allocate(failure(imageCount)[
*])
839 failure(
this_image())
= failed
841 do iid
= 1, imageCount
842 failedParallelism
= failure(iid)[iid]
843 if (failedParallelism)
return
846 logical,
allocatable :: failure(:)
847 integer :: imageCount
853 call mpi_comm_size(
mpi_comm_world, imageCount, ierrMPI)
855 call mpi_allgather (
logical(failed)
&
873 failedParallelism
= logical(
any(failure), LK)
877 mv_failed
= mv_failed
.or. failed
880 failedParallelism
= mv_failed
882 mv_failed
= .false._LK
885 failedParallelism
= failed
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Return the predicted parallel Fork-Join speedup scaling behavior for simulations whose image contribu...
Generate and return the conversion of the input value to an output Fortran string,...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
This module contains procedures and generic interfaces for facilitating parallel computations or comp...
integer(IK) integer(IK) function getImageCountOMP()
Generate and return the number of available processes in the current OpenMP-parallel world communicat...
integer(IK) integer(IK) function getImageCountMPI()
Generate and return the number of available processes in the current MPI-parallel world communication...
type(image_type) function image_typer()
Generate and return an object of class image_type containing information and statistics of the parall...
subroutine setImageCount(count)
Set the number the parallel threads for an OpenMP-enabled application.
character(*, SK), parameter MODULE_NAME
logical(LK) function isFailedImage(failed)
Broadcast the error condition from all images/processes to all images/processes.
character(*, SK), parameter PARALLELIZATION_MODE
The scalar constant of type character of default kind SK, whose value is set to parallel if the ParaM...
subroutine setImageSynced()
Synchronize all existing parallel images and return nothing.
integer(IK) function getImageID()
Generate and return the ID of the current Coarray image / MPI process / OpenMP thread,...
integer(IK) function getImageCount()
Generate and return the number of available processes in the current parallel world communication.
subroutine setImageFinalized()
Finalize the current parallel simulation and return nothing.
This module contains the generic procedures for converting values of different types and kinds to For...
This is the image_type type for generating objects that contain information about the current image/p...
This is the imageis_type type for generating objects with components of type logical of default kind ...