MPI related modules¶
parallel.f90
¶
Description
This module contains the blocking information
Quick access
- Variables
chunksize
,ierr
,n_procs
,nr_per_rank
,nthreads
,rank_bn
,rank_with_l1m0
,rank_with_r_lcr
- Routines
check_mpi_error()
,get_openmp_blocks()
,getblocks()
,mpiio_setup()
,parallel()
Needed modules
mpi
omp_lib
Types
-
type
parallel_mod/
unknown_type
¶ - Type fields
%
n_per_rank
[integer ]%
nstart
[integer ]%
nstop
[integer ]
-
type
Variables
-
parallel_mod/
chunksize
[integer]¶
-
-
parallel_mod/
ierr
[integer]¶
-
-
parallel_mod/
n_procs
[integer]¶
-
-
parallel_mod/
nr_per_rank
[integer]¶
-
-
parallel_mod/
nthreads
[integer]¶
-
-
parallel_mod/
rank_bn
[integer]¶
-
-
parallel_mod/
rank_with_l1m0
[integer]¶
-
-
parallel_mod/
rank_with_r_lcr
[integer]¶
-
Subroutines and functions
-
subroutine
parallel_mod/
check_mpi_error
(code)¶ - Parameters
code [integer ,in]
-
subroutine
parallel_mod/
getblocks
(bal, n_points, n_procs)¶ - Parameters
bal (1) [load ,inout]
n_points [integer ,in]
n_procs [integer ,in]
- Called from
initialize_blocking()
,initialize_geos()
,initialize_radial_data()
,readstartfields_mpi()
-
subroutine
parallel_mod/
get_openmp_blocks
(nstart, nstop)¶ - Parameters
nstart [integer ,inout]
nstop [integer ,inout]
- Called from
parallel_solve_phase()
,parallel_solve()
,get_nl_rms()
,set_imex_rhs()
,rotate_imex()
,get_dtblmfinish()
,fft_many()
,ifft_many()
,get_nl()
,lo2r_redist_wait()
,r2lo_redist_start()
,prepare_mat_5()
,transform_to_lm_space()
,get_dr_rloc()
,get_ddr_rloc()
,native_qst_to_spat()
,native_sphtor_to_spat()
,native_sph_to_spat()
,native_sph_to_grad_spat()
,native_spat_to_sph()
,native_spat_to_sph_tor()
,prepareb_fd()
,fill_ghosts_b()
,updateb_fd()
,finish_exp_mag()
,finish_exp_mag_rdist()
,get_mag_ic_rhs_imp()
,assemble_mag_rloc()
,get_mag_rhs_imp()
,get_mag_rhs_imp_ghost()
,preparephase_fd()
,fill_ghosts_phi()
,updatephase_fd()
,get_phase_rhs_imp()
,get_phase_rhs_imp_ghost()
,assemble_phase_rloc()
,prepares_fd()
,fill_ghosts_s()
,updates_fd()
,finish_exp_entropy()
,finish_exp_entropy_rdist()
,get_entropy_rhs_imp()
,get_entropy_rhs_imp_ghost()
,assemble_entropy_rloc()
,preparew_fd()
,fill_ghosts_w()
,updatew_fd()
,finish_exp_pol()
,finish_exp_pol_rdist()
,get_pol_rhs_imp()
,get_pol_rhs_imp_ghost()
,assemble_pol()
,assemble_pol_rloc()
,finish_exp_smat()
,finish_exp_smat_rdist()
,assemble_single()
,get_single_rhs_imp()
,preparexi_fd()
,fill_ghosts_xi()
,updatexi_fd()
,finish_exp_comp()
,finish_exp_comp_rdist()
,get_comp_rhs_imp()
,get_comp_rhs_imp_ghost()
,assemble_comp_rloc()
,preparez_fd()
,fill_ghosts_z()
,updatez_fd()
,get_tor_rhs_imp()
,get_tor_rhs_imp_ghost()
,assemble_tor_rloc()
-
subroutine
parallel_mod/
mpiio_setup
(info)¶ This routine set ups the default MPI-IO configuration. This is based on recommandations from IDRIS “Best practices for parallel IO and MPI-IO hints”
- Parameters
info [integer ,out]
- Called from
write_pot_mpi()
,open_graph_file()
,readstartfields_mpi()
,store_mpi()
radial_data.f90
¶
Description
This module defines the MPI decomposition in the radial direction.
Quick access
- Variables
n_r_cmb
,n_r_icb
,nrstart
,nrstartmag
,nrstop
,nrstopmag
,radial_balance
- Routines
Needed modules
parallel_mod
(rank()
,n_procs()
,nr_per_rank()
,load()
,getblocks()
): This module contains the blocking informationlogic
(l_mag()
,lverbose()
,l_finite_diff()
): Module containing the logicals that control the run
Variables
-
radial_data/
n_r_cmb
[integer,public/protected]¶
-
-
radial_data/
n_r_icb
[integer,public/protected]¶
-
-
radial_data/
nrstart
[integer,public/protected]¶
-
-
radial_data/
nrstartmag
[integer,public/protected]¶
-
-
radial_data/
nrstop
[integer,public/protected]¶
-
-
radial_data/
nrstopmag
[integer,public/protected]¶
-
-
radial_data/
radial_balance
(*) [load,allocatable/public]¶
-
Subroutines and functions
-
subroutine
radial_data/
initialize_radial_data
(n_r_max)¶ This subroutine is used to set up the MPI decomposition in the radial direction
- Parameters
n_r_max [integer ,in] :: Number of radial grid points
- Called from
- Call to
communications.f90
¶
Description
This module contains the different MPI communicators used in MagIC.
Quick access
- Variables
get_global_sum
,gt_cheb
,gt_ic
,gt_oc
,reduce_radial
,send_lm_pair_to_master
,temp_gather_lo
- Routines
allgather_from_rloc()
,create_gather_type()
,destroy_gather_type()
,finalize_communications()
,find_faster_block()
,find_faster_comm()
,gather_all_from_lo_to_rank0()
,gather_from_lo_to_rank0()
,gather_from_rloc()
,get_global_sum_cmplx_1d()
,get_global_sum_cmplx_2d()
,get_global_sum_real_2d()
,initialize_communications()
,lm2lo_redist()
,lo2lm_redist()
,myallgather()
,reduce_radial_1d()
,reduce_radial_2d()
,reduce_scalar()
,scatter_from_rank0_to_lo()
,send_lm_pair_to_master_arr()
,send_lm_pair_to_master_scal_cmplx()
,send_lm_pair_to_master_scal_real()
,send_scal_lm_to_master()
Needed modules
mpimod
constants
(zero()
): module containing constants and parameters used in the code.precision_mod
: This module controls the precision used in MagICmem_alloc
(memwrite()
,bytes_allocated()
): This little module is used to estimate the global memory allocation used in MagICparallel_mod
(rank()
,n_procs()
,ierr()
): This module contains the blocking informationtruncation
(l_max()
,lm_max()
,minc()
,n_r_max()
,n_r_ic_max()
,fd_order()
,fd_order_bound()
,m_max()
,m_min()
): This module defines the grid points and the truncationblocking
(st_map()
,lo_map()
,lm_balance()
,llm()
,ulm()
): Module containing blocking informationradial_data
(nrstart()
,nrstop()
,radial_balance()
): This module defines the MPI decomposition in the radial direction.logic
(l_mag()
,l_conv()
,l_heat()
,l_chemical_conv()
,l_finite_diff()
,l_mag_kin()
,l_double_curl()
,l_save_out()
,l_packed_transp()
,l_parallel_solve()
,l_mag_par_solve()
,l_phase_field()
): Module containing the logicals that control the runuseful
(abortrun()
): This module contains several useful routines.output_data
(n_log_file()
,log_file()
): This module contains the parameters for output controliso_fortran_env
(output_unit()
)mpi_ptop_mod
(type_mpiptop()
): This module contains the implementation of MPI_Isend/MPI_Irecv global transposempi_alltoall_mod
(type_mpiatoav()
,type_mpiatoaw()
,type_mpiatoap()
): This module contains the implementation of all-to-all global communicatorscharmanip
(capitalize()
): This module contains several useful routines to manipule character stringsnum_param
(mpi_transp()
,mpi_packing()
): Module containing numerical and control parametersmpi_transp_mod
(type_mpitransp()
): This is an abstract class that will be used to define MPI transposers The actual implementation is deferred to either point-to-point (MPI_Isend and MPI_IRecv) communications or all-to-all (MPI_AlltoAll)
Types
-
type
communications/
unknown_type
¶ - Type fields
%
dim2
[integer ]%
gather_mpi_type
(*) [integer ,allocatable]
-
type
Variables
-
communications/
get_global_sum
[public]¶
-
-
communications/
gt_cheb
[gather_type,public]¶
-
-
communications/
gt_ic
[gather_type,public]¶
-
-
communications/
gt_oc
[gather_type,public]¶
-
-
communications/
reduce_radial
[public]¶
-
-
communications/
send_lm_pair_to_master
[public]¶
-
-
communications/
temp_gather_lo
(*) [complex,private/allocatable]¶
-
Subroutines and functions
-
subroutine
communications/
initialize_communications
()¶ - Called from
- Call to
create_gather_type()
,capitalize()
,find_faster_comm()
,memwrite()
-
subroutine
communications/
finalize_communications
()¶ - Called from
- Call to
-
function
communications/
get_global_sum_cmplx_2d
(dwdt_local)¶ - Parameters
dwdt_local (*,*) [complex ,in]
- Return
global_sum [real ]
-
function
communications/
get_global_sum_real_2d
(dwdt_local)¶ - Parameters
dwdt_local (*,*) [real ,in]
- Return
global_sum [real ]
-
function
communications/
get_global_sum_cmplx_1d
(arr_local)¶ Kahan summation algorithm
function KahanSum(input) var sum = 0.0 var c = 0.0 //A running compensation for lost low-order bits. for i = 1 to input.length do y = input[i] - c //So far, so good: c is zero. t = sum + y //Alas, sum is big, y small, //so low-order digits of y are lost. c = (t - sum) - y //(t - sum) recovers the high-order part of y; //subtracting y recovers -(low part of y) sum = t //Algebraically, c should always be zero. //Beware eagerly optimising compilers! //Next time around, the lost low part will be added to y in a fresh attempt. return sum
- Parameters
arr_local (*) [complex ,in] :: So far, so good: c is zero.
- Return
global_sum [real ]
-
subroutine
communications/
gather_all_from_lo_to_rank0
(self, arr_lo, arr_full)¶ - Parameters
- Called from
dtb_gather_lo_on_rank0()
,fields_average()
,write_pot_mpi()
,write_pot()
,write_one_field()
,write_movie_frame()
,output()
,store()
-
subroutine
communications/
create_gather_type
(self, dim2)¶ Define the datatypes for gather_all_from_lo_to_rank0 the sending array has dimension (llm:ulm,1:dim2) receiving array has dimension (1:lm_max,1:dim2)
- Parameters
self [gather_type ]
dim2 [integer ]
- Called from
-
subroutine
communications/
destroy_gather_type
(self)¶ - Parameters
self [gather_type ]
- Called from
-
subroutine
communications/
gather_from_lo_to_rank0
(arr_lo, arr_full)¶ - Parameters
- Called from
fields_average()
,get_e_mag()
,get_onset()
,write_bcmb()
,write_coeff_r()
,output()
,store_mpi()
-
subroutine
communications/
scatter_from_rank0_to_lo
(arr_full, arr_lo)¶ - Parameters
- Called from
readstartfields_old()
,readstartfields()
,read_map_one_field()
,readstartfields_mpi()
,step_time()
-
subroutine
communications/
send_lm_pair_to_master_arr
(b, l, m, blm0)¶
-
subroutine
communications/
send_lm_pair_to_master_scal_real
(b, l, m, blm0)¶
-
subroutine
communications/
send_lm_pair_to_master_scal_cmplx
(b, l, m, blm0)¶
-
subroutine
communications/
send_scal_lm_to_master
(blm0, l, m)¶ - Parameters
blm0 [real ,inout]
l [integer ,in]
m [integer ,in]
-
subroutine
communications/
lm2lo_redist
(arr_lmloc, arr_lo)¶
-
subroutine
communications/
lo2lm_redist
(arr_lo, arr_lmloc)¶
-
subroutine
communications/
gather_from_rloc
(arr_rloc, arr_glob, irank)¶ This subroutine gather a r-distributed array on rank=irank
- Parameters
- Called from
outhemi()
,outhelicity()
,outphase()
,outpar()
,outperppar()
,get_power()
-
subroutine
communications/
allgather_from_rloc
(arr_rloc, arr_glob)¶ This subroutine gather a r-distributed array
- Parameters
- Called from
-
subroutine
communications/
reduce_radial_2d
(arr_dist, arr_glob, irank)¶ - Parameters
arr_dist (*,*) [real ,in]
arr_glob (*,*) [real ,out]
irank [integer ,in]
-
subroutine
communications/
reduce_radial_1d
(arr_dist, arr_glob, irank)¶ - Parameters
arr_dist (*) [real ,in]
arr_glob (*) [real ,inout]
irank [integer ,in]
-
subroutine
communications/
reduce_scalar
(scal_dist, scal_glob, irank)¶ - Parameters
scal_dist [real ,in]
scal_glob [real ,inout]
irank [integer ,in]
- Called from
-
subroutine
communications/
find_faster_block
(idx_type)¶ - Parameters
idx_type [integer ,in]
- Called from
-
subroutine
communications/
find_faster_comm
(idx, mintime, n_fields)¶ This subroutine tests the two MPI transposition strategies and selects the fastest one.
- Parameters
idx [integer ,out]
mintime [real ,out]
n_fields [integer ,in]
- Called from
-
subroutine
communications/
myallgather
(arr, dim1, dim2)¶ - Parameters
arr (dim1,dim2) [complex ,inout]
dim1 [integer ,in,]
dim2 [integer ,in,]
- Use
mpi_transpose.f90
¶
Description
This is an abstract class that will be used to define MPI transposers The actual implementation is deferred to either point-to-point (MPI_Isend and MPI_IRecv) communications or all-to-all (MPI_AlltoAll)
Needed modules
precision_mod
: This module controls the precision used in MagICtruncation
(lm_max()
,n_r_max()
): This module defines the grid points and the truncationradial_data
(nrstart()
,nrstop()
): This module defines the MPI decomposition in the radial direction.blocking
(llm()
,ulm()
): Module containing blocking information
Types
-
type
mpi_transp_mod/
unknown_type
¶ - Type fields
%
n_fields
[integer ]
-
type
parallel_solvers.f90
¶
Description
This module contains the routines that are used to solve linear banded problems with R-distributed arrays.
Quick access
- Routines
finalize_3()
,finalize_5()
,initialize_3()
,initialize_5()
,prepare_mat_3()
,prepare_mat_5()
,solver_dn_3()
,solver_dn_5()
,solver_finish_3()
,solver_finish_5()
,solver_single()
,solver_up_3()
,solver_up_5()
Needed modules
precision_mod
: This module controls the precision used in MagICparallel_mod
: This module contains the blocking informationradial_data
(n_r_cmb()
,n_r_icb()
): This module defines the MPI decomposition in the radial direction.mem_alloc
(bytes_allocated()
): This little module is used to estimate the global memory allocation used in MagICconstants
(one()
): module containing constants and parameters used in the code.truncation
(lm_max()
): This module defines the grid points and the truncation
Types
-
type
parallel_solvers/
unknown_type
¶ - Type fields
%
diag
(*,*) [real ,allocatable]%
lmax
[integer ]%
lmin
[integer ]%
low
(*,*) [real ,allocatable]%
nrmax
[integer ]%
nrmin
[integer ]%
up
(*,*) [real ,allocatable]
-
type
-
type
parallel_solvers/
unknown_type
- Type fields
%
diag
(*,*) [real ,allocatable]%
lmax
[integer ]%
lmin
[integer ]%
low1
(*,*) [real ,allocatable]%
low2
(*,*) [real ,allocatable]%
nrmax
[integer ]%
nrmin
[integer ]%
up1
(*,*) [real ,allocatable]%
up2
(*,*) [real ,allocatable]
-
type
Subroutines and functions
-
subroutine
parallel_solvers/
initialize_3
(this, nrstart, nrstop, lmin, lmax)¶ Memory allocation of a parallel tridiagonal matrix
- Parameters
this [real ]
nrstart [integer ,in]
nrstop [integer ,in]
lmin [integer ,in]
lmax [integer ,in]
-
subroutine
parallel_solvers/
initialize_5
(this, nrstart, nrstop, lmin, lmax)¶ Memory allocation of a parallel tridiagonal matrix
- Parameters
this [real ]
nrstart [integer ,in]
nrstop [integer ,in]
lmin [integer ,in]
lmax [integer ,in]
-
subroutine
parallel_solvers/
finalize_3
(this)¶ Memory deallocation of a parallel tridiagonal solver
- Parameters
this [real ]
-
subroutine
parallel_solvers/
finalize_5
(this)¶ Memory deallocation of a parallel pentadiagonal solver
- Parameters
this [real ]
-
subroutine
parallel_solvers/
prepare_mat_3
(this)¶ LU factorisation of a tridiagonal matrix: the diagonal is overwritten
- Parameters
this [real ]
-
subroutine
parallel_solvers/
prepare_mat_5
(this)¶ LU factorisation of a pentadiagonal matrix
- Parameters
this [real ]
- Call to
-
subroutine
parallel_solvers/
solver_single
(this, x, nrstart, nrstop)¶ This routine is used to solve a solve one single linear system that does not depend on lm. This is for instance used for z10 when l_z10mat is required.
-
subroutine
parallel_solvers/
solver_up_3
(this, x, lmstart, lmstop, nrstart, nrstop, tag, array_req, req, lms_block, nlm_block)¶ First part of the parallel tridiag solver: forward substitution
- Parameters
this [real ]
lmstart [integer ,in] :: Starting lm (OMP thread dependent)
lmstop [integer ,in] :: Stopping lm (OMP thread dependent)
nrstart [integer ,in] :: Starting nR
nrstop [integer ,in] :: Stopping nR
tag [integer ,in]
array_req (*) [integer ,inout]
req [integer ,inout]
lms_block [integer ,in] :: Starting block-index of lm
nlm_block [integer ,in] :: Size of the block
-
subroutine
parallel_solvers/
solver_up_5
(this, x, lmstart, lmstop, nrstart, nrstop, tag, array_req, req, lms_block, nlm_block)¶ First part of the parallel pentadiag solver: forward substitution
- Parameters
this [real ]
lmstart [integer ,in] :: Starting lm (OMP thread dependent)
lmstop [integer ,in] :: Stopping lm (OMP thread dependent)
nrstart [integer ,in] :: Starting nR
nrstop [integer ,in] :: Stopping nR
tag [integer ,in]
array_req (*) [integer ,inout]
req [integer ,inout]
lms_block [integer ,in] :: Starting block-index of lm
nlm_block [integer ,in] :: Size of the block
-
subroutine
parallel_solvers/
solver_dn_3
(this, x, lmstart, lmstop, nrstart, nrstop, tag, array_req, req, lms_block, nlm_block)¶ Second part of the parallel tridiag solver: backward substitution
- Parameters
this [real ]
lmstart [integer ,in] :: Starting lm (OMP thread dependent)
lmstop [integer ,in] :: Stopping lm (OMP thread dependent)
nrstart [integer ,in] :: Starting nR
nrstop [integer ,in] :: Stopping nR
tag [integer ,in]
array_req (*) [integer ,inout]
req [integer ,inout]
lms_block [integer ,in] :: Starting block-index of lm
nlm_block [integer ,in] :: Size of the block
-
subroutine
parallel_solvers/
solver_dn_5
(this, x, lmstart, lmstop, nrstart, nrstop, tag, array_req, req, lms_block, nlm_block)¶ Second part of the parallel pentadiagonal solver: backward substitution
- Parameters
this [real ]
lmstart [integer ,in] :: Starting lm (OMP thread dependent)
lmstop [integer ,in] :: Stopping lm (OMP thread dependent)
nrstart [integer ,in] :: Starting nR
nrstop [integer ,in] :: Stopping nR
tag [integer ,in]
array_req (*) [integer ,inout]
req [integer ,inout]
lms_block [integer ,in] :: Starting block-index of lm
nlm_block [integer ,in] :: Size of the block
-
subroutine
parallel_solvers/
solver_finish_3
(this, x, lms_block, nlm_block, nrstart, nrstop, tag, array_req, req)¶ Last part of the parallel tridiag solver: halo synchronisation
- Parameters
this [real ]
lms_block [integer ,in] :: Starting block-index of lm
nlm_block [integer ,in] :: Size of the block
nrstart [integer ,in] :: Starting index in radius
nrstop [integer ,in] :: Stopping index in radius
tag [integer ,in]
array_req (*) [integer ,inout]
req [integer ,inout]
-
subroutine
parallel_solvers/
solver_finish_5
(this, x, lms_block, nlm_block, nrstart, nrstop, tag, array_req, req)¶ Last part of the parallel pentadiag solver: halo synchronisation
- Parameters
this [real ]
lms_block [integer ,in] :: Starting block-index of lm
nlm_block [integer ,in] :: Size of the block
nrstart [integer ,in] :: Starting index in radius
nrstop [integer ,in] :: Stopping index in radius
tag [integer ,in]
array_req (*) [integer ,inout] :: MPI requests
req [integer ,inout]