SCALE-RM
scale_prc_icoA.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use mpi
17  use scale_precision
18  use scale_io
19  use scale_prc, only: &
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: prc_icoa_setup
29  public :: prc_icoa_rgn_generate
30 
31  !-----------------------------------------------------------------------------
32  !
33  !++ Public parameters & variables
34  !
35 
36  !====== Information for processes ======
37 
38  integer, public, parameter :: prc_rank_pl = prc_masterrank
39  logical, public :: prc_have_pl
40 
41  !====== Information for processes-region relationship ======
42 
43  integer, public, parameter :: i_l = 1
44  integer, public, parameter :: i_prc = 2
45 
46  integer, public, parameter :: i_rgnid = 1
47  integer, public, parameter :: i_dir = 2
48 
49  ! Identifiers of directions of region edges
50  integer, public, parameter :: i_sw = 1
51  integer, public, parameter :: i_nw = 2
52  integer, public, parameter :: i_ne = 3
53  integer, public, parameter :: i_se = 4
54 
55  ! Identifiers of directions of region vertices
56  integer, public, parameter :: i_w = 1
57  integer, public, parameter :: i_n = 2
58  integer, public, parameter :: i_e = 3
59  integer, public, parameter :: i_s = 4
60 
61  ! Identifier of poles (north pole or south pole)
62  integer, public, parameter :: i_npl = 1
63  integer, public, parameter :: i_spl = 2
64 
65  ! main parameter
66  integer, public :: prc_rgn_level = -1
67  integer, public :: prc_rgn_dmd = 10
68  integer, public :: prc_rgn_vlink = 5
69 
70  ! region
71  integer, public :: prc_rgn_total
72  integer, public :: prc_rgn_local
73  integer, public, parameter :: prc_rgn_total_pl = 2
74 
75  integer, public, parameter :: prc_rgn_local_lim = 2560
76 
77  integer, public, allocatable :: prc_rgn_edge_tab (:,:,:)
78 
79  integer, public, allocatable :: prc_rgn_vert_num (:,:)
80  integer, public, allocatable :: prc_rgn_vert_tab (:,:,:,:)
81  integer, public, allocatable :: prc_rgn_vert_tab_pl(:,:,:)
82 
83  integer, public, allocatable :: prc_rgn_lp2r(:,:)
84  integer, public, allocatable :: prc_rgn_r2lp(:,:)
85  integer, public, allocatable :: prc_rgn_l2r (:)
86 
87  integer, public :: prc_rgn_r2p_pl(prc_rgn_total_pl)
88  integer, public :: prc_rgn_rgn4pl(prc_rgn_total_pl)
89 
90  !====== Information for regions ======
91 
92  logical, public, allocatable :: prc_rgn_have_sgp(:)
93 
94  !-----------------------------------------------------------------------------
95  !
96  !++ Private procedure
97  !
98  private :: prc_icoa_rgn_setup
99  private :: prc_icoa_rgn_input
100  private :: prc_icoa_rgn_output
101  private :: prc_icoa_rgn_vertex_walkaround
102  private :: output_info
103 
104  !-----------------------------------------------------------------------------
105  !
106  !++ Private parameters & variables
107  !
108  character(len=H_LONG), private :: prc_icoa_rgn_in_fname = ''
109  character(len=H_LONG), private :: prc_icoa_rgn_out_fname = ''
110 
111  character(len=2), private, parameter :: prc_rgn_edgename(4) = (/'SW','NW','NE','SE'/)
112  character(len=2), private, parameter :: prc_rgn_vertname(4) = (/'W ','N ','E ','S '/)
113 
114  logical, private :: debug = .false.
115 
116  !-----------------------------------------------------------------------------
117 contains
118  !-----------------------------------------------------------------------------
120  subroutine prc_icoa_setup
121  use scale_prc, only: &
122  prc_abort, &
133  prc_nprocs, &
134  prc_myrank, &
136  implicit none
137 
138  namelist / param_prc_icoa / &
139  prc_rgn_level, &
140  debug
141 
142  integer :: l, rgnid
143 
144  integer :: ierr
145  !---------------------------------------------------------------------------
146 
147  log_newline
148  log_info("PRC_ICOA_setup",*) 'Setup'
149 
150  if ( io_l ) then
151  log_newline
152  log_progress(*) 'start MPI'
153  log_newline
154  log_info("PRC_ICOA_setup",*) 'Process information '
155  log_info_cont('(1x,A,I12)') 'UNIVERSAL_COMM_WORLD : ', prc_universal_comm_world
156  log_info_cont('(1x,A,I12)') 'total process [UNIVERSAL] : ', prc_universal_nprocs
157  log_info_cont('(1x,A,I12)') 'my process ID [UNIVERSAL] : ', prc_universal_myrank
158  log_info_cont('(1x,A,L12)') 'master rank? [UNIVERSAL] : ', prc_universal_ismaster
159  log_info_cont('(1x,A,I12)') 'GLOBAL_COMM_WORLD : ', prc_global_comm_world
160  log_info_cont('(1x,A,I12)') 'total process [GLOBAL] : ', prc_global_nprocs
161  log_info_cont('(1x,A,I12)') 'my process ID [GLOBAL] : ', prc_global_myrank
162  log_info_cont('(1x,A,L12)') 'master rank? [GLOBAL] : ', prc_global_ismaster
163  log_info_cont('(1x,A,I12)') 'LOCAL_COMM_WORLD : ', prc_local_comm_world
164  log_info_cont('(1x,A,I12)') 'total process [LOCAL] : ', prc_nprocs
165  log_info_cont('(1x,A,I12)') 'my process ID [LOCAL] : ', prc_myrank
166  log_info_cont('(1x,A,L12)') 'master rank? [LOCAL] : ', prc_ismaster
167  log_info_cont('(1x,A,I12)') 'ABORT_COMM_WORLD : ', prc_abort_comm_world
168  log_info_cont('(1x,A,I12)') 'master rank ID [each world] : ', prc_masterrank
169  endif
170 
171  !--- read namelist
172  rewind(io_fid_conf)
173  read(io_fid_conf,nml=param_prc_icoa,iostat=ierr)
174  if( ierr < 0 ) then !--- missing
175  log_info("PRC_ICOA_setup",*) 'Not found namelist. Default used.'
176  elseif( ierr > 0 ) then !--- fatal error
177  log_error("PRC_ICOA_setup",*) 'Not appropriate names in namelist PARAM_PRC_ICOA. Check!'
178  call prc_abort
179  endif
180  if( io_nml ) write(io_fid_nml,nml=param_prc_icoa)
181 
182  if ( prc_rgn_level < 0 ) then
183  log_error("PRC_ICOA_setup",*) 'PRC_RGN_level is not appropriate :', prc_rgn_level
184  call prc_abort
185  endif
186 
189 
190  if ( mod(prc_rgn_total,prc_nprocs) /= 0 ) then
191  log_error("PRC_ICOA_setup",*) 'Number of total region must be divisible by the number of process', prc_rgn_total, prc_nprocs
192  call prc_abort
193  endif
194 
195  if ( prc_rgn_local > prc_rgn_local_lim ) then
196  log_error("PRC_ICOA_setup",*) 'limit exceed! local region: ', prc_rgn_local, prc_rgn_local_lim
197  call prc_abort
198  endif
199 
200  call prc_icoa_rgn_setup
201 
202  ! pole region management flag
203  if ( prc_myrank == prc_rank_pl ) then
204  prc_have_pl = .true.
205  else
206  prc_have_pl = .false.
207  endif
208 
209  ! singlar point management flag
210  allocate( prc_rgn_have_sgp(prc_rgn_local) )
211  prc_rgn_have_sgp(:) = .false.
212 
213  do l = 1, prc_rgn_local
214  rgnid = prc_rgn_lp2r(l,prc_myrank)
215  if ( prc_rgn_vert_num(i_w,rgnid) == 3 ) then
216  prc_rgn_have_sgp(l) = .true.
217  endif
218  enddo
219 
220  call output_info
221 
222  return
223  end subroutine prc_icoa_setup
224 
225  !-----------------------------------------------------------------------------
226  subroutine prc_icoa_rgn_setup
227  use scale_prc, only: &
228  prc_abort, &
229  prc_nprocs, &
230  prc_myrank
231  implicit none
232 
233  namelist / param_prc_icoa_rgn / &
234  prc_icoa_rgn_in_fname, &
235  prc_icoa_rgn_out_fname, &
236  debug
237 
238  integer :: l, p, r
239  integer :: ierr
240  !---------------------------------------------------------------------------
241 
242  log_newline
243  log_info("PRC_ICOA_RGN_setup",*) 'Setup'
244 
245  !--- read namelist
246  rewind(io_fid_conf)
247  read(io_fid_conf,nml=param_prc_icoa_rgn,iostat=ierr)
248  if( ierr < 0 ) then !--- missing
249  log_info("PRC_ICOA_RGN_setup",*) 'Not found namelist. Default used.'
250  elseif( ierr > 0 ) then !--- fatal error
251  log_error("PRC_ICOA_RGN_setup",*) 'Not appropriate names in namelist PARAM_PRC_ICOA_RGN. Check!'
252  call prc_abort
253  endif
254  if( io_nml ) write(io_fid_nml,nml=param_prc_icoa_rgn)
255 
256  ! Global information (Each process has all the information)
258  allocate( prc_rgn_lp2r(prc_rgn_local,0:prc_nprocs-1) )
259 
260  if ( prc_icoa_rgn_in_fname /= '' ) then
261  call prc_icoa_rgn_input( prc_icoa_rgn_in_fname, & ! [IN]
262  prc_nprocs, & ! [IN]
263  prc_rgn_total, & ! [IN]
264  prc_rgn_local, & ! [IN]
265  prc_rgn_edge_tab(:,:,:), & ! [OUT]
266  prc_rgn_lp2r(:,:) ) ! [OUT]
267  else
268  log_info("PRC_ICOA_RGN_setup",*) 'input file is not specified.'
269 
270  call prc_icoa_rgn_generate( prc_rgn_level, & ! [IN]
271  prc_nprocs, & ! [IN]
272  prc_rgn_total, & ! [IN]
273  prc_rgn_local, & ! [IN]
274  prc_rgn_edge_tab(:,:,:), & ! [OUT]
275  prc_rgn_lp2r(:,:) ) ! [OUT]
276  endif
277 
278  if ( prc_icoa_rgn_out_fname /= '' ) then
279  call prc_icoa_rgn_output( prc_icoa_rgn_out_fname, & ! [IN]
280  prc_nprocs, & ! [IN]
281  prc_rgn_total, & ! [IN]
282  prc_rgn_local, & ! [IN]
283  prc_rgn_edge_tab(:,:,:), & ! [IN]
284  prc_rgn_lp2r(:,:) ) ! [IN]
285  endif
286 
287  !--- additional table (reversal,local)
288  allocate( prc_rgn_r2lp(2,prc_rgn_total) )
289  allocate( prc_rgn_l2r( prc_rgn_local) )
290 
291  do p = 0, prc_nprocs-1
292  do l = 1, prc_rgn_local
293  prc_rgn_r2lp(i_l, prc_rgn_lp2r(l,p)) = l
295  enddo
296  enddo
297 
298  do l = 1, prc_rgn_local
300  enddo
301 
302  !--- region connection chains around the diamond vertexes
306 
307  call prc_icoa_rgn_vertex_walkaround( prc_rgn_total, & ! [IN]
308  prc_rgn_total_pl, & ! [IN]
309  prc_rgn_vlink, & ! [IN]
310  prc_rgn_edge_tab(:,:,:), & ! [IN]
311  prc_rgn_vert_num(:,:), & ! [OUT]
312  prc_rgn_vert_tab(:,:,:,:), & ! [OUT]
313  prc_rgn_vert_tab_pl(:,:,:) ) ! [OUT]
314 
315  !--- tables for pole
316  do r = 1, prc_rgn_total
317  if ( prc_rgn_vert_num(i_n,r) == prc_rgn_vlink ) then
318  prc_rgn_rgn4pl(i_npl) = r
319  exit
320  endif
321  enddo
322 
323  do r = 1, prc_rgn_total
324  if ( prc_rgn_vert_num(i_s,r) == prc_rgn_vlink ) then
325  prc_rgn_rgn4pl(i_spl) = r
326  exit
327  endif
328  enddo
329 
332 
333  return
334  end subroutine prc_icoa_rgn_setup
335 
336  !-----------------------------------------------------------------------------
338  subroutine prc_icoa_rgn_input( &
339  in_fname, &
340  pall, &
341  rall, &
342  lall, &
343  edge_tab, &
344  lp2r )
345  use scale_prc, only: &
346  prc_abort
347  implicit none
348 
349  character(len=*), intent(in) :: in_fname
350  integer, intent(in) :: pall
351  integer, intent(in) :: rall
352  integer, intent(in) :: lall
353  integer, intent(out) :: edge_tab(2,4,rall)
354  integer, intent(out) :: lp2r(lall,0:pall-1)
355 
356  integer :: num_of_rgn
357 
358  namelist / rgn_info / &
359  num_of_rgn
360 
361  integer :: rgnid
362  integer :: sw(2) = -1
363  integer :: nw(2) = -1
364  integer :: ne(2) = -1
365  integer :: se(2) = -1
366 
367  namelist / rgn_link_info / &
368  rgnid, &
369  sw, &
370  nw, &
371  ne, &
372  se
373 
374  integer :: num_of_proc
375 
376  namelist / proc_info / &
377  num_of_proc
378 
379  integer :: peid
380  integer :: num_of_mng
381  integer :: mng_rgnid(prc_rgn_local_lim)
382 
383  namelist / rgn_mng_info / &
384  peid, &
385  num_of_mng, &
386  mng_rgnid
387 
388  integer :: fid, ierr
389  integer :: r, p, l
390  !---------------------------------------------------------------------------
391 
392  log_info("PRC_ICOA_RGN_input",*) 'input region management information file: ', trim(in_fname)
393 
394  fid = io_get_available_fid()
395  open( unit = fid, &
396  file = trim(in_fname), &
397  form = 'formatted', &
398  status = 'old', &
399  iostat = ierr )
400 
401  ! ERROR if filename are not defined
402  if ( ierr /= 0 ) then
403  log_error("PRC_ICOA_RGN_input",*) 'File is not found!', trim(in_fname)
404  call prc_abort
405  endif
406 
407  read(fid,nml=rgn_info)
408 
409  if ( num_of_rgn /= rall ) then
410  log_error("PRC_ICOA_RGN_input",*) 'Missmatch of region number!'
411  log_error_cont(*) 'rall= ', rall,', num_of_rgn=', num_of_rgn
412  call prc_abort
413  endif
414 
415  do r = 1, rall
416  read(fid,nml=rgn_link_info)
417 
418  edge_tab(i_rgnid:i_dir,i_sw,rgnid) = sw(i_rgnid:i_dir)
419  edge_tab(i_rgnid:i_dir,i_nw,rgnid) = nw(i_rgnid:i_dir)
420  edge_tab(i_rgnid:i_dir,i_ne,rgnid) = ne(i_rgnid:i_dir)
421  edge_tab(i_rgnid:i_dir,i_se,rgnid) = se(i_rgnid:i_dir)
422  enddo
423 
424  read(fid,nml=proc_info)
425  if ( num_of_proc /= pall ) then
426  log_error("PRC_ICOA_RGN_input",*) ' Missmatch of process number!'
427  log_error_cont(*) ' pall= ', pall, ', num_of_proc=', num_of_proc
428  call prc_abort
429  endif
430 
431  do p = 0, pall-1
432  read(fid,nml=rgn_mng_info)
433 
434  if ( p /= peid ) then
435  log_error("PRC_ICOA_RGN_input",*) 'Wrong peid: ', p, peid
436  call prc_abort
437  endif
438 
439  if ( num_of_mng /= lall ) then
440  log_error("PRC_ICOA_RGN_input",*) 'number of local region is not match: ', p, num_of_mng, lall
441  call prc_abort
442  endif
443 
444  lp2r(:,p) = -1 ! initialize
445  do l = 1, lall
446  lp2r(l,p) = mng_rgnid(l)
447  enddo
448  enddo
449 
450  close(fid)
451 
452  return
453  end subroutine prc_icoa_rgn_input
454 
455  !-----------------------------------------------------------------------------
457  subroutine prc_icoa_rgn_output( &
458  out_fname, &
459  pall, &
460  rall, &
461  lall, &
462  edge_tab, &
463  lp2r )
464  use scale_prc, only: &
466  implicit none
467 
468  character(len=*), intent(in) :: out_fname
469  integer, intent(in) :: pall
470  integer, intent(in) :: rall
471  integer, intent(in) :: lall
472  integer, intent(in) :: edge_tab(2,4,rall)
473  integer, intent(in) :: lp2r(lall,pall)
474 
475  integer :: num_of_rgn
476 
477  namelist / rgn_info / &
478  num_of_rgn
479 
480  integer :: rgnid
481  integer :: sw(2) = -1
482  integer :: nw(2) = -1
483  integer :: ne(2) = -1
484  integer :: se(2) = -1
485 
486  namelist / rgn_link_info / &
487  rgnid, &
488  sw, &
489  nw, &
490  ne, &
491  se
492 
493  integer :: num_of_proc
494 
495  namelist / proc_info / &
496  num_of_proc
497 
498  integer :: peid
499  integer :: num_of_mng
500  integer :: mng_rgnid(prc_rgn_local_lim)
501 
502  namelist / rgn_mng_info / &
503  peid, &
504  num_of_mng, &
505  mng_rgnid
506 
507  integer :: fid
508  integer :: r, p, l
509  !---------------------------------------------------------------------------
510 
511  if ( prc_ismaster ) then
512  log_info("PRC_ICOA_RGN_output",*) 'output region management information file: ', trim(out_fname)
513 
514  fid = io_get_available_fid()
515  open( unit = fid, &
516  file = trim(out_fname), &
517  form = 'formatted' )
518 
519  num_of_rgn = rall
520  write(fid,nml=rgn_info)
521 
522  do r = 1, rall
523  rgnid = r
524  sw(i_rgnid:i_dir) = edge_tab(i_rgnid:i_dir,i_sw,rgnid)
525  nw(i_rgnid:i_dir) = edge_tab(i_rgnid:i_dir,i_nw,rgnid)
526  ne(i_rgnid:i_dir) = edge_tab(i_rgnid:i_dir,i_ne,rgnid)
527  se(i_rgnid:i_dir) = edge_tab(i_rgnid:i_dir,i_se,rgnid)
528 
529  write(fid,nml=rgn_link_info)
530  enddo
531 
532  num_of_proc = pall
533  write(fid,nml=proc_info)
534 
535  do p = 0, pall-1
536  peid = p
537  num_of_mng = lall
538 
539  mng_rgnid(:) = -1
540  do l = 1, lall
541  mng_rgnid(l) = lp2r(l,p)
542  enddo
543 
544  write(fid,nml=rgn_mng_info)
545  enddo
546 
547  close(fid)
548  else
549  log_info("PRC_ICOA_RGN_output",*) 'output region management information file at the master process'
550  endif
551 
552  return
553  end subroutine prc_icoa_rgn_output
554 
555  !-----------------------------------------------------------------------------
557  Subroutine prc_icoa_rgn_generate( &
558  rlevel, &
559  pall, &
560  rall, &
561  lall, &
562  edge_tab, &
563  lp2r )
564  use scale_prc, only: &
565  prc_abort
566  implicit none
567 
568  integer, intent(in) :: rlevel
569  integer, intent(in) :: pall
570  integer, intent(in) :: rall
571  integer, intent(in) :: lall
572  integer, intent(out) :: edge_tab(2,4,rall)
573  integer, intent(out) :: lp2r(lall,0:pall-1)
574 
575  integer, parameter :: nmax_dmd = 10
576  integer :: dmd_data(4,nmax_dmd)
577  integer :: rall_1d, rall_1dmd
578 
579  integer :: d_nb, i_nb, j_nb, rgnid_nb, direction
580  integer :: d, i, j, rgnid
581  integer :: l, p
582  !---------------------------------------------------------------------------
583 
584  log_info("PRC_ICOA_RGN_generate",*) 'generate region management information file'
585  log_info_cont(*) 'Topology: icosahedral'
586 
587  dmd_data(:, 1) = (/ 6, 5, 2,10 /)
588  dmd_data(:, 2) = (/ 10, 1, 3, 9 /)
589  dmd_data(:, 3) = (/ 9, 2, 4, 8 /)
590  dmd_data(:, 4) = (/ 8, 3, 5, 7 /)
591  dmd_data(:, 5) = (/ 7, 4, 1, 6 /)
592  dmd_data(:, 6) = (/ 7, 5, 1,10 /)
593  dmd_data(:, 7) = (/ 8, 4, 5, 6 /)
594  dmd_data(:, 8) = (/ 9, 3, 4, 7 /)
595  dmd_data(:, 9) = (/ 10, 2, 3, 8 /)
596  dmd_data(:,10) = (/ 6, 1, 2, 9 /)
597 
598  rall_1d = 2**rlevel
599  rall_1dmd = rall_1d*rall_1d
600 
601  !--- make region link table
602  do d = 1, nmax_dmd
603  do i = 1, rall_1d
604  do j = 1, rall_1d
605  rgnid = (d-1)*rall_1dmd + (j-1)*rall_1d + i
606 
607  !--- I_SW
608  if ( j == 1 ) then
609  if ( d <= 5 ) then
610  i_nb = i
611  j_nb = rall_1d
612  d_nb = dmd_data(i_sw,d)
613  direction = i_ne
614  else
615  i_nb = rall_1d
616  j_nb = rall_1d+1-i
617  d_nb = dmd_data(i_sw,d)
618  direction = i_se
619  endif
620  else
621  i_nb = i
622  j_nb = j-1
623  d_nb = d
624  direction = i_ne
625  endif
626  rgnid_nb = (d_nb-1)*rall_1dmd + (j_nb-1)*rall_1d + i_nb
627 
628  edge_tab(i_rgnid,i_sw,rgnid) = rgnid_nb
629  edge_tab(i_dir, i_sw,rgnid) = direction
630 
631  !--- I_NW
632  if ( i == 1 ) then
633  if ( d <= 5 ) then
634  i_nb = rall_1d+1-j
635  j_nb = rall_1d
636  d_nb = dmd_data(i_nw,d)
637  direction = i_ne
638  else
639  i_nb = rall_1d
640  j_nb = j
641  d_nb = dmd_data(i_nw,d)
642  direction = i_se
643  endif
644  else
645  i_nb = i-1
646  j_nb = j
647  d_nb = d
648  direction = i_se
649  endif
650  rgnid_nb = (d_nb-1)*rall_1dmd + (j_nb-1)*rall_1d + i_nb
651 
652  edge_tab(i_rgnid,i_nw,rgnid) = rgnid_nb
653  edge_tab(i_dir, i_nw,rgnid) = direction
654 
655  !--- I_NE
656  if ( j == rall_1d ) then
657  if ( d <= 5 ) then
658  i_nb = 1
659  j_nb = rall_1d+1-i
660  d_nb = dmd_data(i_ne,d)
661  direction = i_nw
662  else
663  i_nb = i
664  j_nb = 1
665  d_nb = dmd_data(i_ne,d)
666  direction = i_sw
667  endif
668  else
669  i_nb = i
670  j_nb = j+1
671  d_nb = d
672  direction = i_sw
673  endif
674  rgnid_nb = (d_nb-1)*rall_1dmd + (j_nb-1)*rall_1d + i_nb
675 
676  edge_tab(i_rgnid,i_ne,rgnid) = rgnid_nb
677  edge_tab(i_dir, i_ne,rgnid) = direction
678 
679  !--- I_SE
680  if ( i == rall_1d ) then
681  if ( d <= 5 ) then
682  i_nb = 1
683  j_nb = j
684  d_nb = dmd_data(i_se,d)
685  direction = i_nw
686  else
687  i_nb = rall_1d+1-j
688  j_nb = 1
689  d_nb = dmd_data(i_se,d)
690  direction = i_sw
691  endif
692  else
693  i_nb = i+1
694  j_nb = j
695  d_nb = d
696  direction = i_nw
697  endif
698  rgnid_nb = (d_nb-1)*rall_1dmd + (j_nb-1)*rall_1d + i_nb
699 
700  edge_tab(i_rgnid,i_se,rgnid) = rgnid_nb
701  edge_tab(i_dir, i_se,rgnid) = direction
702  enddo
703  enddo
704  enddo
705 
706  !--- make region-pe relationship
707  lp2r(:,:) = -1
708 
709  rgnid = 0
710  do p = 0, pall-1
711  do l = 1, lall
712  rgnid = rgnid + 1
713 
714  lp2r(l,p) = rgnid
715  enddo
716  enddo
717 
718  return
719  end Subroutine prc_icoa_rgn_generate
720 
721  !-----------------------------------------------------------------------------
723  subroutine prc_icoa_rgn_vertex_walkaround( &
724  rall, &
725  rall_pl, &
726  vlink, &
727  edge_tab, &
728  vert_num, &
729  vert_tab, &
730  vert_tab_pl )
731  implicit none
732 
733  integer, intent(in) :: rall
734  integer, intent(in) :: rall_pl
735  integer, intent(in) :: vlink
736  integer, intent(in) :: edge_tab(2,4,rall)
737 
738  integer, intent(out) :: vert_num (4,rall)
739  integer, intent(out) :: vert_tab (2,4,rall ,vlink)
740  integer, intent(out) :: vert_tab_pl(2 ,rall_pl,vlink)
741 
742  integer :: rgnid, dir
743  integer :: rgnid_next, dir_next
744 
745  integer :: r, d, v
746  !---------------------------------------------------------------------------
747 
748  vert_num(:,:) = -1
749  vert_tab(:,:,:,:) = -1
750 
751  do r = 1, rall
752  do d = i_w, i_s
753 
754  rgnid = r
755  select case(d)
756  case(i_w)
757  dir = i_sw
758  case(i_n)
759  dir = i_nw
760  case(i_e)
761  dir = i_ne
762  case(i_s)
763  dir = i_se
764  end select
765 
766  v = 0
767  do
768  rgnid_next = edge_tab(i_rgnid,dir,rgnid)
769  dir_next = edge_tab(i_dir, dir,rgnid) - 1
770 
771  if( dir_next == 0 ) dir_next = 4
772  v = v + 1
773  vert_tab(i_rgnid,d,r,v) = rgnid
774  vert_tab(i_dir, d,r,v) = dir
775 
776  rgnid = rgnid_next
777  dir = dir_next
778 
779  if( rgnid == r ) exit
780  enddo
781  vert_num(d,r) = v
782 
783  enddo
784  enddo
785 
786  vert_tab_pl(:,:,:) = -1
787 
788  do r = 1, rall
789  if ( vert_num(i_n,r) == vlink ) then
790  do v = 1, vlink
791  vert_tab_pl(i_rgnid,i_npl,v) = vert_tab(i_rgnid,i_n,r,v)
792  vert_tab_pl(i_dir, i_npl,v) = vert_tab(i_dir, i_n,r,v)
793  enddo
794  exit
795  endif
796  enddo
797 
798  do r = 1, rall
799  if ( vert_num(i_s,r) == vlink ) then
800  do v = 1, vlink
801  vert_tab_pl(i_rgnid,i_spl,v) = vert_tab(i_rgnid,i_s,r,v)
802  vert_tab_pl(i_dir, i_spl,v) = vert_tab(i_dir, i_s,r,v)
803  enddo
804  exit
805  endif
806  enddo
807 
808  return
809  end subroutine prc_icoa_rgn_vertex_walkaround
810 
811  !-----------------------------------------------------------------------------
812  subroutine output_info
813  use scale_prc, only: &
814  prc_myrank
815  implicit none
816 
817  integer :: rgnid, rgnid_next
818  character(len=2) :: dstr, dstr_next
819 
820  integer :: l, d, v
821  !---------------------------------------------------------------------------
822 
823  log_newline
824  log_info("PRC_ICOA_RGN_setup",'(1x,A)') 'Region management information'
825  log_info_cont('(1x,A,A)' ) 'Grid sysytem : Icosahedral'
826  log_info_cont('(1x,A,I7)') 'number of diamond : ', prc_rgn_dmd
827  log_info_cont('(1x,A,I7)') 'maximum number of vertex linkage : ', prc_rgn_vlink
828  log_newline
829  log_info_cont('(1x,A,I7)') 'Region division level (RL) : ', prc_rgn_level
830  log_info_cont('(1x,A,I7,3(A,I4),A)') 'Total number of regular region : ', prc_rgn_total, &
831  ' (', 2**prc_rgn_level, ' x', 2**prc_rgn_level, ' x', prc_rgn_dmd, ' )'
832  log_info_cont('(1x,A,I7)') '# of region per process : ', prc_rgn_local
833  log_info_cont('(1x,A)' ) 'ID of region in my process : '
834  log_info_cont(*) prc_rgn_lp2r(:,prc_myrank)
835  log_info_cont('(1x,A,I7)') 'Region ID, containing north pole : ', prc_rgn_rgn4pl(i_npl)
836  log_info_cont('(1x,A,I7)') 'Region ID, containing south pole : ', prc_rgn_rgn4pl(i_spl)
837  log_info_cont('(1x,A,I7)') 'Process rank, managing north pole : ', prc_rgn_r2p_pl(i_npl)
838  log_info_cont('(1x,A,I7)') 'Process rank, managing south pole : ', prc_rgn_r2p_pl(i_spl)
839 
840  if ( debug ) then
841  log_newline
842  log_info("PRC_ICOA_RGN_setup",'(1x,A)') 'Detailed region management information'
843 
844  log_info_cont(*)'--- (l,myrank) => (rgn)'
845  do l = 1, prc_rgn_local
846  rgnid = prc_rgn_l2r(l)
847  log_info_cont('(1x,A,I4,A,I6,A,I6,A)') '--- (',l,',',prc_myrank,') => (',rgnid,') '
848  enddo
849 
850  log_newline
851  log_info_cont(*)'--- Link information'
852  do l = 1, prc_rgn_local
853  rgnid = prc_rgn_l2r(l)
854 
855  log_newline
856  log_info_cont(*)'--- edge link: (rgn,direction)'
857  do d = i_sw, i_se
858  rgnid_next = prc_rgn_edge_tab(i_rgnid,d,rgnid)
859  dstr = prc_rgn_edgename(d)
860  dstr_next = prc_rgn_edgename(prc_rgn_edge_tab(i_dir,d,rgnid))
861  log_info_cont('(5x,A,I6,A,A,A,I6,A,A,A)') '(',rgnid,',',dstr,') -> (', rgnid_next,',', dstr_next,')'
862  enddo
863 
864  log_info_cont(*)'--- vertex link: (rgn)'
865  do d = i_w, i_s
866  dstr = prc_rgn_vertname(d)
867  log_info_contna('(5x,A,I6,A,A,A)') '(',rgnid,',',dstr,')'
868  do v = 2, prc_rgn_vert_num(d,rgnid)
869  dstr = prc_rgn_vertname(prc_rgn_vert_tab(i_dir,d,rgnid,v))
870  log_info_contna('(A,I6,A,A,A)') ' -> (',prc_rgn_vert_tab(i_rgnid,d,rgnid,v),',',dstr,')'
871  enddo
872  log_newline
873  enddo
874  enddo
875 
876  log_newline
877  log_info_cont(*)'--- Pole information (in the global scope)'
878 
879  log_newline
880  log_info_cont(*)'--- Region ID, containing north pole data : ', prc_rgn_rgn4pl(i_npl)
881  log_info_cont(*)'--- vertex link: (north pole)'
882  do v = 2, prc_rgn_vlink
884  dstr = prc_rgn_vertname(prc_rgn_vert_tab_pl(i_dir,i_npl,v))
885  log_info_contna('(A,I6,A,A,A)') ' -> (',rgnid,',',dstr,')'
886  enddo
888  dstr = prc_rgn_vertname(prc_rgn_vert_tab_pl(i_dir,i_npl,1))
889  log_info_contna('(A,I6,A,A,A)') ' -> (',rgnid,',',dstr,')'
890  log_newline
891  log_info_cont(*)'--- process, managing north pole : ', prc_rgn_r2p_pl(i_npl)
892 
893  log_newline
894  log_info_cont(*)'--- Region ID, containing south pole data : ', prc_rgn_rgn4pl(i_spl)
895  log_info_cont(*)'--- vertex link: (south pole)'
896  do v = 2, prc_rgn_vlink
898  dstr = prc_rgn_vertname(prc_rgn_vert_tab_pl(i_dir,i_spl,v))
899  log_info_contna('(A,I6,A,A,A)') ' -> (',rgnid,',',dstr,')'
900  enddo
902  dstr = prc_rgn_vertname(prc_rgn_vert_tab_pl(i_dir,i_spl,1))
903  log_info_contna('(A,I6,A,A,A)') ' -> (',rgnid,',',dstr,')'
904  log_newline
905  log_info_cont(*)'--- process, managing south pole : ', prc_rgn_r2p_pl(i_spl)
906  endif
907 
908  return
909  end subroutine output_info
910 
911 end module scale_prc_icoa
integer, public prc_universal_nprocs
process num in universal communicator
Definition: scale_prc.F90:72
integer, parameter, public i_prc
process
integer, dimension(:,:,:), allocatable, public prc_rgn_edge_tab
region link information (for 4 edges)
integer, parameter, public i_e
east
integer, parameter, public i_l
local region
logical, public prc_global_ismaster
master process in global communicator?
Definition: scale_prc.F90:81
integer, parameter, public i_se
south east
integer, parameter, public i_spl
south pole
integer, public io_fid_nml
Log file ID (only for output namelist)
Definition: scale_io.F90:57
integer, public prc_global_myrank
myrank in global communicator
Definition: scale_prc.F90:79
integer, dimension(:,:,:,:), allocatable, public prc_rgn_vert_tab
region link information (for 4 vertexes)
integer, dimension(:), allocatable, public prc_rgn_l2r
l,prc_me => rgn
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
integer, parameter, public prc_rgn_local_lim
maximum number of regular region (local)
integer, public prc_universal_comm_world
original communicator
Definition: scale_prc.F90:70
integer, public prc_nprocs
myrank in local communicator
Definition: scale_prc.F90:88
subroutine, public prc_icoa_setup
Setup Processor topology.
integer, public prc_universal_myrank
myrank in universal communicator
Definition: scale_prc.F90:71
integer, parameter, public i_dir
direction
subroutine, public prc_icoa_rgn_generate(rlevel, pall, rall, lall, edge_tab, lp2r)
Generate region management info.
integer, parameter, public i_sw
south west
logical, public io_l
output log or not? (this process)
Definition: scale_io.F90:61
integer, parameter, public i_w
west
integer, parameter, public i_nw
north west
integer, dimension(:,:), allocatable, public prc_rgn_vert_num
number of region around the vertex (4 vertexes)
integer function, public io_get_available_fid()
search & get available file ID
Definition: scale_io.F90:313
module PROCESS
Definition: scale_prc.F90:11
integer, public prc_rgn_level
region division level
integer, dimension(prc_rgn_total_pl), public prc_rgn_r2p_pl
process ID which have the pole regions
integer, dimension(prc_rgn_total_pl), public prc_rgn_rgn4pl
region, having pole data in the halo
integer, parameter, public prc_masterrank
master process in each communicator
Definition: scale_prc.F90:65
integer, public prc_rgn_total
number of regular region (global total)
integer, public prc_global_nprocs
process num in global communicator
Definition: scale_prc.F90:80
integer, parameter, public prc_rank_pl
process ID which manages the pole regions
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:89
integer, parameter, public i_ne
north east
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
integer, public prc_rgn_local
number of regular region (local)
integer, parameter, public i_npl
north pole
integer, parameter, public i_rgnid
region id
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:90
integer, public prc_global_comm_world
global communicator
Definition: scale_prc.F90:78
integer, dimension(:,:,:), allocatable, public prc_rgn_vert_tab_pl
region link information (for 4 vertexes)
integer, parameter, public prc_rgn_total_pl
number of pole region
integer, dimension(:,:), allocatable, public prc_rgn_lp2r
l,prc => rgn
module PRECISION
integer, parameter, public i_s
south
logical, public prc_universal_ismaster
master process in universal communicator?
Definition: scale_prc.F90:73
integer, public prc_local_comm_world
local communicator
Definition: scale_prc.F90:87
integer, public prc_abort_comm_world
communicator for aborting
Definition: scale_prc.F90:95
logical, public io_nml
output log or not? (for namelist, this process)
Definition: scale_io.F90:62
module STDIO
Definition: scale_io.F90:10
integer, parameter, public i_n
north
integer, public prc_rgn_dmd
number of diamonds
logical, public prc_have_pl
this ID manages pole region?
integer, public prc_rgn_vlink
maximum number of vertex linkage, ICO:5
module process / icoA
logical, dimension(:), allocatable, public prc_rgn_have_sgp
region have singlar point?
integer, dimension(:,:), allocatable, public prc_rgn_r2lp
rgn => l,prc