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)
 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_timereorder (rapnlimit, rapnmax, rapttot, rapname)
 reorder rap time 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...
 
integer, public prc_intercomm_parent = MPI_COMM_NULL
 communicator between this rank and parent domain More...
 
integer, public prc_intercomm_child = MPI_COMM_NULL
 communicator between this rank and child domain 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 128 of file scale_prc.F90.

128  implicit none
129 
130  integer, intent(out) :: comm ! communicator
131 
132  integer :: ierr
133  !---------------------------------------------------------------------------
134 
135  call mpi_init(ierr)
136 
137  prc_mpi_alive = .true.
138 ! PRC_UNIVERSAL_handler = MPI_ERRHANDLER_NULL
139 ! call MPI_COMM_CREATE_ERRHANDLER( PRC_MPI_errorhandler, PRC_UNIVERSAL_handler, ierr )
140 
141  comm = mpi_comm_world
142  prc_abort_comm_world = comm
143 
144  return

References prc_abort_comm_world, and prc_mpi_alive.

Referenced by mod_launcher::launcher(), and scale::scale_init().

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 154 of file scale_prc.F90.

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

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

Referenced by mod_launcher::launcher().

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 194 of file scale_prc.F90.

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

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

Referenced by mod_launcher::launcher().

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 229 of file scale_prc.F90.

229  implicit none
230 
231  integer, intent(in) :: comm ! communicator
232  integer, intent(out) :: myrank ! myrank in this communicator
233  logical, intent(out) :: ismaster ! master process in this communicator?
234 
235  integer :: ierr
236  !---------------------------------------------------------------------------
237 
238  prc_local_comm_world = comm
239 
240  call mpi_comm_rank(prc_local_comm_world,prc_myrank,ierr)
241  call mpi_comm_size(prc_local_comm_world,prc_nprocs,ierr)
242  !$acc update device(PRC_myrank)
243 
244  if ( prc_myrank == prc_masterrank ) then
245  prc_ismaster = .true.
246  else
247  prc_ismaster = .false.
248  endif
249 
250  myrank = prc_myrank
251  ismaster = prc_ismaster
252 
253  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 263 of file scale_prc.F90.

263  implicit none
264 
265  integer, intent(in) :: comm ! communicator
266  integer, intent(out) :: nprocs ! number of procs
267  integer, intent(out) :: myrank ! myrank
268  logical, intent(out) :: ismaster ! master process?
269 
270  integer :: ierr
271  !---------------------------------------------------------------------------
272 
273  call mpi_comm_size(comm,nprocs,ierr)
274  call mpi_comm_rank(comm,myrank,ierr)
275 
276  if ( myrank == prc_masterrank ) then
277  ismaster = .true.
278  else
279  ismaster = .false.
280  endif
281 
282  prc_universal_comm_world = comm
283  prc_universal_nprocs = nprocs
284  prc_universal_myrank = myrank
285  prc_universal_ismaster = ismaster
286 
287  prc_global_comm_world = comm
288  prc_global_nprocs = nprocs
289  prc_global_myrank = myrank
290  prc_global_ismaster = ismaster
291 
292  prc_local_comm_world = comm
293  prc_nprocs = nprocs
294  prc_myrank = myrank
295  prc_ismaster = ismaster
296  !$acc update device(PRC_myrank)
297 
298 
299 
300  prc_abort_comm_world = comm
301 
302 ! call MPI_Comm_set_errhandler(PRC_ABORT_COMM_WORLD,PRC_UNIVERSAL_handler,ierr)
303 ! call MPI_Comm_get_errhandler(PRC_ABORT_COMM_WORLD,PRC_ABORT_handler ,ierr)
304 
305  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 313 of file scale_prc.F90.

313  implicit none
314 
315  logical, intent(in) :: use_fpm ! fpm switch
316  logical, intent(in) :: master ! master flag
317 
318  integer :: ierr
319  !---------------------------------------------------------------------------
320 
321  call mpi_comm_create_errhandler(prc_mpi_errorhandler,prc_universal_handler,ierr)
322 
323  call mpi_comm_set_errhandler(prc_abort_comm_world,prc_universal_handler,ierr)
324  call mpi_comm_get_errhandler(prc_abort_comm_world,prc_abort_handler ,ierr)
325 
326  if ( prc_universal_handler /= prc_abort_handler ) then
327  if( prc_universal_ismaster ) write(*,*) ""
328  if( prc_universal_ismaster ) write(*,*) "ERROR: MPI HANDLER is INCONSISTENT"
329  if( prc_universal_ismaster ) write(*,*) " PRC_UNIVERSAL_handler = ", prc_universal_handler
330  if( prc_universal_ismaster ) write(*,*) " PRC_ABORT_handler = ", prc_abort_handler
331  call prc_abort
332  endif
333 
334  if ( use_fpm ) then
335  call sigvars_get_all( master )
336  call signal( sigint, prc_abort )
337 ! call signal( SIGQUIT, PRC_abort )
338 ! call signal( SIGABRT, PRC_abort )
339 ! call signal( SIGFPE, PRC_abort )
340 ! call signal( SIGSEGV, PRC_abort )
341 ! call signal( SIGTERM, PRC_abort )
342  endif
343 
344  return

References prc_abort(), prc_abort_comm_world, prc_abort_handler, prc_universal_handler, prc_universal_ismaster, scale_sigvars::sigint, and scale_sigvars::sigvars_get_all().

Referenced by mod_launcher::launcher(), and scale::scale_init().

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 350 of file scale_prc.F90.

350  implicit none
351 
352  integer :: ierr
353  !---------------------------------------------------------------------------
354 
355  if ( prc_mpi_alive ) then
356  ! tentative approach; input "PRC_UNIVERSAL_COMM_WORLD".
357  call mpi_comm_call_errhandler(prc_universal_comm_world,prc_abort_code,ierr)
358  endif
359 
360  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_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(), 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_3d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_0d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_3d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_2d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_atmos_rev_3d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_3d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_1d(), scale_atmos_hydrostatic::atmos_hydrostatic_buildrho_bytemp_atmos_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_finalize(), mod_atmos_phy_bl_driver::atmos_phy_bl_driver_tracer_setup(), scale_atmos_phy_bl_mynn_jmapplib::atmos_phy_bl_mynn_jmapplib_setup(), scale_atmos_phy_bl_mynn::atmos_phy_bl_mynn_mkinit(), 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_getvar(), scale_atmos_phy_cp_kf_jmapplib::atmos_phy_cp_kf_jmapplib_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_finalize(), 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_calc3d(), scale_atmos_refstate::atmos_refstate_read(), scale_atmos_refstate::atmos_refstate_setup(), scale_atmos_refstate::atmos_refstate_write(), 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_atmos_phy_lt_sato2019::back_sub_ilu(), 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_realinput_netcdf::check_filetype(), 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::comm_bcast_4d_dp(), scale_comm_cartesc::comm_bcast_4d_sp(), scale_comm_cartesc_nest::comm_cartesc_nest_domain_regist_file(), 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_recv(), scale_comm_cartesc_nest::comm_cartesc_nest_nestdown_send(), scale_comm_cartesc_nest::comm_cartesc_nest_parent_info(), scale_comm_cartesc_nest::comm_cartesc_nest_recv_cancel_recv(), scale_comm_cartesc_nest::comm_cartesc_nest_recv_cancel_send(), scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue_recv(), scale_comm_cartesc_nest::comm_cartesc_nest_recvwait_issue_send(), scale_comm_cartesc_nest::comm_cartesc_nest_setup(), scale_comm_cartesc_nest::comm_cartesc_nest_test_recv(), scale_comm_cartesc_nest::comm_cartesc_nest_test_send(), 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_ensemble::comm_ensemble_setup(), scale_comm_cartesc::comm_regist(), scale_comm_icoa::comm_setup(), scale_comm_cartesc::comm_setup(), scale_comm_cartesc::comm_vars8_2d(), scale_comm_cartesc::comm_vars8_3d(), scale_comm_cartesc::comm_vars8_init(), scale_comm_cartesc::comm_vars_2d(), scale_comm_cartesc::comm_vars_3d(), scale_comm_cartesc::comm_vars_init(), scale_comm_cartesc::comm_wait_2d(), scale_comm_cartesc::comm_wait_3d(), 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(), mod_da_admin::da_admin_setup(), mod_da_driver::da_driver_setup(), mod_da_driver::da_driver_update(), mod_da_param_estimation::da_param_estimation_setup(), mod_da_param_estimation::da_param_estimation_update(), mod_da_vars::da_vars_setup(), scale_file::file_add_associatedvariable(), 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_create(), scale_file::file_def_axis(), scale_file_external_input_cartesc::file_external_input_cartesc_finalize(), scale_file_external_input::file_external_input_finalize(), scale_file_external_input::file_external_input_query(), scale_file_external_input::file_external_input_regist_file(), scale_file_external_input::file_external_input_setup(), scale_file::file_get_dimlength(), 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_varcheck(), 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_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(), mod_realinput::get_ijrange(), scale_atmos_phy_lt_sato2019::ilu_decomp(), 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(), mod_launcher::launcher(), scale_letkf::letkf_add_inflation_setup(), scale_letkf::letkf_obs_initialize(), scale_letkf::letkf_obs_readfile(), scale_letkf::letkf_setup(), scale_letkf::letkf_system(), 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(), scale_matrix::matrix_setup(), scale_matrix::matrix_solver_eigenvalue_decomposition(), scale_matrix::matrix_solver_tridiagonal_1d_pcr(), scale_matrix::matrix_solver_tridiagonal_2d_block(), scale_atmos_phy_mp_sn14::mixed_phase_collection_bin(), 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(), scale_atmos_phy_mp_sn14::nucleation_ice_hom(), 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_netcdf::parentatmosinputnetcdf(), mod_realinput_grads::parentatmossetupgrads(), mod_realinput_grads::parentlandinputgrads(), mod_realinput_grads::parentlandsetupgrads(), mod_realinput_netcdf::parentlandsetupnetcdf(), mod_realinput_grads::parentoceaninputgrads(), mod_realinput_grads::parentoceansetupgrads(), mod_realinput_netcdf::parentoceansetupnetcdf(), 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_rapreport(), scale_prof::prof_setup(), scale_random::random_setup(), mod_realinput_netcdf::read1d(), mod_realinput_grads::read2d(), mod_realinput_netcdf::read2d(), mod_realinput_grads::read3d(), mod_realinput_netcdf::read3d(), 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(), scale_file_cartesc::set_dimension(), scale_atmos_dyn_tstep_short_fvm_hivi::solve_bicgstab(), 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_finalize(), 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 366 of file scale_prc.F90.

366  implicit none
367 
368  integer :: ierr
369  logical :: sign_status
370  logical :: sign_exit
371  !---------------------------------------------------------------------------
372 
373  call mpi_barrier(prc_global_comm_world, ierr)
374  if ( prc_intercomm_child /= mpi_comm_null ) &
375  call mpi_comm_free(prc_intercomm_child, ierr)
376  if ( prc_intercomm_parent /= mpi_comm_null ) &
377  call mpi_comm_free(prc_intercomm_parent, ierr)
378 
379  ! FPM polling
380  if ( fpm_alive ) then
381  sign_status = .false.
382  sign_exit = .false.
383  do while ( .NOT. sign_exit )
384  call fpm_polling( sign_status, sign_exit )
385  enddo
386  endif
387 
388  if (prc_universal_handler .NE. mpi_errhandler_null) then
389  call mpi_errhandler_free(prc_universal_handler, ierr)
390  endif
391  if (prc_abort_handler .NE. mpi_errhandler_null) then
392  call mpi_errhandler_free(prc_abort_handler, ierr)
393  endif
394 
395  ! Stop MPI
396  if ( prc_mpi_alive ) then
397  ! free splitted communicator
398  if ( prc_local_comm_world /= prc_global_comm_world ) then
399  call mpi_comm_free(prc_local_comm_world,ierr)
400  endif
401 
402  call mpi_barrier(prc_universal_comm_world,ierr)
403 
404  call mpi_finalize(ierr)
405  endif
406 
407  return

References scale_fpm::fpm_alive, scale_fpm::fpm_polling(), prc_abort_handler, prc_global_comm_world, prc_intercomm_child, prc_intercomm_parent, prc_local_comm_world, prc_mpi_alive, prc_universal_comm_world, and prc_universal_handler.

Referenced by mod_launcher::launcher(), scale::scale_finalize(), 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 419 of file scale_prc.F90.

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

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

Referenced by mod_launcher::launcher().

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 
)

MPI Communicator Split (nesting)

Definition at line 522 of file scale_prc.F90.

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

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

Referenced by mod_launcher::launcher().

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 828 of file scale_prc.F90.

828  implicit none
829 
830  integer :: ierr
831  !---------------------------------------------------------------------------
832 
833  if ( prc_mpi_alive ) then
834  call mpi_barrier(prc_local_comm_world,ierr)
835  endif
836 

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 843 of file scale_prc.F90.

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

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 865 of file scale_prc.F90.

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

References prc_local_comm_world, prc_myrank, and prc_nprocs.

Referenced by scale_prof::prof_rapreport().

Here is the caller graph for this function:

◆ prc_timereorder()

subroutine, public scale_prc::prc_timereorder ( integer, intent(in)  rapnlimit,
integer, intent(inout)  rapnmax,
real(dp), dimension(rapnlimit), intent(inout)  rapttot,
character(len=h_short), dimension(rapnlimit), intent(in)  rapname 
)

reorder rap time

Definition at line 926 of file scale_prc.F90.

926  implicit none
927  integer, intent(in) :: rapnlimit
928  integer, intent(inout) :: rapnmax
929  real(DP), intent(inout) :: rapttot(rapnlimit)
930  character(len=H_SHORT), intent(in) :: rapname(rapnlimit)
931 
932  integer :: rapnmax0
933  character(len=H_SHORT) :: rapname0
934  real(DP) :: rapttot0(rapnlimit)
935 
936  integer :: ierr
937  integer :: i, j
938 
939  if ( prc_ismaster ) rapnmax0 = rapnmax
940  call mpi_bcast( rapnmax0, 1, mpi_integer, prc_masterrank, prc_local_comm_world, ierr )
941 
942  do i = 1, rapnmax0
943  if ( prc_ismaster ) rapname0 = rapname(i)
944  call mpi_bcast( rapname0, h_short, mpi_character, prc_masterrank, prc_local_comm_world, ierr )
945  if ( .not. prc_ismaster ) then
946  do j = 1, rapnmax
947  if ( rapname(j) == rapname0 ) then
948  rapttot0(i) = rapttot(j)
949  exit
950  end if
951  end do
952  if ( j > rapnmax ) then
953  rapttot0(i) = 0.0_dp
954  end if
955  end if
956  end do
957  if ( .not. prc_ismaster ) then
958  rapnmax = rapnmax0
959  do i = 1, rapnmax
960  rapttot(i) = rapttot0(i)
961  end do
962  end if

References scale_fpm::fpm_alive, scale_fpm::fpm_polling(), scale_io::h_short, 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_masterrank, 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 1048 of file scale_prc.F90.

1048  procedure(closer) :: routine
1049 
1050  prc_file_closer => routine
1051 
1052  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 68 of file scale_prc.F90.

68  integer, public, parameter :: PRC_DOMAIN_nlim = 10000

Referenced by mod_launcher::launcher().

◆ prc_comm_null

integer, parameter, public scale_prc::prc_comm_null = MPI_COMM_NULL

Definition at line 69 of file scale_prc.F90.

69  integer, public, parameter :: PRC_COMM_NULL = mpi_comm_null

Referenced by scale_file::file_get_aggregate(), and mod_launcher::launcher().

◆ prc_universal_comm_world

integer, public scale_prc::prc_universal_comm_world = -1

original communicator

Definition at line 72 of file scale_prc.F90.

72  integer, public :: PRC_UNIVERSAL_COMM_WORLD = -1

Referenced by scale_comm_ensemble::comm_ensemble_setup(), 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 73 of file scale_prc.F90.

73  integer, public :: PRC_UNIVERSAL_myrank = -1

Referenced by scale_prc_cartesc::prc_cartesc_setup(), scale_prc_icoa::prc_icoa_setup(), prc_singlecom_setup(), prc_timereorder(), prc_universal_setup(), and scale_random::random_setup().

◆ prc_universal_nprocs

integer, public scale_prc::prc_universal_nprocs = -1

process num in universal communicator

Definition at line 74 of file scale_prc.F90.

74  integer, public :: PRC_UNIVERSAL_nprocs = -1

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

◆ prc_universal_ismaster

logical, public scale_prc::prc_universal_ismaster = .false.

master process in universal communicator?

Definition at line 75 of file scale_prc.F90.

75  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 77 of file scale_prc.F90.

77  integer, public :: PRC_UNIVERSAL_jobID = 0

Referenced by prc_mpisplit_bulk(), and prc_timereorder().

◆ prc_global_comm_world

integer, public scale_prc::prc_global_comm_world = -1

global communicator

Definition at line 80 of file scale_prc.F90.

80  integer, public :: PRC_GLOBAL_COMM_WORLD = -1

Referenced by 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 81 of file scale_prc.F90.

81  integer, public :: PRC_GLOBAL_myrank = -1

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

◆ prc_global_nprocs

integer, public scale_prc::prc_global_nprocs = -1

process num in global communicator

Definition at line 82 of file scale_prc.F90.

82  integer, public :: PRC_GLOBAL_nprocs = -1

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

◆ prc_global_ismaster

logical, public scale_prc::prc_global_ismaster = .false.

master process in global communicator?

Definition at line 83 of file scale_prc.F90.

83  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 85 of file scale_prc.F90.

85  integer, public :: PRC_GLOBAL_domainID = 0

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

◆ prc_global_root

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

root processes in global members

Definition at line 86 of file scale_prc.F90.

86  integer, public :: PRC_GLOBAL_ROOT(PRC_DOMAIN_nlim)

Referenced by mod_launcher::launcher(), and prc_mpisplit_bulk().

◆ 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 91 of file scale_prc.F90.

91  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(), com_gamma(), 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_ensemble::comm_ensemble_setup(), scale_comm_icoa::comm_setup(), mod_copytopo::copytopo(), scale_cpl_phy_sfc_skin::cpl_phy_sfc_skin(), mod_da_driver::da_driver_setup(), 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_query(), scale_file_external_input::file_external_input_regist_file(), 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_finalize(), scale_land_grid_icoa::land_grid_icoa_finalize(), scale_monitor::monitor_finalize(), scale_monitor::monitor_write(), scale_ocean_grid_cartesc::ocean_grid_cartesc_finalize(), 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(), prc_timereorder(), 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_finalize(), scale_urban_grid_icoa::urban_grid_icoa_finalize(), scale_debug::valcheck_1d(), scale_debug::valcheck_2d(), and scale_debug::valcheck_3d().

◆ prc_ismaster

logical, public scale_prc::prc_ismaster = .false.

◆ prc_intercomm_parent

integer, public scale_prc::prc_intercomm_parent = MPI_COMM_NULL

◆ prc_intercomm_child

integer, public scale_prc::prc_intercomm_child = MPI_COMM_NULL

◆ prc_mpi_alive

logical, public scale_prc::prc_mpi_alive = .false.

MPI is alive?

Definition at line 100 of file scale_prc.F90.

100  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_timereorder().

◆ prc_universal_handler

integer, public scale_prc::prc_universal_handler

error handler in universal communicator

Definition at line 101 of file scale_prc.F90.

101  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 102 of file scale_prc.F90.

102  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_singlecom_setup(), and prc_timereorder().

◆ prc_abort_handler

integer, public scale_prc::prc_abort_handler

error handler communicator for aborting

Definition at line 103 of file scale_prc.F90.

103  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:350
scale_fpm::fpm_polling
subroutine, public fpm_polling(run_stat, stop_signal)
Main system of FPM.
Definition: scale_fpm.F90:176
scale_precision::dp
integer, parameter, public dp
Definition: scale_precision.F90:32
scale_fpm::fpm_alive
logical, public fpm_alive
Definition: scale_fpm.F90:32