SCALE-RM
Data Types | Functions/Subroutines | Variables
scale_prc Module Reference

module PROCESS More...

Functions/Subroutines

subroutine, public prc_mpistart (comm)
 Start MPI. More...
 
subroutine, public prc_universal_setup (comm, nprocs, myrank, ismaster)
 setup MPI in universal communicator More...
 
subroutine, public prc_global_setup (abortall, comm)
 setup MPI in global communicator More...
 
subroutine, public prc_local_setup (comm, myrank, ismaster)
 Setup MPI in local communicator. More...
 
subroutine, public prc_singlecom_setup (comm, nprocs, myrank, ismaster)
 Setup MPI single communicator (not use universal-global-local setting) More...
 
subroutine, public prc_errhandler_setup (use_fpm, master)
 Setup MPI error handler. More...
 
subroutine, public prc_abort
 Abort Process. More...
 
subroutine, public prc_mpifinish
 Stop MPI peacefully. More...
 
subroutine, public prc_mpisplit_bulk (ORG_COMM_WORLD, NUM_BULKJOB, PRC_BULKJOB, debug, SUB_COMM_WORLD, ID_BULKJOB)
 MPI Communicator Split (bulk job) More...
 
subroutine, public prc_mpisplit_nest (ORG_COMM_WORLD, NUM_DOMAIN, PRC_DOMAIN, debug, color_reorder, SUB_COMM_WORLD, ID_DOMAIN, INTERCOMM_parent, INTERCOMM_child)
 MPI Communicator Split (nesting) More...
 
subroutine, public prc_mpibarrier
 Barrier MPI. More...
 
real(dp) function, public prc_mpitime ()
 Get MPI time. More...
 
subroutine, public prc_mpitimestat (avgvar, maxvar, minvar, maxidx, minidx, var)
 Calc global statistics for timer. More...
 
subroutine, public prc_set_file_closer (routine)
 

Variables

integer, parameter, public prc_masterrank = 0
 master process in each communicator More...
 
integer, parameter, public prc_domain_nlim = 10000
 max depth of domains More...
 
integer, parameter, public prc_comm_null = MPI_COMM_NULL
 
integer, public prc_universal_comm_world = -1
 original communicator More...
 
integer, public prc_universal_myrank = -1
 myrank in universal communicator More...
 
integer, public prc_universal_nprocs = -1
 process num in universal communicator More...
 
logical, public prc_universal_ismaster = .false.
 master process in universal communicator? More...
 
integer, public prc_universal_jobid = 0
 my job ID in universal communicator More...
 
integer, public prc_global_comm_world = -1
 global communicator More...
 
integer, public prc_global_myrank = -1
 myrank in global communicator More...
 
integer, public prc_global_nprocs = -1
 process num in global communicator More...
 
logical, public prc_global_ismaster = .false.
 master process in global communicator? More...
 
integer, public prc_global_domainid = 0
 my domain ID in global communicator More...
 
integer, dimension(prc_domain_nlim), public prc_global_root
 root processes in global members More...
 
integer, public prc_local_comm_world = -1
 local communicator More...
 
integer, public prc_nprocs = 1
 myrank in local communicator More...
 
integer, public prc_myrank = 0
 process num in local communicator More...
 
logical, public prc_ismaster = .false.
 master process in local communicator? More...
 
logical, public prc_mpi_alive = .false.
 MPI is alive? More...
 
integer, public prc_universal_handler
 error handler in universal communicator More...
 
integer, public prc_abort_comm_world
 communicator for aborting More...
 
integer, public prc_abort_handler
 error handler communicator for aborting More...
 

Detailed Description

module PROCESS

Description
MPI/non-MPI management module
Author
Team SCALE

Function/Subroutine Documentation

◆ prc_mpistart()

subroutine, public scale_prc::prc_mpistart ( integer, intent(out)  comm)

Start MPI.

Definition at line 122 of file scale_prc.F90.

122  implicit none
123 
124  integer, intent(out) :: comm ! communicator
125 
126  integer :: ierr
127  !---------------------------------------------------------------------------
128 
129  call mpi_init(ierr)
130 
131  prc_mpi_alive = .true.
132 ! PRC_UNIVERSAL_handler = MPI_ERRHANDLER_NULL
133 ! call MPI_COMM_CREATE_ERRHANDLER( PRC_MPI_errorhandler, PRC_UNIVERSAL_handler, ierr )
134 
135  comm = mpi_comm_world
136  prc_abort_comm_world = comm
137 
138  return

References prc_abort_comm_world, and prc_mpi_alive.

Referenced by scale::scale_init(), scalerm(), scalerm_init(), and scalerm_pp().

Here is the caller graph for this function:

◆ prc_universal_setup()

subroutine, public scale_prc::prc_universal_setup ( integer, intent(in)  comm,
integer, intent(out)  nprocs,
integer, intent(out)  myrank,
logical, intent(out)  ismaster 
)

setup MPI in universal communicator

Definition at line 148 of file scale_prc.F90.

148  implicit none
149 
150  integer, intent(in) :: comm ! communicator
151  integer, intent(out) :: nprocs ! number of procs in this communicator
152  integer, intent(out) :: myrank ! myrank in this communicator
153  logical, intent(out) :: ismaster ! master process in this communicator?
154 
155  integer :: ierr
156  !---------------------------------------------------------------------------
157 
158  prc_universal_comm_world = comm
159 
160  call mpi_comm_size(prc_universal_comm_world,prc_universal_nprocs,ierr)
161  call mpi_comm_rank(prc_universal_comm_world,prc_universal_myrank,ierr)
162 
163  if ( prc_universal_myrank == prc_masterrank ) then
164  prc_universal_ismaster = .true.
165  else
166  prc_universal_ismaster = .false.
167  endif
168 
169  nprocs = prc_universal_nprocs
170  myrank = prc_universal_myrank
171  ismaster = prc_universal_ismaster
172 
173 
174 
175 ! PRC_ABORT_COMM_WORLD = PRC_UNIVERSAL_COMM_WORLD
176 !
177 ! call MPI_Comm_set_errhandler(PRC_ABORT_COMM_WORLD,PRC_UNIVERSAL_handler,ierr)
178 ! call MPI_Comm_get_errhandler(PRC_ABORT_COMM_WORLD,PRC_ABORT_handler ,ierr)
179 
180  return

References prc_masterrank, prc_universal_comm_world, prc_universal_ismaster, prc_universal_myrank, and prc_universal_nprocs.

Referenced by scalerm(), scalerm_init(), and scalerm_pp().

Here is the caller graph for this function:

◆ prc_global_setup()

subroutine, public scale_prc::prc_global_setup ( logical, intent(in)  abortall,
integer, intent(in)  comm 
)

setup MPI in global communicator

Definition at line 188 of file scale_prc.F90.

188  implicit none
189 
190  logical, intent(in) :: abortall ! abort all jobs?
191  integer, intent(in) :: comm ! communicator
192 
193  integer :: ierr
194  !---------------------------------------------------------------------------
195 
196  prc_global_comm_world = comm
197 
198  call mpi_comm_size(prc_global_comm_world,prc_global_nprocs,ierr)
199  call mpi_comm_rank(prc_global_comm_world,prc_global_myrank,ierr)
200 
201  if ( prc_global_myrank == prc_masterrank ) then
202  prc_global_ismaster = .true.
203  else
204  prc_global_ismaster = .false.
205  endif
206 
207 ! if ( .NOT. abortall ) then
208 ! PRC_ABORT_COMM_WORLD = PRC_GLOBAL_COMM_WORLD
209 !
210 ! call MPI_COMM_SET_ERRHANDLER(PRC_ABORT_COMM_WORLD,PRC_UNIVERSAL_handler,ierr)
211 ! call MPI_COMM_GET_ERRHANDLER(PRC_ABORT_COMM_WORLD,PRC_ABORT_handler ,ierr)
212 ! endif
213 
214  return

References prc_global_comm_world, prc_global_ismaster, prc_global_myrank, prc_global_nprocs, and prc_masterrank.

Referenced by scalerm(), scalerm_init(), and scalerm_pp().

Here is the caller graph for this function:

◆ prc_local_setup()

subroutine, public scale_prc::prc_local_setup ( integer, intent(in)  comm,
integer, intent(out)  myrank,
logical, intent(out)  ismaster 
)

Setup MPI in local communicator.

Definition at line 223 of file scale_prc.F90.

223  implicit none
224 
225  integer, intent(in) :: comm ! communicator
226  integer, intent(out) :: myrank ! myrank in this communicator
227  logical, intent(out) :: ismaster ! master process in this communicator?
228 
229  integer :: ierr
230  !---------------------------------------------------------------------------
231 
232  prc_local_comm_world = comm
233 
234  call mpi_comm_rank(prc_local_comm_world,prc_myrank,ierr)
235  call mpi_comm_size(prc_local_comm_world,prc_nprocs,ierr)
236 
237  if ( prc_myrank == prc_masterrank ) then
238  prc_ismaster = .true.
239  else
240  prc_ismaster = .false.
241  endif
242 
243  myrank = prc_myrank
244  ismaster = prc_ismaster
245 
246  return

References prc_ismaster, prc_local_comm_world, prc_masterrank, prc_myrank, and prc_nprocs.

Referenced by mod_rm_driver::rm_driver(), and mod_rm_prep::rm_prep().

Here is the caller graph for this function:

◆ prc_singlecom_setup()

subroutine, public scale_prc::prc_singlecom_setup ( integer, intent(in)  comm,
integer, intent(out)  nprocs,
integer, intent(out)  myrank,
logical, intent(out)  ismaster 
)

Setup MPI single communicator (not use universal-global-local setting)

Definition at line 256 of file scale_prc.F90.

256  implicit none
257 
258  integer, intent(in) :: comm ! communicator
259  integer, intent(out) :: nprocs ! number of procs
260  integer, intent(out) :: myrank ! myrank
261  logical, intent(out) :: ismaster ! master process?
262 
263  integer :: ierr
264  !---------------------------------------------------------------------------
265 
266  call mpi_comm_size(comm,nprocs,ierr)
267  call mpi_comm_rank(comm,myrank,ierr)
268 
269  if ( myrank == prc_masterrank ) then
270  ismaster = .true.
271  else
272  ismaster = .false.
273  endif
274 
275  prc_universal_comm_world = comm
276  prc_universal_nprocs = nprocs
277  prc_universal_myrank = myrank
278  prc_universal_ismaster = ismaster
279 
280  prc_global_comm_world = comm
281  prc_global_nprocs = nprocs
282  prc_global_myrank = myrank
283  prc_global_ismaster = ismaster
284 
285  prc_local_comm_world = comm
286  prc_nprocs = nprocs
287  prc_myrank = myrank
288  prc_ismaster = ismaster
289 
290 
291 
292  prc_abort_comm_world = comm
293 
294 ! call MPI_Comm_set_errhandler(PRC_ABORT_COMM_WORLD,PRC_UNIVERSAL_handler,ierr)
295 ! call MPI_Comm_get_errhandler(PRC_ABORT_COMM_WORLD,PRC_ABORT_handler ,ierr)
296 
297  return

References prc_abort_comm_world, prc_global_comm_world, prc_global_ismaster, prc_global_myrank, prc_global_nprocs, prc_ismaster, prc_local_comm_world, prc_masterrank, prc_myrank, prc_nprocs, prc_universal_comm_world, prc_universal_ismaster, prc_universal_myrank, and prc_universal_nprocs.

Referenced by scale::scale_init().

Here is the caller graph for this function:

◆ prc_errhandler_setup()

subroutine, public scale_prc::prc_errhandler_setup ( logical, intent(in)  use_fpm,
logical, intent(in)  master 
)

Setup MPI error handler.

Definition at line 305 of file scale_prc.F90.

305  implicit none
306 
307  logical, intent(in) :: use_fpm ! fpm switch
308  logical, intent(in) :: master ! master flag
309 
310  integer :: ierr
311  !---------------------------------------------------------------------------
312 
313  call mpi_comm_create_errhandler(prc_mpi_errorhandler,prc_universal_handler,ierr)
314 
315  call mpi_comm_set_errhandler(prc_abort_comm_world,prc_universal_handler,ierr)
316  call mpi_comm_get_errhandler(prc_abort_comm_world,prc_abort_handler ,ierr)
317 
318  if ( prc_universal_handler /= prc_abort_handler ) then
319  if( prc_universal_ismaster ) write(*,*) ""
320  if( prc_universal_ismaster ) write(*,*) "ERROR: MPI HANDLER is INCONSISTENT"
321  if( prc_universal_ismaster ) write(*,*) " PRC_UNIVERSAL_handler = ", prc_universal_handler
322  if( prc_universal_ismaster ) write(*,*) " PRC_ABORT_handler = ", prc_abort_handler
323  call prc_abort
324  endif
325 
326  if ( use_fpm ) then
327  call sigvars_get_all( master )
328  call signal( sigint, prc_abort )
329  call signal( sigquit, prc_abort )
330  call signal( sigabrt, prc_abort )
331  call signal( sigfpe, prc_abort )
332  call signal( sigsegv, prc_abort )
333  call signal( sigterm, prc_abort )
334  endif
335 
336  return

References prc_abort(), prc_abort_comm_world, prc_abort_handler, prc_universal_handler, prc_universal_ismaster, scale_sigvars::sigabrt, scale_sigvars::sigfpe, scale_sigvars::sigint, scale_sigvars::sigquit, scale_sigvars::sigsegv, scale_sigvars::sigterm, and scale_sigvars::sigvars_get_all().

Referenced by scale::scale_init(), scalerm(), scalerm_init(), and scalerm_pp().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ prc_abort()

subroutine, public scale_prc::prc_abort

Abort Process.

Definition at line 342 of file scale_prc.F90.

342  implicit none
343 
344  integer :: ierr
345  !---------------------------------------------------------------------------
346 
347  if ( prc_mpi_alive ) then
348  ! tentative approach; input "PRC_UNIVERSAL_COMM_WORLD".
349  call mpi_comm_call_errhandler(prc_universal_comm_world,prc_abort_code,ierr)
350  endif
351 
352  stop

References prc_mpi_alive, and prc_universal_comm_world.

Referenced by mod_admin_restart::admin_restart_setup(), mod_admin_time::admin_time_setup(), mod_admin_versioncheck::admin_versioncheck(), scale_atmos_adiabat::atmos_adiabat_cape_3d(), scale_atmos_adiabat::atmos_adiabat_liftparcel_3d(), mod_atmos_admin::atmos_admin_getscheme(), mod_atmos_admin::atmos_admin_setup(), mod_atmos_bnd_driver::atmos_boundary_driver_send(), mod_atmos_bnd_driver::atmos_boundary_driver_set(), mod_atmos_bnd_driver::atmos_boundary_driver_setup(), mod_atmos_bnd_driver::atmos_boundary_driver_update(), mod_atmos_bnd_driver::atmos_boundary_set_file(), mod_atmos_bnd_driver::atmos_boundary_set_online(), mod_atmos_dyn_driver::atmos_dyn_driver_setup(), scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_apply_numfilter(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_setup(), scale_atmos_dyn_fvm_numfilter::atmos_dyn_fvm_numfilter_setup(), scale_atmos_dyn::atmos_dyn_setup(), scale_atmos_dyn_tinteg_large_euler::atmos_dyn_tinteg_large_euler_setup(), scale_atmos_dyn_tinteg_large_rk3::atmos_dyn_tinteg_large_rk3_setup(), scale_atmos_dyn_tinteg_large::atmos_dyn_tinteg_large_setup(), scale_atmos_dyn_tinteg_short_rk11s8o::atmos_dyn_tinteg_short_rk11s8o_setup(), scale_atmos_dyn_tinteg_short_rk3::atmos_dyn_tinteg_short_rk3_setup(), scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4_setup(), scale_atmos_dyn_tinteg_short_rk7s6o::atmos_dyn_tinteg_short_rk7s6o_setup(), scale_atmos_dyn_tinteg_short::atmos_dyn_tinteg_short_setup(), scale_atmos_dyn_tinteg_tracer_euler::atmos_dyn_tinteg_tracer_euler_setup(), scale_atmos_dyn_tinteg_tracer_linrk::atmos_dyn_tinteg_tracer_linrk_setup(), scale_atmos_dyn_tinteg_tracer_rk3::atmos_dyn_tinteg_tracer_rk3_setup(), scale_atmos_dyn_tinteg_tracer::atmos_dyn_tinteg_tracer_setup(), scale_atmos_dyn_tstep_large_fvm_heve::atmos_dyn_tstep_large_fvm_heve_setup(), scale_atmos_dyn_tstep_large::atmos_dyn_tstep_large_setup(), scale_atmos_dyn_tstep_short_fvm_heve::atmos_dyn_tstep_short_fvm_heve_regist(), scale_atmos_dyn_tstep_short_fvm_hevi::atmos_dyn_tstep_short_fvm_hevi_regist(), scale_atmos_dyn_tstep_short_fvm_hivi::atmos_dyn_tstep_short_fvm_hivi_regist(), scale_atmos_dyn_tstep_short_fvm_hivi::atmos_dyn_tstep_short_fvm_hivi_setup(), scale_atmos_dyn_tstep_short::atmos_dyn_tstep_short_regist(), scale_atmos_dyn_tstep_tracer_fvm_heve::atmos_dyn_tstep_tracer_fvm_heve_setup(), scale_atmos_dyn_tstep_tracer::atmos_dyn_tstep_tracer_setup(), mod_atmos_dyn_vars::atmos_dyn_vars_restart_read(), mod_atmos_dyn_vars::atmos_dyn_vars_setup(), scale_atmos_grid_cartesc::atmos_grid_cartesc_generate(), scale_atmos_grid_cartesc_index::atmos_grid_cartesc_index_setup_main(), scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_rotcoef(), scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_setup(), scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_setup(), scale_atmos_grid_cartesc::atmos_grid_cartesc_setup(), scale_atmos_grid_icoa_index::atmos_grid_icoa_index_setup(), scale_atmos_hydrometeor::atmos_hydrometeor_regist(), scale_atmos_hydrometeor::atmos_hydrometeor_setup(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_0d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_3d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_rev_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_setup(), mod_atmos_phy_ae_driver::atmos_phy_ae_driver_calc_tendency(), mod_atmos_phy_ae_driver::atmos_phy_ae_driver_setup(), mod_atmos_phy_ae_driver::atmos_phy_ae_driver_tracer_setup(), scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_setup(), scale_atmos_phy_ae_kajino13::atmos_phy_ae_kajino13_tracer_setup(), scale_atmos_phy_ae_offline::atmos_phy_ae_offline_qtrc2qaero(), scale_atmos_phy_ae_offline::atmos_phy_ae_offline_setup(), scale_atmos_phy_ae_offline::atmos_phy_ae_offline_tendency(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_setup(), mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_setup(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tendency(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_tracer_setup(), mod_atmos_phy_bl_vars::atmos_phy_bl_vars_setup(), mod_atmos_phy_ch_driver::atmos_phy_ch_driver_tracer_setup(), scale_atmos_phy_ch_rn222::atmos_phy_ch_rn222_setup(), mod_atmos_phy_ch_vars::atmos_phy_ch_vars_setup(), scale_atmos_phy_cp_common::atmos_phy_cp_common_setup(), mod_atmos_phy_cp_driver::atmos_phy_cp_driver_setup(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_setup(), scale_atmos_phy_cp_kf::atmos_phy_cp_kf_tendency(), mod_atmos_phy_cp_vars::atmos_phy_cp_vars_setup(), mod_atmos_phy_lt_driver::atmos_phy_lt_driver_setup(), mod_atmos_phy_lt_driver::atmos_phy_lt_driver_tracer_setup(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_electric_field(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_f2013(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_mg2001(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_adjustment(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_setup(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_solve_bicgstab(), mod_atmos_phy_lt_vars::atmos_phy_lt_vars_setup(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_qhyd2qtrc(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_setup(), mod_atmos_phy_mp_driver::atmos_phy_mp_driver_tracer_setup(), scale_atmos_phy_mp_kessler::atmos_phy_mp_kessler_setup(), scale_atmos_phy_mp_common::atmos_phy_mp_negative_fixer(), scale_atmos_phy_mp_common::atmos_phy_mp_saturation_adjustment_3d(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_setup(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_terminal_velocity(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_crg_qtrc2qhyd(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_setup(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_tracer_setup(), scale_atmos_phy_mp_tomita08::atmos_phy_mp_tomita08_setup(), mod_atmos_phy_mp_vars::atmos_phy_mp_vars_setup(), mod_atmos_phy_rd_driver::atmos_phy_rd_driver_setup(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_flux(), scale_atmos_phy_rd_mstrnx::atmos_phy_rd_mstrnx_setup(), scale_atmos_phy_rd_offline::atmos_phy_rd_offline_flux(), scale_atmos_phy_rd_offline::atmos_phy_rd_offline_setup(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup(), scale_atmos_phy_rd_profile::atmos_phy_rd_profile_setup_zgrid(), mod_atmos_phy_rd_vars::atmos_phy_rd_vars_setup(), scale_atmos_phy_sf_bulk::atmos_phy_sf_bulk_setup(), scale_atmos_phy_sf_const::atmos_phy_sf_const_setup(), mod_atmos_phy_sf_driver::atmos_phy_sf_driver_setup(), mod_atmos_phy_sf_vars::atmos_phy_sf_vars_setup(), scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980_config(), scale_atmos_phy_tb_d1980::atmos_phy_tb_d1980_setup(), scale_atmos_phy_tb_dns::atmos_phy_tb_dns_config(), scale_atmos_phy_tb_dns::atmos_phy_tb_dns_setup(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver_tracer_setup(), scale_atmos_phy_tb_smg::atmos_phy_tb_smg_setup(), mod_atmos_phy_tb_vars::atmos_phy_tb_vars_setup(), scale_atmos_refstate::atmos_refstate_read(), scale_atmos_refstate::atmos_refstate_setup(), scale_atmos_saturation::atmos_saturation_moist_conversion_dens_all_0d(), scale_atmos_saturation::atmos_saturation_pote_1d(), scale_atmos_saturation::atmos_saturation_setup(), scale_atmos_saturation::atmos_saturation_tdew_liq_1d(), scale_atmos_sfc_ch_rn222::atmos_sfc_ch_rn222_land_flux(), scale_atmos_sfc_ch_rn222::atmos_sfc_ch_rn222_setup(), scale_atmos_solarins::atmos_solarins_setup(), mod_atmos_vars::atmos_vars_check(), mod_atmos_vars::atmos_vars_get_diagnostic_1d(), mod_atmos_vars::atmos_vars_get_diagnostic_2d(), mod_atmos_vars::atmos_vars_get_diagnostic_3d(), mod_atmos_vars::atmos_vars_restart_open(), mod_atmos_vars::atmos_vars_restart_read(), mod_atmos_vars::atmos_vars_setup(), scale_bulkflux::bulkflux_setup(), scale_land_phy_snow_ky90::calculationmo(), scale_calendar::calendar_cfunits2sec(), scale_calendar::calendar_sec2unit(), scale_calendar::calendar_setup(), scale_calendar::calendar_unit2sec(), scale_debug::check(), mod_cnv2d::cnv2d_exec(), mod_cnv2d::cnv2d_setup(), mod_cnv2d::cnv2d_tile_init(), mod_cnvlanduse::cnvlanduse(), mod_cnvlanduse::cnvlanduse_setup(), mod_cnvtopo::cnvtopo(), mod_cnvtopo::cnvtopo_setup(), mod_cnvuser::cnvuser_setup(), scale_comm_cartesc_nest::comm_cartesc_nest_domain_relate(), scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape(), scale_comm_cartesc_nest::comm_cartesc_nest_intercomm_nestdown_3d(), scale_comm_cartesc_nest::comm_cartesc_nest_issuer_of_receive_3d(), scale_comm_cartesc_nest::comm_cartesc_nest_issuer_of_wait_3d(), scale_comm_cartesc_nest::comm_cartesc_nest_nestdown(), scale_comm_cartesc_nest::comm_cartesc_nest_recv_cancel(), scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue(), scale_comm_cartesc_nest::comm_cartesc_nest_setup(), scale_comm_cartesc_nest::comm_cartesc_nest_test(), scale_comm_cartesc_nest::comm_cartesc_nest_waitall(), scale_comm_icoa::comm_data_transfer_dp(), scale_comm_icoa::comm_data_transfer_nopl(), scale_comm_icoa::comm_data_transfer_sp(), scale_comm_icoa::comm_setup(), scale_comm_cartesc::comm_vars8_init(), scale_comm_cartesc::comm_vars_init(), scale_const::const_setup(), mod_convert::convert(), mod_convert::convert_setup(), mod_copytopo::copytopo(), scale_coriolis::coriolis_setup(), mod_cpl_admin::cpl_admin_setup(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin_setup(), mod_cpl_vars::cpl_vars_setup(), scale_file::file_add_associatedvariable(), scale_file::file_attach_buffer(), scale_file_cartesc::file_cartesc_create(), scale_file_cartesc::file_cartesc_def_var(), scale_file_cartesc::file_cartesc_read_1d(), scale_file_cartesc::file_cartesc_read_auto_2d(), scale_file_cartesc::file_cartesc_read_auto_3d(), scale_file_cartesc::file_cartesc_read_var_2d(), scale_file_cartesc::file_cartesc_read_var_3d(), scale_file_cartesc::file_cartesc_read_var_4d(), scale_file_cartesc::file_cartesc_setup(), scale_file_cartesc::file_cartesc_write_var_1d(), scale_file_cartesc::file_cartesc_write_var_2d(), scale_file_cartesc::file_cartesc_write_var_3d(), scale_file_cartesc::file_cartesc_write_var_3d_t(), scale_file_cartesc::file_cartesc_write_var_4d(), scale_file::file_close(), scale_file::file_create(), scale_file::file_def_associatedcoordinate(), scale_file::file_def_axis(), scale_file::file_def_variable(), scale_file::file_detach_buffer(), scale_file::file_enddef(), scale_file_external_input_cartesc::file_external_input_cartesc_setup(), scale_file_external_input::file_external_input_regist(), scale_file_external_input::file_external_input_setup(), scale_file_external_input::file_external_input_update_3d(), scale_file::file_flush(), scale_file::file_get_dimlength(), scale_file::file_get_stepsize(), scale_file_grads::file_grads_close(), scale_file_grads::file_grads_get_shape_id(), scale_file_grads::file_grads_get_shape_name(), scale_file_grads::file_grads_isoned(), scale_file_grads::file_grads_open(), scale_file_grads::file_grads_read_1d_id(), scale_file_grads::file_grads_read_2d_id(), scale_file_grads::file_grads_varid(), scale_file_history_cartesc::file_history_cartesc_truncate_1d(), scale_file_history_cartesc::file_history_cartesc_truncate_3d(), scale_file_history::file_history_finalize(), scale_file_history::file_history_in_0d(), scale_file_history::file_history_in_1d(), scale_file_history::file_history_in_2d(), scale_file_history::file_history_in_3d(), scale_file_history::file_history_reg(), scale_file_history::file_history_set_associatedcoordinate_1d(), scale_file_history::file_history_set_associatedcoordinate_2d(), scale_file_history::file_history_set_associatedcoordinate_3d(), scale_file_history::file_history_set_attribute_double(), scale_file_history::file_history_set_attribute_float(), scale_file_history::file_history_set_attribute_int(), scale_file_history::file_history_set_attribute_int_ary(), scale_file_history::file_history_set_attribute_text(), scale_file_history::file_history_set_axis(), scale_file_history::file_history_set_dim(), scale_file_history::file_history_setup(), scale_file::file_make_fname(), scale_file::file_redef(), scale_file::file_set_option(), scale_file_tiledata::file_tiledata_get_data_int1(), scale_file_tiledata::file_tiledata_get_info(), scale_file_tiledata::file_tiledata_get_latlon(), scale_file_tiledata::file_tiledata_read_catalog_file(), scale_file_tiledata::file_tiledata_read_data_int1_int(), scale_file_tiledata::file_tiledata_read_data_int2_int(), scale_file_tiledata::file_tiledata_read_data_int2_real(), scale_file_tiledata::file_tiledata_read_data_int4_int(), scale_file_tiledata::file_tiledata_read_data_int4_real(), scale_file_tiledata::file_tiledata_read_data_real4_int(), scale_file_tiledata::file_tiledata_read_data_real4_real(), scale_file_tiledata::file_tiledata_read_data_real8_real(), scale_filter::filter_hyperdiff_3d(), mod_mkinit::flux_setup(), scale_interp::interp_div_block(), scale_interp::interp_domain_compatibility(), scale_interp::interp_factor1d(), scale_interp::interp_factor2d_linear_xy(), scale_interp::interp_factor2d_weight(), scale_interp::interp_factor3d_weight(), scale_interp::interp_setup(), mod_lake_admin::lake_admin_setup(), mod_land_admin::land_admin_setup(), mod_land_driver::land_driver_setup(), scale_land_dyn_bucket::land_dyn_bucket(), scale_land_dyn_bucket::land_dyn_bucket_setup(), scale_land_grid_cartesc_index::land_grid_cartesc_index_setup(), scale_land_grid_cartesc::land_grid_cartesc_setup(), scale_land_grid_icoa_index::land_grid_icoa_index_setup(), scale_land_grid_icoa::land_grid_icoa_setup(), mod_realinput::land_interporation(), scale_land_phy_matsiro::land_phy_matsiro_setup(), scale_land_phy_snow_ky90::land_phy_snow_ky90(), scale_land_phy_snow_ky90::land_phy_snow_ky90_setup(), mod_mkinit::land_setup(), mod_land_vars::land_vars_check(), mod_land_vars::land_vars_restart_read(), mod_land_vars::land_vars_setup(), scale_landuse::landuse_fillhalo(), scale_landuse::landuse_setup(), scale_atmos_dyn_tstep_short_fvm_hivi::make_matrix(), scale_mapprojection::mapprojection_get_param(), scale_mapprojection::mapprojection_get_param_equidistantcylindrical(), scale_mapprojection::mapprojection_get_param_lambertconformal(), scale_mapprojection::mapprojection_get_param_mercator(), scale_mapprojection::mapprojection_lonlat2xy_0d_initialized(), scale_mapprojection::mapprojection_lonlat2xy_2d_param(), scale_mapprojection::mapprojection_mapfactor_initialized(), scale_mapprojection::mapprojection_rotcoef_initialized(), scale_mapprojection::mapprojection_setup(), scale_mapprojection::mapprojection_xy2lonlat_0d_initialized(), scale_mapprojection::mapprojection_xy2lonlat_2d_param(), scale_mapprojection::mapprojection_xy2lonlat_equidistantcylindrical(), mod_mkinit::mkinit(), mod_mkinit::mkinit_setup(), mod_mktopo::mktopo(), mod_mktopo::mktopo_setup(), scale_monitor::monitor_reg(), scale_monitor::monitor_setup(), scale_monitor::monitor_write(), scale_atmos_phy_mp_sn14::nucleation(), mod_ocean_admin::ocean_admin_setup(), mod_ocean_driver::ocean_driver_setup(), scale_ocean_dyn_offline::ocean_dyn_offline(), scale_ocean_dyn_offline::ocean_dyn_offline_setup(), scale_ocean_dyn_slab::ocean_dyn_slab(), scale_ocean_dyn_slab::ocean_dyn_slab_setup(), scale_ocean_grid_cartesc_index::ocean_grid_cartesc_index_setup(), scale_ocean_grid_cartesc::ocean_grid_cartesc_setup(), scale_ocean_grid_icoa_index::ocean_grid_icoa_index_setup(), scale_ocean_grid_icoa::ocean_grid_icoa_setup(), mod_realinput::ocean_interporation(), scale_ocean_phy_albedo::ocean_phy_albedo_const_setup(), scale_ocean_phy_albedo::ocean_phy_albedo_seaice_setup(), scale_ocean_phy_ice_simple::ocean_phy_ice_setup(), scale_ocean_phy_ice_simple::ocean_phy_ice_simple(), scale_ocean_phy_roughness::ocean_phy_roughness_const_setup(), scale_ocean_phy_roughness_miller92::ocean_phy_roughness_miller92_setup(), scale_ocean_phy_roughness_moon07::ocean_phy_roughness_moon07_setup(), scale_ocean_phy_roughness::ocean_phy_roughness_seaice_setup(), scale_ocean_phy_roughness::ocean_phy_roughness_setup(), scale_ocean_phy_tc::ocean_phy_tc_seaice_setup(), mod_mkinit::ocean_setup(), mod_ocean_vars::ocean_vars_restart_read(), mod_ocean_vars::ocean_vars_setup(), mod_realinput_grads::parentatmosinputgrads(), mod_realinput_scale::parentatmosopenscale(), mod_realinput_grads::parentatmossetupgrads(), mod_realinput_wrfarw::parentatmossetupwrfarw(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_grads::parentlandsetupgrads(), mod_realinput_wrfarw::parentlandsetupwrfarw(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_grads::parentoceansetupgrads(), mod_realinput_wrfarw::parentoceansetupwrfarw(), scale_prc_cartesc::prc_cartesc_setup(), prc_errhandler_setup(), scale_prc_icoa::prc_icoa_rgn_generate(), scale_prc_icoa::prc_icoa_setup(), prc_mpisplit_bulk(), prc_mpisplit_nest(), scale_prof::prof_setup(), scale_random::random_setup(), scale_file_grads::read_data_int1(), mod_mkinit::read_sounding(), mod_realinput::realinput_atmos(), mod_realinput::realinput_surface(), mod_mkinit::rect_setup(), mod_realinput::replace_misval_map(), mod_rm_driver::rm_driver(), scalerm(), scalerm_init(), scalerm_pp(), scale_file_cartesc::set_dimension(), scale_atmos_dyn_tstep_short_fvm_hivi::solve_bicgstab(), scale_atmos_dyn_tstep_short_fvm_hevi::solve_direct(), scale_spnudge::spnudge_setup(), scale_statistics::statistics_setup(), scale_statistics::statistics_total_2d(), scale_statistics::statistics_total_3d(), mod_mkinit::tke_setup(), scale_topography::topography_fillhalo(), scale_topography::topography_setup(), scale_tracer::tracer_regist(), mod_urban_admin::urban_admin_setup(), mod_urban_driver::urban_driver_setup(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01_setup(), scale_urban_grid_cartesc_index::urban_grid_cartesc_index_setup(), scale_urban_grid_cartesc::urban_grid_cartesc_setup(), scale_urban_grid_icoa_index::urban_grid_icoa_index_setup(), scale_urban_grid_icoa::urban_grid_icoa_setup(), mod_mkinit::urban_setup(), mod_urban_vars::urban_vars_restart_read(), mod_urban_vars::urban_vars_setup(), mod_user::user_setup(), mod_user::user_update(), scale_debug::valcheck_1d(), scale_debug::valcheck_2d(), scale_debug::valcheck_3d(), scale_comm_cartesc::vars8_2d_mpi(), scale_comm_cartesc::vars8_3d_mpi(), scale_comm_cartesc::vars_2d_mpi(), scale_comm_cartesc::vars_3d_mpi(), and scale_comm_cartesc::vars_3d_mpi_pc().

◆ prc_mpifinish()

subroutine, public scale_prc::prc_mpifinish

Stop MPI peacefully.

Definition at line 358 of file scale_prc.F90.

358  implicit none
359 
360  integer :: ierr
361  logical :: sign_status
362  logical :: sign_exit
363  !---------------------------------------------------------------------------
364 
365  ! FPM polling
366  if ( fpm_alive ) then
367  sign_status = .false.
368  sign_exit = .false.
369  do while ( .NOT. sign_exit )
370  call fpm_polling( sign_status, sign_exit )
371  enddo
372  endif
373 
374  if (prc_universal_handler .NE. mpi_errhandler_null) then
375  call mpi_errhandler_free(prc_universal_handler, ierr)
376  endif
377  if (prc_abort_handler .NE. mpi_errhandler_null) then
378  call mpi_errhandler_free(prc_abort_handler, ierr)
379  endif
380 
381  ! Stop MPI
382  if ( prc_mpi_alive ) then
383  log_newline
384  log_progress(*) 'finalize MPI...'
385 
386  ! free splitted communicator
387  if ( prc_local_comm_world /= prc_global_comm_world ) then
388  call mpi_comm_free(prc_local_comm_world,ierr)
389  endif
390 
391  call mpi_barrier(prc_universal_comm_world,ierr)
392 
393  call mpi_finalize(ierr)
394  log_progress(*) 'MPI is peacefully finalized'
395  endif
396 
397  ! Close logfile, configfile
398  if ( io_l ) then
399  if( io_fid_log /= io_fid_stdout ) close(io_fid_log)
400  endif
401  close(io_fid_conf)
402 
403  return

References scale_fpm::fpm_alive, scale_fpm::fpm_polling(), scale_io::io_fid_conf, scale_io::io_fid_log, scale_io::io_fid_stdout, scale_io::io_l, prc_abort_handler, prc_global_comm_world, prc_local_comm_world, prc_mpi_alive, prc_universal_comm_world, and prc_universal_handler.

Referenced by scale::scale_finalize(), scalerm(), scalerm_init(), scalerm_pp(), and scale_comm_icoa::suf().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ prc_mpisplit_bulk()

subroutine, public scale_prc::prc_mpisplit_bulk ( integer, intent(in)  ORG_COMM_WORLD,
integer, intent(in)  NUM_BULKJOB,
integer, dimension(:), intent(in)  PRC_BULKJOB,
logical, intent(in)  debug,
integer, intent(out)  SUB_COMM_WORLD,
integer, intent(out)  ID_BULKJOB 
)

MPI Communicator Split (bulk job)

Definition at line 415 of file scale_prc.F90.

415  implicit none
416 
417  integer, intent(in) :: ORG_COMM_WORLD ! communicator (original group)
418  integer, intent(in) :: NUM_BULKJOB ! number of bulk jobs
419  integer, intent(in) :: PRC_BULKJOB(:) ! number of ranks in subgroup communicator
420  logical, intent(in) :: debug ! log-output for mpi splitting?
421  integer, intent(out) :: SUB_COMM_WORLD ! communicator (new subgroup)
422  integer, intent(out) :: ID_BULKJOB ! bulk job id
423 
424  integer :: ORG_myrank ! my rank in the original communicator
425  integer :: ORG_nrank ! number of ranks in the original communicator
426  integer :: sum_nrank
427 
428  integer, allocatable :: prc2color (:) ! color id for each process
429  integer, allocatable :: prc2key (:) ! local rank id for each process
430  integer :: COL_domain(0:PRC_DOMAIN_nlim) ! domain id of this color
431  integer :: COL_master(0:PRC_DOMAIN_nlim) ! master rank of this color
432  integer :: COL_parent(0:PRC_DOMAIN_nlim) ! parent color of this color
433  integer :: COL_child (0:PRC_DOMAIN_nlim) ! child color of this color
434 
435  logical :: color_reorder = .false. ! dummy
436 
437  integer :: i, color
438  integer :: itag, ierr
439  !---------------------------------------------------------------------------
440 
441 
442  if ( num_bulkjob == 1 ) then ! single domain run
443 
444  sub_comm_world = org_comm_world
445  id_bulkjob = 0
446  prc_universal_jobid = 0
447 
448  elseif( num_bulkjob > 1 ) then ! multi domain run
449 
450  call mpi_comm_rank(org_comm_world,org_myrank,ierr)
451  call mpi_comm_size(org_comm_world,org_nrank ,ierr)
452 
453  sum_nrank = sum(prc_bulkjob(1:num_bulkjob))
454 
455  if ( sum_nrank /= org_nrank ) then
456  if ( prc_universal_ismaster ) then
457  log_error("PRC_MPIsplit",*) "MPI PROCESS NUMBER is INCONSISTENT"
458  log_error_cont(*) " REQUESTED NPROCS = ", sum_nrank, " LAUNCHED NPROCS = ", org_nrank
459  endif
460  call prc_abort
461  endif
462 
463  allocate( prc2color(0:org_nrank-1) )
464  allocate( prc2key(0:org_nrank-1) )
465 
466  call prc_mpicoloring( org_comm_world, & ! [IN]
467  org_nrank, & ! [IN]
468  num_bulkjob, & ! [IN]
469  prc_bulkjob(:), & ! [IN]
470  color_reorder, & ! [IN]
471  .true., & ! [IN]
472  debug, & ! [IN]
473  prc2color(:), & ! [OUT]
474  prc2key(:), & ! [OUT]
475  col_domain(:), & ! [OUT]
476  col_master(:), & ! [OUT]
477  col_parent(:), & ! [OUT]
478  col_child(:) ) ! [OUT]
479 
480  ! split comm_world
481  call mpi_comm_split( org_comm_world, & ! [IN]
482  prc2color(org_myrank), & ! [IN]
483  prc2key(org_myrank), & ! [IN]
484  sub_comm_world, & ! [OUT]
485  ierr ) ! [OUT]
486 
487  do i = 1, num_bulkjob
488  color = i-1
489  prc_global_root(i) = col_master(color)
490  enddo
491 
492  id_bulkjob = prc2color(org_myrank) ! color = bulk id
493  prc_universal_jobid = prc2color(org_myrank) ! color = bulk id
494 
495  deallocate( prc2color )
496  deallocate( prc2key )
497 
498  else
499  if ( prc_universal_ismaster ) then
500  log_error("PRC_MPIsplit",*) "REQUESTED DOMAIN NUMBER IS NOT ACCEPTABLE"
501  endif
502  call prc_abort
503  endif
504 
505  return

References prc_abort(), prc_global_root, prc_universal_ismaster, and prc_universal_jobid.

Referenced by scalerm(), scalerm_init(), and scalerm_pp().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ prc_mpisplit_nest()

subroutine, public scale_prc::prc_mpisplit_nest ( integer, intent(in)  ORG_COMM_WORLD,
integer, intent(in)  NUM_DOMAIN,
integer, dimension(:), intent(in)  PRC_DOMAIN,
logical, intent(in)  debug,
logical, intent(in)  color_reorder,
integer, intent(out)  SUB_COMM_WORLD,
integer, intent(out)  ID_DOMAIN,
integer, intent(out)  INTERCOMM_parent,
integer, intent(out)  INTERCOMM_child 
)

MPI Communicator Split (nesting)

Definition at line 520 of file scale_prc.F90.

520  implicit none
521 
522  integer, intent(in) :: ORG_COMM_WORLD ! communicator (original group)
523  integer, intent(in) :: NUM_DOMAIN ! number of bulk jobs
524  integer, intent(in) :: PRC_DOMAIN(:) ! number of ranks in subgroup communicator
525  logical, intent(in) :: debug ! log-output for mpi splitting?
526  logical, intent(in) :: color_reorder ! reorder
527  integer, intent(out) :: SUB_COMM_WORLD ! communicator (new subgroup)
528  integer, intent(out) :: ID_DOMAIN ! domain id
529  integer, intent(out) :: INTERCOMM_parent ! communicator between this rank and parent domain
530  integer, intent(out) :: INTERCOMM_child ! communicator between this rank and child domain
531 
532  integer :: ORG_myrank ! my rank in the original communicator
533  integer :: ORG_nrank ! number of ranks in the original communicator
534  integer :: sum_nrank
535 
536  integer, allocatable :: prc2color (:) ! color id for each process
537  integer, allocatable :: prc2key (:) ! local rank id for each process
538  integer :: COL_domain(0:PRC_DOMAIN_nlim) ! domain id of this color
539  integer :: COL_master(0:PRC_DOMAIN_nlim) ! master rank of this color
540  integer :: COL_parent(0:PRC_DOMAIN_nlim) ! parent color of this color
541  integer :: COL_child (0:PRC_DOMAIN_nlim) ! child color of this color
542 
543  integer :: i, color
544  integer :: itag, ierr
545  !---------------------------------------------------------------------------
546 
547  intercomm_parent = mpi_comm_null
548  intercomm_child = mpi_comm_null
549 
550  if ( num_domain == 1 ) then ! single domain run
551 
552  sub_comm_world = org_comm_world
553  id_domain = 1
554  prc_global_domainid = 1
555 
556  elseif( num_domain > 1 ) then ! multi domain run
557 
558  call mpi_comm_rank(org_comm_world,org_myrank,ierr)
559  call mpi_comm_size(org_comm_world,org_nrank ,ierr)
560 
561  sum_nrank = sum(prc_domain(1:num_domain))
562 
563  if ( sum_nrank /= org_nrank ) then
564  if ( prc_universal_ismaster ) then
565  log_error("PRC_MPIsplit",*) "MPI PROCESS NUMBER is INCONSISTENT"
566  log_error_cont(*) " REQUESTED NPROCS = ", sum_nrank, " LAUNCHED NPROCS = ", org_nrank
567  endif
568  call prc_abort
569  endif
570 
571  allocate( prc2color(0:org_nrank-1) )
572  allocate( prc2key(0:org_nrank-1) )
573 
574  call prc_mpicoloring( org_comm_world, & ! [IN]
575  org_nrank, & ! [IN]
576  num_domain, & ! [IN]
577  prc_domain(:), & ! [IN]
578  color_reorder, & ! [IN]
579  .false., & ! [IN]
580  debug, & ! [IN]
581  prc2color(:), & ! [OUT]
582  prc2key(:), & ! [OUT]
583  col_domain(:), & ! [OUT]
584  col_master(:), & ! [OUT]
585  col_parent(:), & ! [OUT]
586  col_child(:) ) ! [OUT]
587 
588  ! split comm_world
589  call mpi_comm_split( org_comm_world, & ! [IN]
590  prc2color(org_myrank), & ! [IN]
591  prc2key(org_myrank), & ! [IN]
592  sub_comm_world, & ! [OUT]
593  ierr ) ! [OUT]
594 
595  id_domain = col_domain(prc2color(org_myrank)) ! color /= domain id
596  prc_global_domainid = col_domain(prc2color(org_myrank)) ! color /= domain id
597 
598  if ( prc_universal_ismaster ) then
599  write(*,*)
600  write(*,*) "INFO [PRC_MPIsplit] Inter-domain relationship information"
601  do i = 1, num_domain
602  color = i-1
603  if ( col_parent(color) >= 0 ) then
604  write(*,'(5x,A,I2.2)') "Relationship No.", i
605  write(*,'(5x,A,I2.2,A,A,I2.2)') "Parent color: ", col_parent(color), " <=> ", &
606  "Child color: ", color
607  endif
608  enddo
609  endif
610 
611  ! create inter-communicator
612  do i = 1, num_domain
613  color = i-1
614  if ( col_parent(color) >= 0 ) then
615  itag = i*100
616 
617  if ( prc2color(org_myrank) == col_parent(color) ) then ! as a parent
618 
619  call mpi_intercomm_create( sub_comm_world, prc_masterrank, &
620  org_comm_world, col_master(color), &
621  itag, intercomm_child, ierr )
622 
623  elseif( prc2color(org_myrank) == color ) then ! as a child
624 
625  call mpi_intercomm_create( sub_comm_world, prc_masterrank, &
626  org_comm_world, col_master(col_parent(color)), &
627  itag, intercomm_parent, ierr )
628 
629  endif
630 
631  call mpi_barrier(org_comm_world,ierr)
632  endif
633  enddo
634 
635  deallocate( prc2color )
636  deallocate( prc2key )
637 
638  else
639  if ( prc_universal_ismaster ) then
640  write(*,*)"ERROR [RPC_MPIsplit] REQUESTED DOMAIN NUMBER IS NOT ACCEPTABLE"
641  endif
642  call prc_abort
643  endif
644 
645  return

References prc_abort(), prc_global_domainid, prc_masterrank, and prc_universal_ismaster.

Referenced by scalerm(), scalerm_init(), and scalerm_pp().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ prc_mpibarrier()

subroutine, public scale_prc::prc_mpibarrier

Barrier MPI.

Definition at line 831 of file scale_prc.F90.

831  implicit none
832 
833  integer :: ierr
834  !---------------------------------------------------------------------------
835 
836  if ( prc_mpi_alive ) then
837  call mpi_barrier(prc_local_comm_world,ierr)
838  endif
839 

References prc_local_comm_world, and prc_mpi_alive.

Referenced by scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_f2013(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_mg2001(), scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Here is the caller graph for this function:

◆ prc_mpitime()

real(dp) function, public scale_prc::prc_mpitime

Get MPI time.

Returns
time

Definition at line 846 of file scale_prc.F90.

846  implicit none
847 
848  real(DP) :: time
849  !---------------------------------------------------------------------------
850 
851  if ( prc_mpi_alive ) then
852  time = real(mpi_wtime(), kind=dp)
853  else
854  call cpu_time(time)
855  endif
856 

References scale_precision::dp, and prc_mpi_alive.

Referenced by mod_admin_time::admin_time_advance(), mod_admin_time::admin_time_checkstate(), mod_admin_time::admin_time_setup(), scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

Here is the caller graph for this function:

◆ prc_mpitimestat()

subroutine, public scale_prc::prc_mpitimestat ( real(dp), dimension(:), intent(out)  avgvar,
real(dp), dimension(:), intent(out)  maxvar,
real(dp), dimension(:), intent(out)  minvar,
integer, dimension(:), intent(out)  maxidx,
integer, dimension(:), intent(out)  minidx,
real(dp), dimension(:), intent(in)  var 
)

Calc global statistics for timer.

Parameters
[out]avgvaraverage
[out]maxvarmaximum
[out]minvarminimum
[out]maxidxindex of maximum
[out]minidxindex of minimum
[in]varvalues for statistics

Definition at line 868 of file scale_prc.F90.

868  implicit none
869 
870  real(DP), intent(out) :: avgvar(:)
871  real(DP), intent(out) :: maxvar(:)
872  real(DP), intent(out) :: minvar(:)
873  integer, intent(out) :: maxidx(:)
874  integer, intent(out) :: minidx(:)
875  real(DP), intent(in) :: var(:)
876 
877  real(DP), allocatable :: statval(:,:)
878  integer :: vsize
879 
880  real(DP) :: totalvar
881  integer :: ierr
882  integer :: v, p
883  !---------------------------------------------------------------------------
884 
885  vsize = size(var(:))
886 
887  allocate( statval(vsize,0:prc_nprocs-1) )
888  statval(:,:) = 0.0_dp
889 
890  do v = 1, vsize
891  statval(v,prc_myrank) = var(v)
892  enddo
893 
894  ! MPI broadcast
895  do p = 0, prc_nprocs-1
896  call mpi_bcast( statval(1,p), &
897  vsize, &
898  mpi_double_precision, &
899  p, &
900  prc_local_comm_world, &
901  ierr )
902  enddo
903 
904  do v = 1, vsize
905  totalvar = 0.0_dp
906  do p = 0, prc_nprocs-1
907  totalvar = totalvar + statval(v,p)
908  enddo
909  avgvar(v) = totalvar / prc_nprocs
910 
911  maxvar(v) = maxval(statval(v,:))
912  minvar(v) = minval(statval(v,:))
913  maxidx(v:v) = maxloc(statval(v,:))
914  minidx(v:v) = minloc(statval(v,:))
915  enddo
916 
917  deallocate( statval )
918 
919  return

References scale_fpm::fpm_alive, scale_fpm::fpm_polling(), scale_io::io_fid_conf, scale_io::io_fid_log, scale_io::io_fid_stdout, scale_io::io_l, prc_abort_comm_world, prc_global_domainid, prc_global_myrank, prc_global_nprocs, prc_ismaster, prc_local_comm_world, prc_mpi_alive, prc_myrank, prc_nprocs, prc_universal_jobid, prc_universal_myrank, and prc_universal_nprocs.

Referenced by scale_prof::prof_rapreport().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ prc_set_file_closer()

subroutine, public scale_prc::prc_set_file_closer ( procedure(closer)  routine)

Definition at line 1005 of file scale_prc.F90.

1005  procedure(closer) :: routine
1006 
1007  prc_file_closer => routine
1008 
1009  return

Referenced by scale_file::file_setup().

Here is the caller graph for this function:

Variable Documentation

◆ prc_masterrank

integer, parameter, public scale_prc::prc_masterrank = 0

◆ prc_domain_nlim

integer, parameter, public scale_prc::prc_domain_nlim = 10000

max depth of domains

Definition at line 67 of file scale_prc.F90.

67  integer, public, parameter :: PRC_DOMAIN_nlim = 10000

Referenced by scalerm(), scalerm_init(), and scalerm_pp().

◆ prc_comm_null

integer, parameter, public scale_prc::prc_comm_null = MPI_COMM_NULL

Definition at line 68 of file scale_prc.F90.

68  integer, public, parameter :: PRC_COMM_NULL = mpi_comm_null

Referenced by scale_file::file_get_aggregate(), scalerm(), scalerm_init(), and scalerm_pp().

◆ prc_universal_comm_world

integer, public scale_prc::prc_universal_comm_world = -1

original communicator

Definition at line 71 of file scale_prc.F90.

71  integer, public :: PRC_UNIVERSAL_COMM_WORLD = -1

Referenced by prc_abort(), scale_prc_cartesc::prc_cartesc_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpifinish(), prc_singlecom_setup(), and prc_universal_setup().

◆ prc_universal_myrank

integer, public scale_prc::prc_universal_myrank = -1

myrank in universal communicator

Definition at line 72 of file scale_prc.F90.

72  integer, public :: PRC_UNIVERSAL_myrank = -1

Referenced by scale_prc_cartesc::prc_cartesc_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpitimestat(), prc_singlecom_setup(), and prc_universal_setup().

◆ prc_universal_nprocs

integer, public scale_prc::prc_universal_nprocs = -1

process num in universal communicator

Definition at line 73 of file scale_prc.F90.

73  integer, public :: PRC_UNIVERSAL_nprocs = -1

Referenced by scale_prc_cartesc::prc_cartesc_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpitimestat(), prc_singlecom_setup(), and prc_universal_setup().

◆ prc_universal_ismaster

logical, public scale_prc::prc_universal_ismaster = .false.

master process in universal communicator?

Definition at line 74 of file scale_prc.F90.

74  logical, public :: PRC_UNIVERSAL_IsMaster = .false.

Referenced by mod_admin_time::admin_time_checkstate(), scale_prc_cartesc::prc_cartesc_setup(), prc_errhandler_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpisplit_bulk(), prc_mpisplit_nest(), prc_singlecom_setup(), and prc_universal_setup().

◆ prc_universal_jobid

integer, public scale_prc::prc_universal_jobid = 0

my job ID in universal communicator

Definition at line 76 of file scale_prc.F90.

76  integer, public :: PRC_UNIVERSAL_jobID = 0

Referenced by prc_mpisplit_bulk(), and prc_mpitimestat().

◆ prc_global_comm_world

integer, public scale_prc::prc_global_comm_world = -1

global communicator

Definition at line 79 of file scale_prc.F90.

79  integer, public :: PRC_GLOBAL_COMM_WORLD = -1

Referenced by scale_comm_cartesc_nest::comm_cartesc_nest_disconnect(), scale_prc_cartesc::prc_cartesc_setup(), prc_global_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpifinish(), and prc_singlecom_setup().

◆ prc_global_myrank

integer, public scale_prc::prc_global_myrank = -1

myrank in global communicator

Definition at line 80 of file scale_prc.F90.

80  integer, public :: PRC_GLOBAL_myrank = -1

Referenced by scale_prc_cartesc::prc_cartesc_setup(), prc_global_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpitimestat(), and prc_singlecom_setup().

◆ prc_global_nprocs

integer, public scale_prc::prc_global_nprocs = -1

process num in global communicator

Definition at line 81 of file scale_prc.F90.

81  integer, public :: PRC_GLOBAL_nprocs = -1

Referenced by scale_prc_cartesc::prc_cartesc_setup(), prc_global_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpitimestat(), and prc_singlecom_setup().

◆ prc_global_ismaster

logical, public scale_prc::prc_global_ismaster = .false.

master process in global communicator?

Definition at line 82 of file scale_prc.F90.

82  logical, public :: PRC_GLOBAL_IsMaster = .false.

Referenced by scale_prc_cartesc::prc_cartesc_setup(), prc_global_setup(), scale_prc_icoa::prc_icoa_setup(), and prc_singlecom_setup().

◆ prc_global_domainid

integer, public scale_prc::prc_global_domainid = 0

my domain ID in global communicator

Definition at line 84 of file scale_prc.F90.

84  integer, public :: PRC_GLOBAL_domainID = 0

Referenced by scale_comm_cartesc_nest::comm_cartesc_nest_setup(), prc_mpisplit_nest(), and prc_mpitimestat().

◆ prc_global_root

integer, dimension(prc_domain_nlim), public scale_prc::prc_global_root

root processes in global members

Definition at line 85 of file scale_prc.F90.

85  integer, public :: PRC_GLOBAL_ROOT(PRC_DOMAIN_nlim)

Referenced by prc_mpisplit_bulk(), scalerm(), scalerm_init(), and scalerm_pp().

◆ prc_local_comm_world

integer, public scale_prc::prc_local_comm_world = -1

◆ prc_nprocs

integer, public scale_prc::prc_nprocs = 1

◆ prc_myrank

integer, public scale_prc::prc_myrank = 0

process num in local communicator

Definition at line 90 of file scale_prc.F90.

90  integer, public :: PRC_myrank = 0

Referenced by mod_admin_time::admin_time_setup(), scale_atmos_grid_cartesc::atmos_grid_cartesc_allocate(), scale_atmos_grid_cartesc::atmos_grid_cartesc_generate(), scale_atmos_grid_cartesc_index::atmos_grid_cartesc_index_setup_main(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_electric_field(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_f2013(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_neutralization_mg2001(), scale_atmos_phy_lt_sato2019::atmos_phy_lt_sato2019_setup(), mod_atmos_vars::atmos_vars_check(), mod_atmos_vars::atmos_vars_restart_check(), scale_comm_cartesc_nest::comm_cartesc_nest_domain_relate(), scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape(), scale_comm_cartesc_nest::comm_cartesc_nest_issuer_of_receive_3d(), scale_comm_icoa::comm_setup(), mod_copytopo::copytopo(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), scale_atmos_phy_mp_sn14::debug_tem(), scale_dft::dft_setup(), scale_file_cartesc::file_cartesc_create(), scale_file_cartesc::file_cartesc_open(), scale_file_cartesc::file_cartesc_setup(), scale_file_cartesc::file_cartesc_write_axes(), scale_file_cartesc::file_cartesc_write_var_1d(), scale_file_cartesc::file_cartesc_write_var_2d(), scale_file_cartesc::file_cartesc_write_var_3d(), scale_file_cartesc::file_cartesc_write_var_3d_t(), scale_file_cartesc::file_cartesc_write_var_4d(), scale_file_external_input::file_external_input_regist(), scale_file_external_input::file_external_input_update_3d(), scale_file_history_cartesc::file_history_cartesc_set_pres(), scale_file_history_cartesc::file_history_cartesc_setup(), scale_file_history_cartesc::file_history_cartesc_truncate_3d(), scale_land_grid_cartesc::land_grid_cartesc_setup(), scale_land_grid_icoa::land_grid_icoa_setup(), scale_monitor::monitor_finalize(), scale_monitor::monitor_write(), scale_ocean_grid_cartesc::ocean_grid_cartesc_setup(), scale_ocean_grid_icoa::ocean_grid_icoa_setup(), scale_prc_cartesc::prc_cartesc_setup(), scale_prc_icoa::prc_icoa_rgn_generate(), scale_prc_icoa::prc_icoa_setup(), prc_local_setup(), prc_mpitimestat(), prc_singlecom_setup(), scale_random::random_setup(), scale_statistics::statistics_detail_2d(), scale_statistics::statistics_total_2d(), scale_statistics::statistics_total_3d(), scale_comm_icoa::suf(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01(), scale_urban_dyn_kusaka01::urban_dyn_kusaka01_setup(), scale_urban_grid_cartesc::urban_grid_cartesc_setup(), scale_urban_grid_icoa::urban_grid_icoa_setup(), scale_debug::valcheck_1d(), scale_debug::valcheck_2d(), and scale_debug::valcheck_3d().

◆ prc_ismaster

logical, public scale_prc::prc_ismaster = .false.

◆ prc_mpi_alive

logical, public scale_prc::prc_mpi_alive = .false.

MPI is alive?

Definition at line 94 of file scale_prc.F90.

94  logical, public :: PRC_mpi_alive = .false.

Referenced by prc_abort(), scale_prc_cartesc::prc_cartesc_setup(), prc_mpibarrier(), prc_mpifinish(), prc_mpistart(), prc_mpitime(), and prc_mpitimestat().

◆ prc_universal_handler

integer, public scale_prc::prc_universal_handler

error handler in universal communicator

Definition at line 95 of file scale_prc.F90.

95  integer, public :: PRC_UNIVERSAL_handler

Referenced by prc_errhandler_setup(), and prc_mpifinish().

◆ prc_abort_comm_world

integer, public scale_prc::prc_abort_comm_world

communicator for aborting

Definition at line 96 of file scale_prc.F90.

96  integer, public :: PRC_ABORT_COMM_WORLD

Referenced by scale_prc_cartesc::prc_cartesc_setup(), prc_errhandler_setup(), scale_prc_icoa::prc_icoa_setup(), prc_mpistart(), prc_mpitimestat(), and prc_singlecom_setup().

◆ prc_abort_handler

integer, public scale_prc::prc_abort_handler

error handler communicator for aborting

Definition at line 97 of file scale_prc.F90.

97  integer, public :: PRC_ABORT_handler

Referenced by prc_errhandler_setup(), and prc_mpifinish().

scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:342
scale_fpm::fpm_polling
subroutine, public fpm_polling(run_stat, stop_signal)
Main system of FPM.
Definition: scale_fpm.F90:176