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