SCALE-RM
mod_copytopo.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
20  use scale_tracer
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
28  public :: copytopo
29 
30  !-----------------------------------------------------------------------------
31  !
32  !++ Public parameters & variables
33  !
34  !-----------------------------------------------------------------------------
35  !
36  !++ Private procedure
37  !
38  private :: copytopo_transgrid
39  private :: copytopo_setalpha
40  private :: copytopo_input_data
41  private :: copytopo_mix_data
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private parameters & variables
46  !
47  integer, private, parameter :: handle = 1
48 
49  character(len=H_SHORT), private :: COPYTOPO_IN_FILETYPE = ''
50  character(len=H_LONG), private :: COPYTOPO_IN_BASENAME = ''
51  integer, private :: COPYTOPO_IN_PRC_NUM_X = 0 ! for old files
52  integer, private :: COPYTOPO_IN_PRC_NUM_Y = 0 ! for old files
53  character(len=H_LONG), private :: COPYTOPO_LATLON_CATALOGUE = '' ! for old files
54  integer, private :: COPYTOPO_DOMID
55 
56  character(len=H_LONG), private :: COPYTOPO_IN_POSTFIX = ''
57  character(len=H_SHORT), private :: COPYTOPO_IN_VARNAME = 'topo'
58  character(len=H_SHORT), private :: COPYTOPO_IN_LONNAME = 'lon'
59  character(len=H_SHORT), private :: COPYTOPO_IN_LATNAME = 'lat'
60  real(RP), private :: COPYTOPO_TRANSITION_DX = -1.0_rp
61  real(RP), private :: COPYTOPO_TRANSITION_DY = -1.0_rp
62  real(RP), private :: COPYTOPO_FRACX = 1.0_rp
63  real(RP), private :: COPYTOPO_FRACY = 1.0_rp
64  real(RP), private :: COPYTOPO_taux = 1.0_rp
65  real(RP), private :: COPYTOPO_tauy = 1.0_rp
66  logical, private :: COPYTOPO_ENTIRE_REGION = .false.
67  logical, private :: COPYTOPO_LINEAR_H = .true.
68  real(RP), private :: COPYTOPO_EXP_H = 2.0_rp
69  character(len=H_SHORT), private :: COPYTOPO_INTRP_TYPE = "LINEAR"
70  integer, private :: COPYTOPO_hypdiff_ORDER = 4
71  integer, private :: COPYTOPO_hypdiff_NITER = 20
72 
73  !-----------------------------------------------------------------------------
74 contains
75  !-----------------------------------------------------------------------------
77  subroutine copytopo( &
78  TOPO_child )
79  use scale_prc, only: &
80  prc_abort
81  implicit none
82 
83  real(rp), intent(inout) :: topo_child(:,:)
84 
85  real(rp) :: ctrx (ia)
86  real(rp) :: ctry (ja)
87  real(rp) :: alpha (ia,ja)
88  real(rp) :: topo_parent(ia,ja)
89 
90  namelist / param_copytopo / &
91  copytopo_in_filetype, &
92  copytopo_in_basename, &
93  copytopo_in_prc_num_x, & ! for old files
94  copytopo_in_prc_num_y, & ! for old files
95  copytopo_latlon_catalogue, & ! for old files
96  copytopo_in_postfix, &
97  copytopo_in_varname, &
98  copytopo_in_lonname, &
99  copytopo_in_latname, &
100  copytopo_transition_dx, &
101  copytopo_transition_dy, &
102  copytopo_fracx, &
103  copytopo_fracy, &
104  copytopo_taux, &
105  copytopo_tauy, &
106  copytopo_entire_region, &
107  copytopo_linear_h, &
108  copytopo_exp_h, &
109  copytopo_hypdiff_order, &
110  copytopo_hypdiff_niter
111 
112  integer :: ierr
113  !---------------------------------------------------------------------------
114 
115  log_newline
116  log_info("COPYTOPO",*) 'Setup'
117 
118  !--- read namelist
119  rewind(io_fid_conf)
120  read(io_fid_conf,nml=param_copytopo,iostat=ierr)
121  if( ierr < 0 ) then !--- missing
122  log_info("COPYTOPO",*) 'Not found namelist. Default used.'
123  elseif( ierr > 0 ) then !--- fatal error
124  log_error("COPYTOPO",*) 'Not appropriate names in namelist PARAM_COPYTOPO. Check!'
125  call prc_abort
126  endif
127  log_nml(param_copytopo)
128 
129  ! copy topography from parent domain to transition region
130 
131  call copytopo_transgrid( ctrx(:), ctry(:) ) ! [OUT]
132 
133  call copytopo_setalpha ( ctrx(:), ctry(:), & ! [IN]
134  alpha(:,:) ) ! [OUT]
135 
136  call copytopo_input_data( topo_parent(:,:) ) ! [OUT]
137 
138  call copytopo_mix_data( topo_parent(:,:), & ! [IN]
139  alpha(:,:), & ! [IN]
140  topo_child(:,:) ) ! [INOUT]
141 
142  return
143  end subroutine copytopo
144 
145  !-----------------------------------------------------------------------------
147  subroutine copytopo_transgrid( &
148  CTRX, &
149  CTRY )
150  use scale_prc, only: &
151  prc_abort, &
152  prc_myrank
153  use scale_prc_cartesc, only: &
154  prc_2drank
155  use scale_const, only: &
156  eps => const_eps
157  use scale_atmos_grid_cartesc, only: &
158  cxg => atmos_grid_cartesc_cxg, &
159  fxg => atmos_grid_cartesc_fxg, &
160  cyg => atmos_grid_cartesc_cyg, &
161  fyg => atmos_grid_cartesc_fyg, &
162  cbfxg => atmos_grid_cartesc_cbfxg, &
163  cbfyg => atmos_grid_cartesc_cbfyg
164  implicit none
165 
166  real(rp), intent(out) :: ctrx(ia)
167  real(rp), intent(out) :: ctry(ja)
168 
169  real(rp) :: ctrxg ( iag)
170  real(rp) :: ctryg ( jag)
171 
172  real(rp) :: bufftotx, transtotx
173  real(rp) :: bufftoty, transtoty
174  integer :: imain, ibuff, itrans
175  integer :: jmain, jbuff, jtrans
176  integer :: copy_is, copy_ie
177  integer :: copy_js, copy_je
178 
179  integer :: i, j, ii, jj
180  !---------------------------------------------------------------------------
181 
182  ! X-direction
183  ! calculate buffer grid size
184 
185  do i = ihalo+1, iag
186  if( abs(cbfxg(i)) < eps ) exit
187  enddo
188  ibuff = i - 1 - ihalo
189  if ( ibuff == 0 ) then
190  bufftotx = 0.0_rp
191  else
192  bufftotx = fxg(ibuff+ihalo) - fxg(ihalo)
193  end if
194 
195  if ( copytopo_transition_dx < 0.0_rp ) then
196  copytopo_transition_dx = bufftotx
197  end if
198  do i = ibuff+ihalo+1, iag
199  if( cxg(i) - bufftotx - fxg(ihalo) >= copytopo_transition_dx ) exit
200  enddo
201  itrans = i - 1 - ihalo - ibuff
202  if ( itrans == 0 ) then
203  transtotx = 0.0_rp
204  else
205  transtotx = fxg(itrans+ibuff+ihalo) - fxg(ihalo) - bufftotx
206  end if
207 
208  imain = iag - 2*ibuff - 2*itrans - 2*ihalo
209 
210  if ( imain < 1 ) then
211  log_error("COPYTOPO_transgrid",*) 'Not appropriate transition width for global domain(X).', copytopo_transition_dx
212  log_error_cont(*) '# of buffer region (one side) = ', ibuff
213  log_error_cont(*) '# of transion region (one side) = ', itrans
214  call prc_abort
215  endif
216 
217  ! calc transition factor (global domaim)
218  ctrxg(:) = 0.0_rp
219 
220  copy_is = 1
221  copy_ie = ihalo+ibuff
222  !$omp parallel do
223  do i = copy_is, copy_ie
224  ctrxg(i) = 1.0_rp
225  enddo
226 
227  if ( itrans > 0 ) then
228  copy_is = ihalo+ibuff+1
229  copy_ie = ihalo+ibuff+itrans
230  !$omp parallel do
231  do i = copy_is, copy_ie
232  ctrxg(i) = (transtotx+bufftotx+fxg(ihalo )-cxg(i)) / transtotx
233  enddo
234  copy_is = iag - ihalo - ibuff - itrans + 1
235  copy_ie = iag - ihalo - ibuff
236  !$omp parallel do
237  do i = copy_is, copy_ie
238  ctrxg(i) = (transtotx+bufftotx-fxg(iag-ihalo)+cxg(i)) / transtotx
239  enddo
240  endif
241 
242  copy_is = iag - ihalo - ibuff + 1
243  copy_ie = iag
244  !$omp parallel do
245  do i = copy_is, copy_ie
246  ctrxg(i) = 1.0_rp
247  enddo
248 
249  ctrxg(:) = max( min( ctrxg(:), 1.0_rp ), 0.0_rp )
250 
251  ! Y-direction
252  ! calculate buffer grid size
253  do j = jhalo+1, jag
254  if( abs(cbfyg(j)) < eps ) exit
255  enddo
256  jbuff = j - 1 - jhalo
257  if ( jbuff == 0 ) then
258  bufftoty = 0.0_rp
259  else
260  bufftoty = fyg(jbuff+jhalo) - fyg(jhalo)
261  end if
262 
263  if ( copytopo_transition_dy < 0.0_rp ) then
264  copytopo_transition_dy = bufftoty
265  end if
266  do j = jbuff+jhalo+1, jag
267  if( cyg(j) - bufftoty - fyg(jhalo) >= copytopo_transition_dy ) exit
268  enddo
269  jtrans = j - 1 - jhalo - jbuff
270  if ( jtrans == 0 ) then
271  transtoty = 0.0_rp
272  else
273  transtoty = fyg(jtrans+jbuff+jhalo) - fyg(jhalo) - bufftoty
274  end if
275 
276  jmain = jag - 2*jbuff - 2*jtrans - 2*jhalo
277 
278  if ( jmain < 1 ) then
279  log_error("COPYTOPO_transgrid",*) 'Not appropriate transition width for global domain(Y).', copytopo_transition_dy
280  log_error_cont(*) '# of buffer region (one side)', jbuff
281  log_error_cont(*) '# of transion region (one side)', jtrans
282  call prc_abort
283  endif
284 
285  ! calc transition factor (global domaim)
286  ctryg(:) = 0.0_rp
287 
288  copy_js = 1
289  copy_je = jhalo+jbuff
290  !$omp parallel do
291  do j = copy_js, copy_je
292  ctryg(j) = 1.0_rp
293  enddo
294 
295  if ( jtrans > 0 ) then
296  copy_js = jhalo+jbuff+1
297  copy_je = jhalo+jbuff+jtrans
298  !$omp parallel do
299  do j = copy_js, copy_je
300  ctryg(j) = (transtoty+bufftoty+fyg(jhalo )-cyg(j)) / transtoty
301  enddo
302  copy_js = jag - jhalo - jbuff - jtrans + 1
303  copy_je = jag - jhalo - jbuff
304  !$omp parallel do
305  do j = copy_js, copy_je
306  ctryg(j) = (transtoty+bufftoty-fyg(jag-jhalo)+cyg(j)) / transtoty
307  enddo
308  endif
309 
310  copy_js = jag - jhalo - jbuff + 1
311  copy_je = jag
312  !$omp parallel do
313  do j = copy_js, copy_je
314  ctryg(j) = 1.0_rp
315  enddo
316 
317  ctryg(:) = max( min( ctryg(:), 1.0_rp ), 0.0_rp )
318 
319 
320 
321  ! horizontal coordinate (local domaim)
322  !$omp parallel do &
323  !$omp private(ii)
324  do i = 1, ia
325  ii = i + prc_2drank(prc_myrank,1) * imax
326  ctrx(i) = ctrxg(ii)
327  enddo
328 
329  !$omp parallel do &
330  !$omp private(jj)
331  do j = 1, ja
332  jj = j + prc_2drank(prc_myrank,2) * jmax
333  ctry(j) = ctryg(jj)
334  enddo
335 
336  log_info("COPYTOPO_transgrid",*) 'COPYTOPO grid Information'
337  log_info_cont(*) 'size of buffer region (x and y; one side) = ', bufftotx, bufftoty
338  log_info_cont(*) '# of buffer region (x and y; one side) = ', ibuff, jbuff
339  log_info_cont(*) 'size of transion region (x and y; one side) = ', transtotx, transtoty
340  log_info_cont(*) '# of transion region (x and y; one side) = ', itrans, jtrans
341 
342 
343  return
344  end subroutine copytopo_transgrid
345 
346 
347  !-----------------------------------------------------------------------------
349  subroutine copytopo_setalpha( &
350  CTRX, &
351  CTRY, &
352  alpha )
353  use scale_const, only: &
354  eps => const_eps
355  use scale_comm_cartesc, only: &
356  comm_vars8, &
357  comm_wait
358  implicit none
359 
360  real(rp), intent(in) :: ctrx (ia)
361  real(rp), intent(in) :: ctry (ja)
362  real(rp), intent(out) :: alpha(ia,ja)
363 
364  real(rp) :: coef_x, alpha_x1
365  real(rp) :: coef_y, alpha_y1
366  real(rp) :: ee1
367 
368  integer :: i, j
369  !---------------------------------------------------------------------------
370 
371  ! check invalid fraction
372  copytopo_fracx = max( min( copytopo_fracx, 1.0_rp ), eps )
373  copytopo_fracy = max( min( copytopo_fracy, 1.0_rp ), eps )
374 
375  if ( copytopo_taux <= 0.0_rp ) then ! invalid tau
376  coef_x = 0.0_rp
377  else
378  coef_x = 1.0_rp / copytopo_taux
379  endif
380 
381  if ( copytopo_tauy <= 0.0_rp ) then ! invalid tau
382  coef_y = 0.0_rp
383  else
384  coef_y = 1.0_rp / copytopo_tauy
385  endif
386 
387  !$omp parallel do &
388  !$omp private(ee1,alpha_x1,alpha_y1)
389  do j = 1, ja
390  do i = 1, ia
391  ee1 = ctrx(i)
392 
393  if ( ee1 <= 1.0_rp - copytopo_fracx ) then
394  ee1 = 0.0_rp
395  else
396  ee1 = ( ee1 - 1.0_rp + copytopo_fracx ) / copytopo_fracx
397  endif
398 
399  if ( copytopo_linear_h ) then
400  alpha_x1 = coef_x * ee1
401  else
402  alpha_x1 = coef_x * ee1 * exp( -(1.0_rp-ee1) * copytopo_exp_h )
403  endif
404 
405  ee1 = ctry(j)
406 
407  if ( ee1 <= 1.0_rp - copytopo_fracy ) then
408  ee1 = 0.0_rp
409  else
410  ee1 = ( ee1 - 1.0_rp + copytopo_fracy ) / copytopo_fracy
411  endif
412 
413  if ( copytopo_linear_h ) then
414  alpha_y1 = coef_y * ee1
415  else
416  alpha_y1 = coef_y * ee1 * exp( -(1.0_rp-ee1) * copytopo_exp_h )
417  endif
418 
419  alpha(i,j) = max( alpha_x1, alpha_y1 )
420  enddo
421  enddo
422 
423  call comm_vars8( alpha(:,:), 1 )
424  call comm_wait ( alpha(:,:), 1 )
425 
426  return
427  end subroutine copytopo_setalpha
428 
429 
430  !-----------------------------------------------------------------------------
432  subroutine copytopo_input_data( &
433  TOPO_parent )
434  use scale_prc, only: &
435  prc_abort, &
437  use scale_const, only: &
438  eps => const_eps, &
439  pi => const_pi
440  use scale_interp, only: &
442  interp_factor2d, &
444  use scale_comm_cartesc, only: &
445  comm_bcast, &
446  comm_vars8, &
447  comm_wait
448  use scale_atmos_grid_cartesc_real, only: &
451  use scale_atmos_grid_cartesc, only: &
452  cx => atmos_grid_cartesc_cx, &
454  use scale_mapprojection, only: &
455  mapprojection_lonlat2xy
456  use scale_filter, only: &
457  filter_hyperdiff
458  use scale_landuse, only: &
460  use scale_comm_cartesc_nest, only: &
461  nest_interp_level => comm_cartesc_nest_interp_level
462  implicit none
463 
464  real(rp), intent(out) :: topo_parent(:,:)
465 
466  integer :: ia_org, ja_org
467 
468  real(rp), allocatable :: lon_org (:,:), lat_org(:,:)
469  real(rp), allocatable :: x_org (:,:), y_org (:,:)
470  real(rp), allocatable :: topo_org(:,:)
471 
472  logical :: single
473 
474  ! for interpolation
475  real(rp) :: dummy(1,1)
476  integer :: idx_i(ia,ja,nest_interp_level)
477  integer :: idx_j(ia,ja,nest_interp_level)
478  real(rp) :: hfact(ia,ja,nest_interp_level)
479  integer :: interp_level
480  logical :: zonal, pole
481 
482  real(rp) :: topo_sign(ia,ja)
483  real(rp) :: ocean_flag
484 
485  integer :: i, j
486  !---------------------------------------------------------------------------
487 
488  select case ( copytopo_in_filetype )
489  case ( 'SCALE' )
490  call copytopo_get_size_scale( ia_org, ja_org ) ! (out)
491  single = .false.
492  case ( 'GrADS' )
493  if ( prc_ismaster ) call copytopo_get_size_grads( ia_org, ja_org ) ! (out)
494  single = .true.
495  case ( 'WRF-ARW' )
496  if ( prc_ismaster ) call copytopo_get_size_wrfarw( ia_org, ja_org ) ! (out)
497  single = .true.
498  case default
499  log_error("COPYTOPO_input_data",*) 'COPYTOPO_IN_FILETYPE must be "SCALE" or "GrADS"'
500  call prc_abort
501  end select
502 
503  if ( single ) then
504  call comm_bcast( ia_org )
505  call comm_bcast( ja_org )
506  end if
507 
508  allocate( lon_org(ia_org,ja_org) )
509  allocate( lat_org(ia_org,ja_org) )
510  allocate( topo_org(ia_org,ja_org) )
511 
512  select case ( copytopo_in_filetype )
513  case ( 'SCALE' )
514  call copytopo_get_data_scale( ia_org, ja_org, & ! (in)
515  lon_org(:,:), lat_org(:,:), topo_org(:,:) ) ! (out)
516  case ( 'GrADS' )
517  if ( prc_ismaster ) &
518  call copytopo_get_data_grads( ia_org, ja_org, & ! (in)
519  lon_org(:,:), lat_org(:,:), topo_org(:,:) ) ! (out)
520  case ( 'WRF-ARW' )
521  if ( prc_ismaster ) &
522  call copytopo_get_data_wrfarw( ia_org, ja_org, & ! (in)
523  lon_org(:,:), lat_org(:,:), topo_org(:,:) ) ! (out)
524  end select
525 
526  if ( single ) then
527  call comm_bcast( ia_org, ja_org, lon_org(:,:) )
528  call comm_bcast( ia_org, ja_org, lat_org(:,:) )
529  call comm_bcast( ia_org, ja_org, topo_org(:,:) )
530  end if
531 
532  call interp_domain_compatibility( lon_org(:,:), lat_org(:,:), dummy(:,:), &
533  lon(:,:), lat(:,:), dummy(:,:), dummy(:,:), &
534  skip_z=.true. )
535 
536  select case ( copytopo_intrp_type )
537  case ( "LINEAR" )
538 
539  allocate( x_org(ia_org,ja_org) )
540  allocate( y_org(ia_org,ja_org) )
541  call mapprojection_lonlat2xy( ia_org, 1, ia_org, ja_org, 1, ja_org, &
542  lon_org(:,:), lat_org(:,:), & ! [IN]
543  x_org(:,:), y_org(:,:) ) ! [OUT]
544 
545  zonal = ( maxval(lon_org) - minval(lon_org) ) > 2.0_rp * pi * 0.9_rp
546  pole = ( maxval(lat_org) > pi * 0.5_rp * 0.9_rp ) .or. ( minval(lat_org) < - pi * 0.5_rp * 0.9_rp )
547  call interp_factor2d( ia_org, ja_org, ia, ja, & ! [IN]
548  x_org(:,:), y_org(:,:), & ! [IN]
549  cx(:), cy(:), & ! [IN]
550  idx_i(:,:,:), idx_j(:,:,:), & ! [OUT]
551  hfact(:,:,:), & ! [OUT]
552  zonal = zonal, pole = pole ) ! [IN]
553  deallocate( x_org, y_org )
554  interp_level = 4
555 
556  case ( "DIST-WEIGHT" )
557 
558  call interp_factor2d( nest_interp_level, & ! [IN]
559  ia_org, ja_org, ia, ja, & ! [IN]
560  lon_org(:,:), lat_org(:,:), & ! [IN]
561  lon(:,:), lat(:,:), & ! [IN]
562  idx_i(:,:,:), idx_j(:,:,:), & ! [OUT]
563  hfact(:,:,:) ) ! [OUT]
564  interp_level = nest_interp_level
565 
566  case default
567  log_error("COPYTOPO_input_data",*) "Unsupported type of COPYTOPO_INTRP_TYPE : ", trim(copytopo_intrp_type)
568  log_error_cont(*) 'It must be "LINEAR" or "DIST-WEIGHT"'
569  call prc_abort
570  end select
571 
572  call interp_interp2d( interp_level, & ! [IN]
573  ia_org, ja_org, & ! [IN]
574  ia, ja, & ! [IN]
575  idx_i(:,:,:), & ! [IN]
576  idx_j(:,:,:), & ! [IN]
577  hfact(:,:,:), & ! [IN]
578  topo_org(:,:), & ! [IN]
579  topo_parent(:,:) ) ! [OUT]
580 
581  deallocate( topo_org, lon_org, lat_org )
582 
583  if ( copytopo_hypdiff_niter > 1 ) then
584  !$omp parallel do &
585  !$omp private(ocean_flag)
586  do j = 1, ja
587  do i = 1, ia
588  ocean_flag = ( 0.5_rp + sign( 0.5_rp, landuse_fact_ocean(i,j) - 1.0_rp + eps ) ) & ! fact_ocean==1
589  * ( 0.5_rp + sign( 0.5_rp, eps - abs(topo_parent(i,j)) ) ) ! |TOPO| > EPS
590  topo_sign(i,j) = sign( 1.0_rp, topo_parent(i,j) ) * ( 1.0_rp - ocean_flag )
591  end do
592  end do
593 
594  call filter_hyperdiff( ia, isb, ieb, ja, jsb, jeb, &
595  topo_parent(:,:), copytopo_hypdiff_order, copytopo_hypdiff_niter, &
596  limiter_sign = topo_sign(:,:) )
597  end if
598 
599 
600  call comm_vars8( topo_parent(:,:), 1 )
601  call comm_wait ( topo_parent(:,:), 1 )
602 
603 
604  return
605  end subroutine copytopo_input_data
606 
607  !-----------------------------------------------------------------------------
609  subroutine copytopo_mix_data( &
610  TOPO_parent, &
611  alpha, &
612  TOPO_child )
613  implicit none
614 
615  real(rp), intent(in) :: topo_parent(:,:) ! topography of parent domain
616  real(rp), intent(in) :: alpha (:,:) ! dumping coefficient (0-1)
617  real(rp), intent(inout) :: topo_child (:,:) ! topography of child domain
618 
619  integer :: i, j
620  !---------------------------------------------------------------------------
621 
622  if ( copytopo_entire_region ) then
623  !$omp parallel do
624  do j = 1, ja
625  do i = 1, ia
626  topo_child(i,j) = topo_parent(i,j)
627  enddo
628  enddo
629  else
630  !$omp parallel do
631  do j = 1, ja
632  do i = 1, ia
633  topo_child(i,j) = ( 1.0_rp-alpha(i,j) ) * topo_child(i,j) &
634  + ( alpha(i,j) ) * topo_parent(i,j)
635  enddo
636  enddo
637  endif
638 
639  return
640  end subroutine copytopo_mix_data
641 
642  !-----------------------------------------------------------------------------
643  subroutine copytopo_get_size_scale( &
644  IA_org, JA_org )
645  use scale_comm_cartesc_nest, only: &
648  implicit none
649  integer, intent(out) :: ia_org
650  integer, intent(out) :: ja_org
651 
653  copytopo_domid, & ! [OUT]
654  copytopo_in_basename, & ! [IN]
655  copytopo_in_prc_num_x, & ! [IN]
656  copytopo_in_prc_num_y, & ! [IN]
657  copytopo_latlon_catalogue ) ! [IN]
658 
660  copytopo_domid, & ! [IN]
661  imaxg=ia_org, jmaxg=ja_org ) ! [OUT]
662 
663  return
664  end subroutine copytopo_get_size_scale
665 
666  !-----------------------------------------------------------------------------
667  subroutine copytopo_get_data_scale( &
668  IA_org, JA_org, &
669  LON_org, LAT_org, &
670  TOPO_org )
671  use scale_const, only: &
672  d2r => const_d2r
673  use scale_comm_cartesc_nest, only: &
676  use scale_file, only: &
677  file_open, &
678  file_read, &
679  file_close
680  implicit none
681  integer, intent(in) :: IA_org, JA_org
682  real(RP), intent(out) :: LON_org (IA_org,JA_org)
683  real(RP), intent(out) :: LAT_org (IA_org,JA_org)
684  real(RP), intent(out) :: TOPO_org(IA_org,JA_org)
685 
686  integer :: num_tile
687  integer, allocatable :: tile_id(:)
688 
689  real(RP), allocatable :: read2D(:,:)
690 
691  integer :: tilei, tilej
692  integer :: cxs, cxe, cys, cye ! for child domain
693  integer :: pxs, pxe, pys, pye ! for parent domain
694  integer :: fid
695 
696  integer :: rank
697  integer :: n
698 
699  call comm_cartesc_nest_parent_info( copytopo_domid, & ! [IN]
700  num_tile = num_tile )
701  allocate( tile_id(num_tile) )
702  call comm_cartesc_nest_parent_info( copytopo_domid, & ! [IN]
703  tile_id = tile_id )
704 
705  do n = 1, num_tile
706  ! read data from split files
707  rank = tile_id(n)
708 
709  call comm_cartesc_nest_domain_shape( tilei, tilej, & ! [OUT]
710  cxs, cxe, & ! [OUT]
711  cys, cye, & ! [OUT]
712  pxs, pxe, & ! [OUT]
713  pys, pye, & ! [OUT]
714  copytopo_domid, & ! [IN]
715  n ) ! [IN]
716 
717  allocate( read2d(tilei,tilej) )
718 
719  call file_open( copytopo_in_basename, & ! [IN]
720  fid, & ! [OUT]
721  aggregate=.false., rankid=rank ) ! [IN]
722 
723  call file_read( fid, "lon", read2d(:,:) )
724  lon_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
725  call file_read( fid, "lat", read2d(:,:) )
726  lat_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
727  call file_read( fid, copytopo_in_varname, read2d(:,:) )
728  topo_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye)
729  deallocate( read2d )
730 
731  call file_close( fid )
732 
733  enddo
734 
735  return
736  end subroutine copytopo_get_data_scale
737 
738  !-----------------------------------------------------------------------------
739  subroutine copytopo_get_size_grads( &
740  IA_org, JA_org )
741  use scale_file_grads, only: &
742  file_grads_open, &
743  file_grads_get_shape
744  implicit none
745  integer, intent(out) :: IA_org, JA_org
746 
747  integer :: file_id
748  integer :: shape(2)
749 
750  call file_grads_open( copytopo_in_basename, & ! (in)
751  file_id ) ! (out)
752 
753  call file_grads_get_shape( file_id, copytopo_in_varname, & ! (in)
754  shape(:) ) ! (out)
755 
756  ia_org = shape(1)
757  ja_org = shape(2)
758 
759  return
760  end subroutine copytopo_get_size_grads
761 
762  !-----------------------------------------------------------------------------
763  subroutine copytopo_get_data_grads( &
764  IA_org, JA_org, &
765  LON_org, LAT_org, &
766  TOPO_org )
767  use scale_const, only: &
768  d2r => const_d2r
769  use scale_file_grads, only: &
770  file_grads_open, &
773  file_grads_read, &
775  implicit none
776  integer, intent(in) :: IA_org, JA_org
777  real(RP), intent(out) :: LON_org (IA_org,JA_org)
778  real(RP), intent(out) :: LAT_org (IA_org,JA_org)
779  real(RP), intent(out) :: TOPO_org(IA_org,JA_org)
780 
781  integer :: file_id, var_id
782  real(RP) :: lon1D(IA_org), lat1D(JA_org)
783 
784  integer :: i, j
785 
786  call file_grads_open( copytopo_in_basename, & ! (in)
787  file_id ) ! (out)
788 
789  ! lon
790  call file_grads_varid( file_id, copytopo_in_lonname, & ! (in)
791  var_id ) ! (out)
792  if ( file_grads_isoned( file_id, var_id ) ) then
793  call file_grads_read( file_id, var_id, & ! (in)
794  lon1d(:) ) ! (out)
795  !$omp parallel do collapse(2)
796  do j = 1, ja_org
797  do i = 1, ia_org
798  lon_org(i,j) = lon1d(i) * d2r
799  end do
800  end do
801  else
802  call file_grads_read( file_id, var_id, & ! (in)
803  lon_org(:,:), & ! (out)
804  postfix = copytopo_in_postfix ) ! (in)
805  !$omp parallel do collapse(2)
806  do j = 1, ja_org
807  do i = 1, ia_org
808  lon_org(i,j) = lon_org(i,j) * d2r
809  end do
810  end do
811  end if
812 
813  ! lat
814  call file_grads_varid( file_id, copytopo_in_latname, & ! (in)
815  var_id ) ! (out)
816  if ( file_grads_isoned( file_id, var_id ) ) then
817  call file_grads_read( file_id, var_id, & ! (in)
818  lat1d(:) ) ! (out)
819  !$omp parallel do collapse(2)
820  do j = 1, ja_org
821  do i = 1, ia_org
822  lat_org(i,j) = lat1d(j) * d2r
823  end do
824  end do
825  else
826  call file_grads_read( file_id, var_id, & ! (in)
827  lat_org(:,:), & ! (out)
828  postfix = copytopo_in_postfix ) ! (in)
829  !$omp parallel do collapse(2)
830  do j = 1, ja_org
831  do i = 1, ia_org
832  lat_org(i,j) = lat_org(i,j) * d2r
833  end do
834  end do
835  end if
836 
837  ! topo
838  call file_grads_varid( file_id, copytopo_in_varname, & ! (in)
839  var_id ) ! (out)
840  call file_grads_read( file_id, var_id, & ! (in)
841  topo_org(:,:), & ! (out)
842  postfix = copytopo_in_postfix ) ! (in)
843 
844  ! close
845  call file_grads_close( file_id )
846 
847  return
848  end subroutine copytopo_get_data_grads
849 
850  !-----------------------------------------------------------------------------
851  subroutine copytopo_get_size_wrfarw( &
852  IA_org, JA_org )
853  use scale_file_grads, only: &
854  file_grads_open, &
855  file_grads_get_shape
856  use scale_file, only: &
857  file_open, &
859  implicit none
860  integer, intent(out) :: IA_org, JA_org
861 
862  integer :: fid
863 
864  call file_open( copytopo_in_basename, & ! (in)
865  fid, & ! (out)
866  postfix = copytopo_in_postfix, & ! (in)
867  single = .true. ) ! (in)
868 
869  call file_get_dimlength( fid, "west_east", ia_org )
870  call file_get_dimlength( fid, "south_north", ja_org )
871 
872  return
873  end subroutine copytopo_get_size_wrfarw
874 
875  subroutine copytopo_get_data_wrfarw( &
876  IA_org, JA_org, &
877  LON_org, LAT_org, &
878  TOPO_org )
879  use scale_const, only: &
880  d2r => const_d2r
881  use scale_file, only: &
882  file_open, &
883  file_read, &
884  file_close
885  implicit none
886  integer, intent(in) :: IA_org, JA_org
887  real(RP), intent(out) :: LON_org (IA_org,JA_org)
888  real(RP), intent(out) :: LAT_org (IA_org,JA_org)
889  real(RP), intent(out) :: TOPO_org(IA_org,JA_org)
890 
891  integer :: fid
892  integer :: i, j
893 
894  call file_open( copytopo_in_basename, & ! (in)
895  fid, & ! (out)
896  postfix = copytopo_in_postfix, & ! (in)
897  single = .true. ) ! (in)
898 
899  call file_read( fid, "XLAT", lat_org(:,:) )
900 
901  call file_read( fid, "XLONG", lon_org(:,:) )
902 
903  !$omp parallel do collapse(2)
904  do j = 1, ja_org
905  do i = 1, ia_org
906  lat_org(i,j) = lat_org(i,j) * d2r
907  lon_org(i,j) = lon_org(i,j) * d2r
908  end do
909  end do
910 
911  call file_read( fid, "HGT", topo_org(:,:) )
912 
913  call file_close( fid )
914 
915  end subroutine copytopo_get_data_wrfarw
916 
917 end module mod_copytopo
scale_interp::interp_interp2d
subroutine, public interp_interp2d(npoints, IA_ref, JA_ref, IA, JA, idx_i, idx_j, hfact, val_ref, val, threshold_undef, wsum, val2)
Definition: scale_interp.F90:1376
scale_comm_cartesc_nest::comm_cartesc_nest_interp_level
integer, public comm_cartesc_nest_interp_level
horizontal interpolation level
Definition: scale_comm_cartesC_nest.F90:76
scale_atmos_grid_cartesc_index::isb
integer, public isb
Definition: scale_atmos_grid_cartesC_index.F90:64
scale_landuse::landuse_fact_ocean
real(rp), dimension(:,:), allocatable, public landuse_fact_ocean
ocean factor
Definition: scale_landuse.F90:45
scale_comm_cartesc_nest::comm_cartesc_nest_domain_shape
subroutine, public comm_cartesc_nest_domain_shape(tilei, tilej, cxs, cxe, cys, cye, pxs, pxe, pys, pye, dom_id, iloc, xstg, ystg)
Return shape of ParentDomain at the specified rank (for offline)
Definition: scale_comm_cartesC_nest.F90:1187
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
scale_file_grads::file_grads_isoned
logical function, public file_grads_isoned(file_id, var_id)
Definition: scale_file_grads.F90:373
scale_file_grads::file_grads_open
subroutine, public file_grads_open(file_name, file_id)
Open.
Definition: scale_file_grads.F90:104
scale_atmos_grid_cartesc_index::ihalo
integer, public ihalo
Definition: scale_atmos_grid_cartesC_index.F90:44
scale_file::file_open
subroutine, public file_open(basename, fid, mode, single, allnodes, aggregate, rankid, postfix)
Definition: scale_file.F90:536
scale_precision
module PRECISION
Definition: scale_precision.F90:14
scale_comm_cartesc_nest
module Communication CartesianC nesting
Definition: scale_comm_cartesC_nest.F90:12
scale_interp
module INTERPOLATION
Definition: scale_interp.F90:12
scale_const::const_eps
real(rp), public const_eps
small number
Definition: scale_const.F90:35
scale_atmos_grid_cartesc::atmos_grid_cartesc_fyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fyg
face coordinate [m]: y, global
Definition: scale_atmos_grid_cartesC.F90:79
scale_prc::prc_myrank
integer, public prc_myrank
process num in local communicator
Definition: scale_prc.F90:91
mod_copytopo::copytopo
subroutine, public copytopo(TOPO_child)
Setup and Main.
Definition: mod_copytopo.F90:79
scale_file_grads::file_grads_varid
subroutine, public file_grads_varid(file_id, var_name, var_id)
Definition: scale_file_grads.F90:312
scale_file::file_get_dimlength
subroutine, public file_get_dimlength(fid, dimname, len, error)
get length of dimension
Definition: scale_file.F90:633
scale_atmos_grid_cartesc_index::imaxg
integer, public imaxg
Definition: scale_atmos_grid_cartesC_index.F90:72
scale_comm_cartesc_nest::comm_cartesc_nest_parent_info
subroutine, public comm_cartesc_nest_parent_info(dom_id, KMAX, LKMAX, IMAXG, JMAXG, num_tile, tile_id)
Return infomation of parent domain (for offline)
Definition: scale_comm_cartesC_nest.F90:1130
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lon
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
Definition: scale_atmos_grid_cartesC_real.F90:49
mod_copytopo::copytopo_get_data_wrfarw
subroutine copytopo_get_data_wrfarw(IA_org, JA_org, LON_org, LAT_org, TOPO_org)
Definition: mod_copytopo.F90:879
mod_copytopo::copytopo_get_data_grads
subroutine copytopo_get_data_grads(IA_org, JA_org, LON_org, LAT_org, TOPO_org)
Definition: mod_copytopo.F90:767
scale_file
module file
Definition: scale_file.F90:15
scale_atmos_grid_cartesc_index::imax
integer, public imax
Definition: scale_atmos_grid_cartesC_index.F90:37
scale_atmos_grid_cartesc_index::jmaxg
integer, public jmaxg
Definition: scale_atmos_grid_cartesC_index.F90:73
scale_atmos_grid_cartesc_index::jeb
integer, public jeb
Definition: scale_atmos_grid_cartesC_index.F90:67
scale_prc
module PROCESS
Definition: scale_prc.F90:11
scale_mapprojection
module Map projection
Definition: scale_mapprojection.F90:12
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc_index::jag
integer, public jag
Definition: scale_atmos_grid_cartesC_index.F90:75
scale_atmos_grid_cartesc_index::iag
integer, public iag
Definition: scale_atmos_grid_cartesC_index.F90:74
scale_interp::interp_domain_compatibility
subroutine, public interp_domain_compatibility(lon_org, lat_org, topc_org, lon_loc, lat_loc, topc_loc, topf_loc, skip_x, skip_y, skip_z)
Definition: scale_interp.F90:150
scale_file::file_close
subroutine, public file_close(fid, abort)
Definition: scale_file.F90:6234
scale_atmos_grid_cartesc::atmos_grid_cartesc_cxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cxg
center coordinate [m]: x, global
Definition: scale_atmos_grid_cartesC.F90:76
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
scale_const
module CONSTANT
Definition: scale_const.F90:11
scale_atmos_grid_cartesc_index::ia
integer, public ia
Definition: scale_atmos_grid_cartesC_index.F90:48
scale_comm_cartesc_nest::comm_cartesc_nest_domain_regist_file
subroutine, public comm_cartesc_nest_domain_regist_file(dom_id, PARENT_BASENAME, PARENT_PRC_NUM_X, PARENT_PRC_NUM_Y, LATLON_CATALOGUE_FNAME)
offline setup
Definition: scale_comm_cartesC_nest.F90:708
scale_prc_cartesc::prc_2drank
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
Definition: scale_prc_cartesC.F90:45
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
scale_prof
module profiler
Definition: scale_prof.F90:11
scale_file_grads
module file_grads
Definition: scale_file_grads.F90:11
scale_atmos_grid_cartesc_index::jhalo
integer, public jhalo
Definition: scale_atmos_grid_cartesC_index.F90:45
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfxg
center buffer factor (0-1): x, global
Definition: scale_atmos_grid_cartesC.F90:86
scale_atmos_grid_cartesc_index::ja
integer, public ja
Definition: scale_atmos_grid_cartesC_index.F90:49
mod_copytopo
module Copy topography
Definition: mod_copytopo.F90:11
scale_const::const_pi
real(rp), parameter, public const_pi
pi
Definition: scale_const.F90:32
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_fxg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fxg
face coordinate [m]: x, global
Definition: scale_atmos_grid_cartesC.F90:78
mod_copytopo::copytopo_get_size_wrfarw
subroutine copytopo_get_size_wrfarw(IA_org, JA_org)
Definition: mod_copytopo.F90:853
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
scale_atmos_grid_cartesc::atmos_grid_cartesc_cyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cyg
center coordinate [m]: y, global
Definition: scale_atmos_grid_cartesC.F90:77
scale_atmos_grid_cartesc_index::jmax
integer, public jmax
Definition: scale_atmos_grid_cartesC_index.F90:38
scale_atmos_grid_cartesc::atmos_grid_cartesc_cy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cy
center coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:57
scale_atmos_grid_cartesc_index::ieb
integer, public ieb
Definition: scale_atmos_grid_cartesC_index.F90:65
scale_filter
module FILTER
Definition: scale_filter.F90:11
mod_copytopo::copytopo_get_size_grads
subroutine copytopo_get_size_grads(IA_org, JA_org)
Definition: mod_copytopo.F90:741
scale_const::const_d2r
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:33
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:53
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_landuse
module LANDUSE
Definition: scale_landuse.F90:19
scale_atmos_grid_cartesc::atmos_grid_cartesc_cbfyg
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cbfyg
center buffer factor (0-1): y, global
Definition: scale_atmos_grid_cartesC.F90:87
scale_file_grads::file_grads_close
subroutine, public file_grads_close(file_id)
Definition: scale_file_grads.F90:713
scale_atmos_grid_cartesc_index::jsb
integer, public jsb
Definition: scale_atmos_grid_cartesC_index.F90:66
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
scale_atmos_grid_cartesc::atmos_grid_cartesc_cx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cx
center coordinate [m]: x, local
Definition: scale_atmos_grid_cartesC.F90:56
scale_prc::prc_ismaster
logical, public prc_ismaster
master process in local communicator?
Definition: scale_prc.F90:92
mod_copytopo::copytopo_get_data_scale
subroutine copytopo_get_data_scale(IA_org, JA_org, LON_org, LAT_org, TOPO_org)
Definition: mod_copytopo.F90:671