SCALE-RM
Functions/Subroutines | Variables
scale_rm_process Module Reference

module RM PROCESS More...

Functions/Subroutines

subroutine, public prc_setup
 Setup Processor topology. More...
 

Variables

integer, parameter, public prc_w = 1
 [node direction] west More...
 
integer, parameter, public prc_n = 2
 [node direction] north More...
 
integer, parameter, public prc_e = 3
 [node direction] east More...
 
integer, parameter, public prc_s = 4
 [node direction] south More...
 
integer, parameter, public prc_nw = 5
 [node direction] northwest More...
 
integer, parameter, public prc_ne = 6
 [node direction] northeast More...
 
integer, parameter, public prc_sw = 7
 [node direction] southwest More...
 
integer, parameter, public prc_se = 8
 [node direction] southeast More...
 
integer, public prc_num_x = 1
 x length of 2D processor topology More...
 
integer, public prc_num_y = 1
 y length of 2D processor topology More...
 
integer, dimension(:,:), allocatable, public prc_2drank
 node index in 2D topology More...
 
integer, dimension(8), public prc_next = -1
 node ID of 8 neighbour process More...
 
logical, public prc_has_w
 
logical, public prc_has_n
 
logical, public prc_has_e
 
logical, public prc_has_s
 
logical, public prc_periodic_x = .true.
 periodic condition or not (X)? More...
 
logical, public prc_periodic_y = .true.
 periodic condition or not (Y)? More...
 

Detailed Description

module RM PROCESS

Description
MPI process management module for regional model
Author
Team SCALE
NAMELIST
  • PARAM_PRC
    nametypedefault valuecomment
    PRC_NUM_X integer 1 x length of 2D processor topology
    PRC_NUM_Y integer 1 y length of 2D processor topology
    PRC_PERIODIC_X logical .true. periodic condition or not (X)?
    PRC_PERIODIC_Y logical .true. periodic condition or not (Y)?
    PRC_CART_REORDER logical .false. flag for rank reordering over the cartesian map

History Output
No history output

Function/Subroutine Documentation

◆ prc_setup()

subroutine, public scale_rm_process::prc_setup ( )

Setup Processor topology.

Definition at line 70 of file scale_rm_process.F90.

References scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_fid_nml, scale_stdio::io_l, scale_stdio::io_nml, prc_2drank, scale_process::prc_abort_comm_world, prc_e, scale_process::prc_global_comm_world, scale_process::prc_global_ismaster, scale_process::prc_global_myrank, scale_process::prc_global_nprocs, prc_has_e, prc_has_n, prc_has_s, prc_has_w, scale_process::prc_ismaster, scale_process::prc_local_comm_world, scale_process::prc_masterrank, scale_process::prc_mpi_alive, scale_process::prc_mpistop(), scale_process::prc_myrank, prc_n, prc_ne, prc_next, scale_process::prc_nprocs, prc_num_x, prc_num_y, prc_nw, prc_periodic_x, prc_periodic_y, prc_s, prc_se, prc_sw, scale_process::prc_universal_comm_world, scale_process::prc_universal_ismaster, scale_process::prc_universal_myrank, scale_process::prc_universal_nprocs, and prc_w.

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

70  use scale_process, only: &
71  prc_mpistop, &
72  prc_masterrank, &
73  prc_mpi_alive, &
74  prc_abort_comm_world, &
75  prc_universal_comm_world, &
76  prc_universal_nprocs, &
77  prc_universal_myrank, &
78  prc_universal_ismaster, &
79  prc_global_comm_world, &
80  prc_global_nprocs, &
81  prc_global_myrank, &
82  prc_global_ismaster, &
83  prc_local_comm_world, &
84  prc_nprocs, &
85  prc_myrank, &
86  prc_ismaster
87  implicit none
88 
89  logical :: PRC_CART_REORDER = .false.
90 
91  namelist / param_prc / &
92  prc_num_x, &
93  prc_num_y, &
94  prc_periodic_x, &
95  prc_periodic_y, &
96  prc_cart_reorder
97 
98  logical :: period(2)
99  integer :: divide(2)
100  integer :: coords_W(2)
101  integer :: coords_N(2)
102  integer :: coords_E(2)
103  integer :: coords_S(2)
104  integer :: next_coords(2)
105  integer :: iptbl
106  integer :: next(8)
107 
108  integer :: ierr
109  integer :: p
110  !---------------------------------------------------------------------------
111 
112  if( io_l ) write(io_fid_log,*)
113  if( io_l ) write(io_fid_log,*) '++++++ Module[PROCESS] / Categ[ATMOS-RM COMM] / Origin[SCALElib]'
114 
115  if ( io_l ) then
116  write(io_fid_log,*) '++++++ Start MPI'
117  write(io_fid_log,*)
118  write(io_fid_log,*) '*** Process information ***'
119  write(io_fid_log,'(1x,A,I12)') '*** UNIVERSAL_COMM_WORLD : ', prc_universal_comm_world
120  write(io_fid_log,'(1x,A,I12)') '*** total process [UNIVERSAL] : ', prc_universal_nprocs
121  write(io_fid_log,'(1x,A,I12)') '*** my process ID [UNIVERSAL] : ', prc_universal_myrank
122  write(io_fid_log,'(1x,A,L12)') '*** master rank? [UNIVERSAL] : ', prc_universal_ismaster
123  write(io_fid_log,'(1x,A,I12)') '*** GLOBAL_COMM_WORLD : ', prc_global_comm_world
124  write(io_fid_log,'(1x,A,I12)') '*** total process [GLOBAL] : ', prc_global_nprocs
125  write(io_fid_log,'(1x,A,I12)') '*** my process ID [GLOBAL] : ', prc_global_myrank
126  write(io_fid_log,'(1x,A,L12)') '*** master rank? [GLOBAL] : ', prc_global_ismaster
127  write(io_fid_log,'(1x,A,I12)') '*** LOCAL_COMM_WORLD : ', prc_local_comm_world
128  write(io_fid_log,'(1x,A,I12)') '*** total process [LOCAL] : ', prc_nprocs
129  write(io_fid_log,'(1x,A,I12)') '*** my process ID [LOCAL] : ', prc_myrank
130  write(io_fid_log,'(1x,A,L12)') '*** master rank? [LOCAL] : ', prc_ismaster
131  write(io_fid_log,'(1x,A,I12)') '*** ABORT_COMM_WORLD : ', prc_abort_comm_world
132  write(io_fid_log,'(1x,A,I12)') '*** master rank ID [each world] : ', prc_masterrank
133  endif
134 
135  !--- read namelist
136  rewind(io_fid_conf)
137  read(io_fid_conf,nml=param_prc,iostat=ierr)
138  if( ierr < 0 ) then !--- missing
139  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
140  elseif( ierr > 0 ) then !--- fatal error
141  write(*,*) 'xxx Not appropriate names in namelist PARAM_PRC. Check!'
142  call prc_mpistop
143  endif
144  if( io_nml ) write(io_fid_nml,nml=param_prc)
145 
146  if( io_l ) write(io_fid_log,*)
147  if( io_l ) write(io_fid_log,*) '*** Process allocation ***'
148  if( io_l ) write(io_fid_log,*) '*** No. of Node :', prc_num_x," x ",prc_num_y
149 
150  if ( prc_num_x*prc_num_y /= prc_nprocs ) then
151  write(*,*) 'xxx total number of node does not match that requested. Check!'
152  call prc_mpistop
153  endif
154 
155  if ( mod(prc_nprocs,prc_num_x) /= 0 ) then
156  write(*,*) 'xxx number of requested node cannot devide to 2D. Check!'
157  call prc_mpistop
158  endif
159 
160  ! set communication topology
161  allocate( prc_2drank(-1:prc_nprocs-1,2) )
162  prc_2drank(:,:) = -1
163 
164  do p = 0, prc_nprocs-1
165  prc_2drank(p,1) = mod(p,prc_num_x)
166  prc_2drank(p,2) = (p-prc_2drank(p,1)) / prc_num_x
167  enddo
168 
169  divide(1) = prc_num_y
170  divide(2) = prc_num_x
171  period(1) = prc_periodic_y
172  period(2) = prc_periodic_x
173  if ( prc_mpi_alive ) then
174  call mpi_cart_create(prc_local_comm_world,2,divide,period,prc_cart_reorder,iptbl,ierr)
175  call mpi_cart_shift(iptbl,0,1,prc_next(prc_s),prc_next(prc_n),ierr) ! next rank search Down/Up
176  call mpi_cart_shift(iptbl,1,1,prc_next(prc_w),prc_next(prc_e),ierr) ! next rank search Left/Right
177 
178  ! get neighbor_coordinates
179  prc_has_w = prc_next(prc_w) /= mpi_proc_null
180  if( prc_has_w ) call mpi_cart_coords(iptbl,prc_next(prc_w),2,coords_w,ierr)
181  prc_has_n = prc_next(prc_n) /= mpi_proc_null
182  if( prc_has_n ) call mpi_cart_coords(iptbl,prc_next(prc_n),2,coords_n,ierr)
183  prc_has_e = prc_next(prc_e) /= mpi_proc_null
184  if( prc_has_e ) call mpi_cart_coords(iptbl,prc_next(prc_e),2,coords_e,ierr)
185  prc_has_s = prc_next(prc_s) /= mpi_proc_null
186  if( prc_has_s ) call mpi_cart_coords(iptbl,prc_next(prc_s),2,coords_s,ierr)
187  ! next rank search NorthWest
188  if ( .NOT. prc_has_n &
189  .OR. .NOT. prc_has_w ) then
190  prc_next(prc_nw) = mpi_proc_null
191  else
192  next_coords(1) = coords_n(1)
193  next_coords(2) = coords_w(2)
194  call mpi_cart_rank(iptbl, next_coords, prc_next(prc_nw), ierr)
195  endif
196  ! next rank search NorthEast
197  if ( .NOT. prc_has_n &
198  .OR. .NOT. prc_has_e ) then
199  prc_next(prc_ne) = mpi_proc_null
200  else
201  next_coords(1) = coords_n(1)
202  next_coords(2) = coords_e(2)
203  call mpi_cart_rank(iptbl, next_coords, prc_next(prc_ne), ierr)
204  endif
205  ! next rank search SouthWest
206  if ( .NOT. prc_has_s &
207  .OR. .NOT. prc_has_w ) then
208  prc_next(prc_sw) = mpi_proc_null
209  else
210  next_coords(1) = coords_s(1)
211  next_coords(2) = coords_w(2)
212  call mpi_cart_rank(iptbl, next_coords, prc_next(prc_sw), ierr)
213  endif
214  ! next rank search SouthEast
215  if ( .NOT. prc_has_s &
216  .OR. .NOT. prc_has_e ) then
217  prc_next(prc_se) = mpi_proc_null
218  else
219  next_coords(1) = coords_s(1)
220  next_coords(2) = coords_e(2)
221  call mpi_cart_rank(iptbl, next_coords, prc_next(prc_se), ierr)
222  endif
223  call mpi_comm_free(iptbl,ierr)
224  endif
225 
226  next(:) = max(prc_next(:),-1) ! avoid if MPI_PROC_NULL < -1
227 
228  if( io_l ) write(io_fid_log,*) '*** Node topology :'
229  if( io_l ) write(io_fid_log,'(1x,A,I5,A,I5,A,I5,A,A,I5,A,I5,A,I5,A,A,I5,A,I5,A,I5,A)') &
230  '*** NW(',next(prc_nw),',',prc_2drank(next(prc_nw),1),',',prc_2drank(next(prc_nw),2),')', &
231  ' - N(',next(prc_n) ,',',prc_2drank(next(prc_n) ,1),',',prc_2drank(next(prc_n) ,2),')', &
232  ' - NE(',next(prc_ne),',',prc_2drank(next(prc_ne),1),',',prc_2drank(next(prc_ne),2),')'
233  if( io_l ) write(io_fid_log,'(1x,A)') '*** | | |'
234  if( io_l ) write(io_fid_log,'(1x,A,I5,A,I5,A,I5,A,A,I5,A,I5,A,I5,A,A,I5,A,I5,A,I5,A)') &
235  '*** W(',next(prc_w),',',prc_2drank(next(prc_w),1),',',prc_2drank(next(prc_w),2),')', &
236  ' - P(',prc_myrank ,',',prc_2drank(prc_myrank, 1),',',prc_2drank(prc_myrank, 2),')', &
237  ' - E(',next(prc_e),',',prc_2drank(next(prc_e),1),',',prc_2drank(next(prc_e),2),')'
238  if( io_l ) write(io_fid_log,'(1x,A)') '*** | | |'
239  if( io_l ) write(io_fid_log,'(1x,A,I5,A,I5,A,I5,A,A,I5,A,I5,A,I5,A,A,I5,A,I5,A,I5,A)') &
240  '*** SW(',next(prc_sw),',',prc_2drank(next(prc_sw),1),',',prc_2drank(next(prc_sw),2),')', &
241  ' - S(',next(prc_s) ,',',prc_2drank(next(prc_s) ,1),',',prc_2drank(next(prc_s) ,2),')', &
242  ' - SE(',next(prc_se),',',prc_2drank(next(prc_se),1),',',prc_2drank(next(prc_se),2),')'
243 
244  return
integer, public prc_num_x
x length of 2D processor topology
subroutine, public prc_mpistop
Abort MPI.
integer, public prc_num_y
y length of 2D processor topology
module PROCESS
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ prc_w

integer, parameter, public scale_rm_process::prc_w = 1

[node direction] west

Definition at line 34 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), scale_comm::vars_2d_mpi(), and scale_comm::vars_init_mpi_pc().

34  integer, public, parameter :: PRC_W = 1

◆ prc_n

integer, parameter, public scale_rm_process::prc_n = 2

[node direction] north

Definition at line 35 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), scale_comm::vars_2d_mpi(), and scale_comm::vars_init_mpi_pc().

35  integer, public, parameter :: PRC_N = 2

◆ prc_e

integer, parameter, public scale_rm_process::prc_e = 3

[node direction] east

Definition at line 36 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), scale_comm::vars_2d_mpi(), and scale_comm::vars_init_mpi_pc().

36  integer, public, parameter :: PRC_E = 3

◆ prc_s

integer, parameter, public scale_rm_process::prc_s = 4

[node direction] south

Definition at line 37 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), scale_comm::vars_2d_mpi(), and scale_comm::vars_init_mpi_pc().

37  integer, public, parameter :: PRC_S = 4

◆ prc_nw

integer, parameter, public scale_rm_process::prc_nw = 5

[node direction] northwest

Definition at line 38 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), and scale_comm::vars_init_mpi_pc().

38  integer, public, parameter :: PRC_NW = 5

◆ prc_ne

integer, parameter, public scale_rm_process::prc_ne = 6

[node direction] northeast

Definition at line 39 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), and scale_comm::vars_init_mpi_pc().

39  integer, public, parameter :: PRC_NE = 6

◆ prc_sw

integer, parameter, public scale_rm_process::prc_sw = 7

[node direction] southwest

Definition at line 40 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), and scale_comm::vars_init_mpi_pc().

40  integer, public, parameter :: PRC_SW = 7

◆ prc_se

integer, parameter, public scale_rm_process::prc_se = 8

[node direction] southeast

Definition at line 41 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), and scale_comm::vars_init_mpi_pc().

41  integer, public, parameter :: PRC_SE = 8

◆ prc_num_x

integer, public scale_rm_process::prc_num_x = 1

◆ prc_num_y

integer, public scale_rm_process::prc_num_y = 1

◆ prc_2drank

integer, dimension(:,:), allocatable, public scale_rm_process::prc_2drank

◆ prc_next

integer, dimension(8), public scale_rm_process::prc_next = -1

node ID of 8 neighbour process

Definition at line 47 of file scale_rm_process.F90.

Referenced by prc_setup(), scale_comm::vars8_2d_mpi(), scale_comm::vars8_3d_mpi(), scale_comm::vars_2d_mpi(), and scale_comm::vars_init_mpi_pc().

47  integer, public :: PRC_next(8) = -1

◆ prc_has_w

logical, public scale_rm_process::prc_has_w

◆ prc_has_n

logical, public scale_rm_process::prc_has_n

◆ prc_has_e

logical, public scale_rm_process::prc_has_e

◆ prc_has_s

logical, public scale_rm_process::prc_has_s

◆ prc_periodic_x

logical, public scale_rm_process::prc_periodic_x = .true.

periodic condition or not (X)?

Definition at line 54 of file scale_rm_process.F90.

Referenced by scale_fileio::fileio_create(), scale_grid_index::grid_index_setup(), and prc_setup().

54  logical, public :: PRC_PERIODIC_X = .true.

◆ prc_periodic_y

logical, public scale_rm_process::prc_periodic_y = .true.

periodic condition or not (Y)?

Definition at line 55 of file scale_rm_process.F90.

Referenced by scale_fileio::fileio_create(), scale_grid_index::grid_index_setup(), and prc_setup().

55  logical, public :: PRC_PERIODIC_Y = .true.