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, 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 (ORG_COMM, NUM_DOMAIN, PRC_DOMAINS, CONF_FILES, LOG_SPLIT, bulk_split, color_reorder, INTRA_COMM, inter_parent, inter_child, fname_local)
 MPI Communicator Split. 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+1), 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 121 of file scale_prc.F90.

References prc_abort_comm_world, and prc_mpi_alive.

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

121  implicit none
122 
123  integer, intent(out) :: comm ! communicator
124 
125  integer :: ierr
126  !---------------------------------------------------------------------------
127 
128  call mpi_init(ierr)
129 
130  prc_mpi_alive = .true.
131 ! PRC_UNIVERSAL_handler = MPI_ERRHANDLER_NULL
132 ! call MPI_COMM_CREATE_ERRHANDLER( PRC_MPI_errorhandler, PRC_UNIVERSAL_handler, ierr )
133 
134  comm = mpi_comm_world
135  prc_abort_comm_world = comm
136 
137  return
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,
logical, intent(out)  ismaster 
)

setup MPI in universal communicator

Definition at line 146 of file scale_prc.F90.

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

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

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

184  implicit none
185 
186  logical, intent(in) :: abortall ! abort all jobs?
187  integer, intent(in) :: comm ! communicator
188 
189  integer :: ierr
190  !---------------------------------------------------------------------------
191 
192  prc_global_comm_world = comm
193 
194  call mpi_comm_size(prc_global_comm_world,prc_global_nprocs,ierr)
195  call mpi_comm_rank(prc_global_comm_world,prc_global_myrank,ierr)
196 
197  if ( prc_global_myrank == prc_masterrank ) then
198  prc_global_ismaster = .true.
199  else
200  prc_global_ismaster = .false.
201  endif
202 
203 ! if ( .NOT. abortall ) then
204 ! PRC_ABORT_COMM_WORLD = PRC_GLOBAL_COMM_WORLD
205 !
206 ! call MPI_COMM_SET_ERRHANDLER(PRC_ABORT_COMM_WORLD,PRC_UNIVERSAL_handler,ierr)
207 ! call MPI_COMM_GET_ERRHANDLER(PRC_ABORT_COMM_WORLD,PRC_ABORT_handler ,ierr)
208 ! endif
209 
210  return
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 219 of file scale_prc.F90.

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

219  implicit none
220 
221  integer, intent(in) :: comm ! communicator
222  integer, intent(out) :: myrank ! myrank in this communicator
223  logical, intent(out) :: ismaster ! master process in this communicator?
224 
225  integer :: ierr
226  !---------------------------------------------------------------------------
227 
228  prc_local_comm_world = comm
229 
230  call mpi_comm_rank(prc_local_comm_world,prc_myrank,ierr)
231  call mpi_comm_size(prc_local_comm_world,prc_nprocs,ierr)
232 
233  if ( prc_myrank == prc_masterrank ) then
234  prc_ismaster = .true.
235  else
236  prc_ismaster = .false.
237  endif
238 
239  myrank = prc_myrank
240  ismaster = prc_ismaster
241 
242  return
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 252 of file scale_prc.F90.

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

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

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

301  implicit none
302 
303  logical, intent(in) :: use_fpm ! fpm switch
304  logical, intent(in) :: master ! master flag
305 
306  integer :: ierr
307  !---------------------------------------------------------------------------
308 
309  call mpi_comm_create_errhandler(prc_mpi_errorhandler,prc_universal_handler,ierr)
310 
311  call mpi_comm_set_errhandler(prc_abort_comm_world,prc_universal_handler,ierr)
312  call mpi_comm_get_errhandler(prc_abort_comm_world,prc_abort_handler ,ierr)
313 
314  if ( prc_universal_handler /= prc_abort_handler ) then
315  if( prc_universal_ismaster ) write(*,*) ""
316  if( prc_universal_ismaster ) write(*,*) "ERROR: MPI HANDLER is INCONSISTENT"
317  if( prc_universal_ismaster ) write(*,*) " PRC_UNIVERSAL_handler = ", prc_universal_handler
318  if( prc_universal_ismaster ) write(*,*) " PRC_ABORT_handler = ", prc_abort_handler
319  call prc_abort
320  endif
321 
322  if ( use_fpm ) then
323  call sigvars_get_all( master )
324  call signal( sigint, prc_abort )
325  call signal( sigquit, prc_abort )
326  call signal( sigabrt, prc_abort )
327  call signal( sigfpe, prc_abort )
328  call signal( sigsegv, prc_abort )
329  call signal( sigterm, prc_abort )
330  endif
331 
332  return
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 338 of file scale_prc.F90.

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_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_common::atmos_dyn_filter_setup(), scale_atmos_dyn_fvm_flux::atmos_dyn_fvm_flux_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_rk3::atmos_dyn_tinteg_short_rk3_setup(), scale_atmos_dyn_tinteg_short_rk4::atmos_dyn_tinteg_short_rk4_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_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(), 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_real::atmos_grid_cartesc_real_update_z(), 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(), mod_atmos_phy_ae_vars::atmos_phy_ae_vars_setup(), mod_atmos_phy_bl_driver::atmos_phy_bl_driver_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_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_qhyd2qtrc(), scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14_setup(), scale_atmos_phy_mp_suzuki10::atmos_phy_mp_suzuki10_qhyd2qtrc(), 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_get_diagnostic_1d(), mod_atmos_vars::atmos_vars_get_diagnostic_2d(), mod_atmos_vars::atmos_vars_get_diagnostic_3d(), mod_atmos_vars::atmos_vars_monitor(), 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(), mod_cnv2d::cnv2d_setup(), mod_cnvlanduse::cnvlanduse(), mod_cnvlanduse::cnvlanduse_setup(), mod_cnvtopo::cnvtopo(), mod_cnvtopo::cnvtopo_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(), 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_1d(), 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_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_query_id(), 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_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_set_option(), scale_file_tiledata::file_tiledata_get_data_int1(), scale_file_tiledata::file_tiledata_get_data_real(), scale_file_tiledata::file_tiledata_get_info(), 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_factor2d(), scale_interp::interp_factor3d(), scale_interp::interp_search_horiz_struct(), 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_restart_read(), mod_land_vars::land_vars_setup(), mod_land_vars::land_vars_total(), scale_landuse::landuse_fillhalo(), scale_landuse::landuse_setup(), scale_atmos_dyn_tstep_short_fvm_hivi::make_matrix(), scale_mapprojection::mapprojection_get_attributes(), scale_mapprojection::mapprojection_setup(), 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_kij(), 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(), 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_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(), scale_prof::prof_setup(), scale_random::random_setup(), mod_realinput_grads::read_namelist(), 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_statistics::statistics_setup(), scale_statistics::statistics_total_2d(), scale_statistics::statistics_total_3d(), mod_mkinit::tke_setup(), scale_topography::topo_fillhalo(), scale_topography::topo_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_pc(), and scale_comm_cartesc::vars_init_mpi_pc().

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

◆ prc_mpifinish()

subroutine, public scale_prc::prc_mpifinish ( )

Stop MPI peacefully.

Definition at line 354 of file scale_prc.F90.

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

354  implicit none
355 
356  integer :: ierr
357  logical :: sign_status
358  logical :: sign_exit
359  !---------------------------------------------------------------------------
360 
361  ! FPM polling
362  if ( fpm_alive ) then
363  sign_status = .false.
364  sign_exit = .false.
365  do while ( .NOT. sign_exit )
366  call fpm_polling( sign_status, sign_exit )
367  enddo
368  endif
369 
370  if (prc_universal_handler .NE. mpi_errhandler_null) then
371  call mpi_errhandler_free(prc_universal_handler, ierr)
372  endif
373  if (prc_abort_handler .NE. mpi_errhandler_null) then
374  call mpi_errhandler_free(prc_abort_handler, ierr)
375  endif
376 
377  ! Stop MPI
378  if ( prc_mpi_alive ) then
379  log_newline
380  log_progress(*) 'finalize MPI...'
381 
382  ! free splitted communicator
383  if ( prc_local_comm_world /= prc_global_comm_world ) then
384  call mpi_comm_free(prc_local_comm_world,ierr)
385  endif
386 
387  call mpi_barrier(prc_universal_comm_world,ierr)
388 
389  call mpi_finalize(ierr)
390  log_progress(*) 'MPI is peacefully finalized'
391  endif
392 
393  ! Close logfile, configfile
394  if ( io_l ) then
395  if( io_fid_log /= io_fid_stdout ) close(io_fid_log)
396  endif
397  close(io_fid_conf)
398 
399  return
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
logical, public fpm_alive
Definition: scale_fpm.F90:32
logical, public io_l
output log or not? (this process)
Definition: scale_io.F90:61
integer, parameter, public io_fid_stdout
Definition: scale_io.F90:54
subroutine, public fpm_polling(run_stat, stop_signal)
Main system of FPM.
Definition: scale_fpm.F90:176
integer, public io_fid_log
Log file ID.
Definition: scale_io.F90:56
Here is the call graph for this function:
Here is the caller graph for this function:

◆ prc_mpisplit()

subroutine, public scale_prc::prc_mpisplit ( integer, intent(in)  ORG_COMM,
integer, intent(in)  NUM_DOMAIN,
integer, dimension(:), intent(in)  PRC_DOMAINS,
character(len=*), dimension(:), intent(in)  CONF_FILES,
logical, intent(in)  LOG_SPLIT,
logical, intent(in)  bulk_split,
logical, intent(in)  color_reorder,
integer, intent(out)  INTRA_COMM,
integer, intent(out)  inter_parent,
integer, intent(out)  inter_child,
character(len=h_long), intent(out)  fname_local 
)

MPI Communicator Split.

Definition at line 416 of file scale_prc.F90.

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

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

416  implicit none
417 
418  integer, intent(in) :: org_comm
419  integer, intent(in) :: num_domain
420  integer, intent(in) :: prc_domains(:)
421  character(len=*), intent(in) :: conf_files(:)
422  logical, intent(in) :: log_split
423  logical, intent(in) :: bulk_split
424  logical, intent(in) :: color_reorder
425  integer, intent(out) :: intra_comm
426  integer, intent(out) :: inter_parent
427  integer, intent(out) :: inter_child
428  character(len=H_LONG), intent(out) :: fname_local
429 
430  integer :: parent_col(prc_domain_nlim) ! parent color number
431  integer :: child_col(prc_domain_nlim) ! child color number
432  integer :: prc_root(0:prc_domain_nlim) ! root process in the color
433  integer, allocatable :: color_list(:) ! member list in each color
434  integer, allocatable :: key_list(:) ! local process number in each color
435 
436  integer :: total_nmax
437  integer :: org_myrank ! my rank number in the original communicator
438  integer :: org_nmax ! total rank number in the original communicator
439 
440  logical :: do_create_p(prc_domain_nlim)
441  logical :: do_create_c(prc_domain_nlim)
442  logical :: reordering
443 
444  character(len=H_LONG) :: col_file(0:prc_domain_nlim)
445  character(len=4) :: col_num
446 
447  integer :: i, ii
448  integer :: itag, ierr
449  !---------------------------------------------------------------------------
450 
451  intra_comm = org_comm
452  inter_parent = mpi_comm_null
453  inter_child = mpi_comm_null
454  fname_local = conf_files(1)
455 
456  if ( num_domain > 1 ) then ! multi domain run
457  call mpi_comm_rank(org_comm,org_myrank,ierr)
458  call mpi_comm_size(org_comm,org_nmax, ierr)
459  allocate( color_list(0:org_nmax-1) )
460  allocate( key_list(0:org_nmax-1) )
461 
462  total_nmax = 0
463  do i = 1, num_domain
464  total_nmax = total_nmax + prc_domains(i)
465  enddo
466  if ( total_nmax /= org_nmax ) then
467  if( prc_universal_ismaster ) then
468  log_error("PRC_MPIsplit",*) "MPI PROCESS NUMBER is INCONSISTENT"
469  log_error_cont(*) " REQUESTED NPROCS = ", total_nmax, " LAUNCHED NPROCS = ", org_nmax
470  end if
471  call prc_abort
472  endif
473 
474  reordering = color_reorder
475  if ( bulk_split ) then
476  reordering = .false.
477  endif
478  call prc_mpicoloring( org_comm, & ! [IN]
479  num_domain, & ! [IN]
480  prc_domains, & ! [IN]
481  conf_files, & ! [IN]
482  reordering, & ! [IN]
483  log_split, & ! [IN]
484  color_list, & ! [OUT]
485  prc_root, & ! [OUT]
486  key_list, & ! [OUT]
487  parent_col, & ! [OUT]
488  child_col, & ! [OUT]
489  col_file ) ! [OUT]
490 
491  if ( bulk_split ) then
492  ii = 1
493  do i=0, prc_domain_nlim
494  prc_global_root(ii) = prc_root(i)
495  ii = ii + 1
496  enddo
497  endif
498 
499  ! split comm_world
500  call mpi_comm_split(org_comm, &
501  color_list(org_myrank), &
502  key_list(org_myrank), &
503  intra_comm, ierr)
504  if ( bulk_split ) then
505  write(col_num,'(I4.4)') color_list(org_myrank)
506  fname_local = col_num
507  prc_universal_jobid = color_list(org_myrank)
508  else
509  fname_local = col_file(color_list(org_myrank))
510  endif
511 
512  ! set parent-child relationship
513  do_create_p(:) = .false.
514  do_create_c(:) = .false.
515  if ( .NOT. bulk_split ) then
516  if ( prc_universal_ismaster ) write(*,*)
517  if ( prc_universal_ismaster ) write(*,*) "INFO [PRC_MPIsplit] Inter-domain relationship information"
518  do i = 1, num_domain-1
519  if ( prc_universal_ismaster ) write(*,'(5x,A,I2.2)') "Relationship No. ", i
520  if ( prc_universal_ismaster ) write(*,'(5x,2(A,I2))') "Parent color = ", parent_col(i), &
521  " <=> child color = ", child_col(i)
522  if ( color_list(org_myrank) == parent_col(i) ) then
523  do_create_p(i) = .true.
524  elseif ( color_list(org_myrank) == child_col(i) ) then
525  do_create_c(i) = .true.
526  endif
527  enddo
528  endif
529 
530  ! create inter-commnunicator
531  inter_parent = mpi_comm_null
532  inter_child = mpi_comm_null
533  if ( .NOT. bulk_split ) then
534  do i = 1, num_domain-1
535  itag = i*100
536  if ( do_create_p(i) ) then ! as a parent
537  call mpi_intercomm_create( intra_comm, prc_masterrank, &
538  org_comm, prc_root(child_col(i)), &
539  itag, inter_child, ierr )
540  elseif( do_create_c(i) ) then ! as a child
541  call mpi_intercomm_create( intra_comm, prc_masterrank, &
542  org_comm, prc_root(parent_col(i)), &
543  itag, inter_parent, ierr )
544  endif
545  call mpi_barrier(org_comm, ierr)
546  enddo
547  endif
548 
549  deallocate( color_list, key_list )
550 
551  elseif ( num_domain == 1 ) then ! single domain run
552  ! if ( PRC_UNIVERSAL_IsMaster ) write (*,*) "INFO [PRC_MPIsplit] a single communicator"
553  else
554  if ( prc_universal_ismaster ) then
555  write(*,*)"ERROR [RPC_MPIsplit] REQUESTED DOMAIN NUMBER IS NOT ACCEPTABLE"
556  end if
557  call prc_abort
558  endif
559 
560  return
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
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 798 of file scale_prc.F90.

References prc_local_comm_world, and prc_mpi_alive.

Referenced by scale_prof::prof_rapend(), and scale_prof::prof_rapstart().

798  implicit none
799 
800  integer :: ierr
801  !---------------------------------------------------------------------------
802 
803  if ( prc_mpi_alive ) then
804  call mpi_barrier(prc_local_comm_world,ierr)
805  endif
806 
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 813 of file scale_prc.F90.

References 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().

813  implicit none
814 
815  real(DP) :: time
816  !---------------------------------------------------------------------------
817 
818  if ( prc_mpi_alive ) then
819  time = real(MPI_WTIME(), kind=dp)
820  else
821  call cpu_time(time)
822  endif
823 
integer, parameter, public dp
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 835 of file scale_prc.F90.

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

835  implicit none
836 
837  real(DP), intent(out) :: avgvar(:)
838  real(DP), intent(out) :: maxvar(:)
839  real(DP), intent(out) :: minvar(:)
840  integer, intent(out) :: maxidx(:)
841  integer, intent(out) :: minidx(:)
842  real(DP), intent(in) :: var(:)
843 
844  real(DP), allocatable :: statval(:,:)
845  integer :: vsize
846 
847  real(DP) :: totalvar
848  integer :: ierr
849  integer :: v, p
850  !---------------------------------------------------------------------------
851 
852  vsize = size(var(:))
853 
854  allocate( statval(vsize,0:prc_nprocs-1) )
855  statval(:,:) = 0.0_dp
856 
857  do v = 1, vsize
858  statval(v,prc_myrank) = var(v)
859  enddo
860 
861  ! MPI broadcast
862  do p = 0, prc_nprocs-1
863  call mpi_bcast( statval(1,p), &
864  vsize, &
865  mpi_double_precision, &
866  p, &
867  prc_local_comm_world, &
868  ierr )
869  enddo
870 
871  do v = 1, vsize
872  totalvar = 0.0_dp
873  do p = 0, prc_nprocs-1
874  totalvar = totalvar + statval(v,p)
875  enddo
876  avgvar(v) = totalvar / prc_nprocs
877 
878  maxvar(v) = maxval(statval(v,:))
879  minvar(v) = minval(statval(v,:))
880  maxidx(v:v) = maxloc(statval(v,:))
881  minidx(v:v) = minloc(statval(v,:))
882  enddo
883 
884  deallocate( statval )
885 
886  return
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 972 of file scale_prc.F90.

Referenced by scale_file::file_setup().

972  procedure(closer) :: routine
973 
974  prc_file_closer => routine
975 
976  return
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 66 of file scale_prc.F90.

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

66  integer, public, parameter :: prc_domain_nlim = 10000

◆ prc_comm_null

integer, parameter, public scale_prc::prc_comm_null = MPI_COMM_NULL

Definition at line 67 of file scale_prc.F90.

Referenced by scale_file::file_get_aggregate().

67  integer, public, parameter :: prc_comm_null = mpi_comm_null

◆ prc_universal_comm_world

integer, public scale_prc::prc_universal_comm_world = -1

original communicator

Definition at line 70 of file scale_prc.F90.

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

70  integer, public :: prc_universal_comm_world = -1

◆ prc_universal_myrank

integer, public scale_prc::prc_universal_myrank = -1

myrank in universal communicator

Definition at line 71 of file scale_prc.F90.

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

71  integer, public :: prc_universal_myrank = -1

◆ prc_universal_nprocs

integer, public scale_prc::prc_universal_nprocs = -1

process num in universal communicator

Definition at line 72 of file scale_prc.F90.

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

72  integer, public :: prc_universal_nprocs = -1

◆ prc_universal_ismaster

logical, public scale_prc::prc_universal_ismaster = .false.

master process in universal communicator?

Definition at line 73 of file scale_prc.F90.

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(), prc_singlecom_setup(), and prc_universal_setup().

73  logical, public :: prc_universal_ismaster = .false.

◆ prc_universal_jobid

integer, public scale_prc::prc_universal_jobid = 0

my job ID in universal communicator

Definition at line 75 of file scale_prc.F90.

Referenced by prc_mpisplit(), and prc_mpitimestat().

75  integer, public :: prc_universal_jobid = 0

◆ prc_global_comm_world

integer, public scale_prc::prc_global_comm_world = -1

global communicator

Definition at line 78 of file scale_prc.F90.

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

78  integer, public :: prc_global_comm_world = -1

◆ prc_global_myrank

integer, public scale_prc::prc_global_myrank = -1

myrank in global communicator

Definition at line 79 of file scale_prc.F90.

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

79  integer, public :: prc_global_myrank = -1

◆ prc_global_nprocs

integer, public scale_prc::prc_global_nprocs = -1

process num in global communicator

Definition at line 80 of file scale_prc.F90.

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

80  integer, public :: prc_global_nprocs = -1

◆ prc_global_ismaster

logical, public scale_prc::prc_global_ismaster = .false.

master process in global communicator?

Definition at line 81 of file scale_prc.F90.

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

81  logical, public :: prc_global_ismaster = .false.

◆ prc_global_domainid

integer, public scale_prc::prc_global_domainid = 0

my domain ID in global communicator

Definition at line 83 of file scale_prc.F90.

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

83  integer, public :: prc_global_domainid = 0

◆ prc_global_root

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

root processes in global members

Definition at line 84 of file scale_prc.F90.

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

84  integer, public :: prc_global_root(prc_domain_nlim+1)

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

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(), mod_atmos_vars::atmos_vars_monitor(), 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_kij(), 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_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().

89  integer, public :: prc_myrank = 0

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

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

93  logical, public :: prc_mpi_alive = .false.

◆ prc_universal_handler

integer, public scale_prc::prc_universal_handler

error handler in universal communicator

Definition at line 94 of file scale_prc.F90.

Referenced by prc_errhandler_setup(), and prc_mpifinish().

94  integer, public :: prc_universal_handler

◆ prc_abort_comm_world

integer, public scale_prc::prc_abort_comm_world

communicator for aborting

Definition at line 95 of file scale_prc.F90.

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

95  integer, public :: prc_abort_comm_world

◆ prc_abort_handler

integer, public scale_prc::prc_abort_handler

error handler communicator for aborting

Definition at line 96 of file scale_prc.F90.

Referenced by prc_errhandler_setup(), and prc_mpifinish().

96  integer, public :: prc_abort_handler