SCALE-RM
mod_copytopo.f90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
9 !-------------------------------------------------------------------------------
11  !-----------------------------------------------------------------------------
12  !
13  !++ used modules
14  !
15  use gtool_file_h
16  use gtool_file, only: &
17  fileread
18  use scale_precision
19  use scale_stdio
20  use scale_prof
22  use scale_tracer
23  !-----------------------------------------------------------------------------
24  implicit none
25  private
26  !-----------------------------------------------------------------------------
27  !
28  !++ Public procedure
29  !
30  public :: copytopo
31 
32  !-----------------------------------------------------------------------------
33  !
34  !++ Public parameters & variables
35  !
36  integer, public :: cnvtopo_type = -1
37 
38  integer, public, parameter :: i_ignore = 0
39  integer, public, parameter :: i_gtopo30 = 1
40  integer, public, parameter :: i_dem50m = 2
41  integer, public, parameter :: i_gmted2010 = 3
42 
43  !-----------------------------------------------------------------------------
44  !
45  !++ Private procedure
46  !
47  private :: copytopo_transgrid
48  private :: copytopo_setalpha
49  private :: copytopo_input_data
50  private :: copytopo_mix_data
51 
52  !-----------------------------------------------------------------------------
53  !
54  !++ Private parameters & variables
55  !
56  integer, private, parameter :: handle = 1
57  integer, private :: itp_nh = 4
58 
59  character(len=H_LONG), private :: copytopo_in_basename = ''
60 
61  real(RP), private :: copytopo_transition_dx = -1.0_rp
62  real(RP), private :: copytopo_transition_dy = -1.0_rp
63  real(RP), private :: copytopo_transfact = -1.0_rp
64  real(RP), private :: copytopo_fracx = 1.0_rp
65  real(RP), private :: copytopo_fracy = 1.0_rp
66  real(RP), private :: copytopo_taux = 1.0_rp
67  real(RP), private :: copytopo_tauy = 1.0_rp
68 
69  logical, private :: copytopo_entire_region = .false.
70  logical, private :: copytopo_linear_h = .true.
71  real(RP), private :: copytopo_exp_h = 2.0_rp
72 
73  real(RP), private, allocatable :: ctrx(:)
74  real(RP), private, allocatable :: ctry(:)
75  real(RP), private, allocatable :: copytopo_alpha(:,:)
76  real(RP), private, allocatable :: topo_pd(:,:)
77 
78  !-----------------------------------------------------------------------------
79 contains
80  !-----------------------------------------------------------------------------
82  subroutine copytopo( &
83  topo_cd ) ! [inout]
84  use scale_process, only: &
86  use scale_grid, only: &
87  dx, &
88  dy, &
89  buffer_dx, &
90  buffer_dy, &
91  bufffact
92  use scale_grid_nest, only: &
94  implicit none
95 
96  real(RP), intent(inout) :: topo_cd(:,:)
97 
98  namelist / param_copytopo / &
99  copytopo_in_basename, &
100  copytopo_transition_dx, &
101  copytopo_transition_dy, &
102  copytopo_transfact, &
103  copytopo_fracx, &
104  copytopo_fracy, &
105  copytopo_taux, &
106  copytopo_tauy, &
107  copytopo_entire_region, &
108  copytopo_linear_h, &
109  copytopo_exp_h
110 
111  integer :: ierr
112  !---------------------------------------------------------------------------
113 
114  if( io_l ) write(io_fid_log,*)
115  if( io_l ) write(io_fid_log,*) '+++ Module[COPYTOPO]/Categ[COPYTOPO]'
116 
117  !--- read namelist
118  rewind(io_fid_conf)
119  read(io_fid_conf,nml=param_copytopo,iostat=ierr)
120  if( ierr < 0 ) then !--- missing
121  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
122  elseif( ierr > 0 ) then !--- fatal error
123  write(*,*) 'xxx Not appropriate names in namelist PARAM_COPYTOPO. Check!'
124  call prc_mpistop
125  endif
126  if( io_lnml ) write(io_fid_log,nml=param_copytopo)
127 
128  if ( copytopo_transition_dx < 0.0_rp ) then
129  copytopo_transition_dx = buffer_dx
130  end if
131  if ( copytopo_transition_dy < 0.0_rp ) then
132  copytopo_transition_dy = buffer_dy
133  end if
134  if ( copytopo_transfact < 0.0_rp ) then
135  copytopo_transfact = bufffact
136  end if
137 
138  allocate( ctrx(ia) )
139  allocate( ctry(ja) )
140  allocate( copytopo_alpha(ia,ja) )
141  allocate( topo_pd(ia,ja) )
142  copytopo_alpha(:,:) = 0.0_rp
143 
144  itp_nh = int( nest_interp_level )
145 
146  ! copy topography from parent domain to transition region
147  call copytopo_transgrid
148 
149  call copytopo_setalpha
150 
151  call copytopo_input_data( topo_pd ) ! (out)
152 
153  call copytopo_mix_data( topo_cd, & ! (inout)
154  topo_pd ) ! (in)
155 
156  return
157  end subroutine copytopo
158 
159  !-----------------------------------------------------------------------------
161  !-----------------------------------------------------------------------------
162 
163  !-----------------------------------------------------------------------------
165  subroutine copytopo_transgrid
166  use scale_process, only: &
167  prc_mpistop, &
168  prc_myrank
169  use scale_rm_process, only: &
170  prc_2drank, &
171  prc_num_x, &
172  prc_num_y
173  use scale_const, only: &
174  radius => const_radius, &
175  eps => const_eps
176  use scale_grid, only: &
177  dx, &
178  dy, &
179  cxg => grid_cxg, &
180  fxg => grid_fxg, &
181  cyg => grid_cyg, &
182  fyg => grid_fyg, &
183  cbfxg => grid_cbfxg, &
184  cbfyg => grid_cbfyg, &
185  bufffact
186  implicit none
187 
188  real(RP), allocatable :: CTRXG(:)
189  real(RP), allocatable :: CTRYG(:)
190 
191  real(RP), allocatable :: buffx(:), buffy(:)
192  real(RP) :: bufftotx, bufftoty
193  real(RP), allocatable :: transx(:), transy(:)
194  real(RP) :: transtotx, transtoty
195 
196  integer :: IAG, JAG
197  integer :: imain, ibuff, itrans
198  integer :: jmain, jbuff, jtrans
199  integer :: copy_is, copy_ie, copy_js, copy_je
200  integer :: i, j, ii, jj
201  !---------------------------------------------------------------------------
202 
203  ! array size (global domain)
204  iag = ihalo + imax*prc_num_x + ihalo
205  jag = jhalo + jmax*prc_num_y + jhalo
206 
207  allocate( buffx(0:iag) )
208  allocate( buffy(0:jag) )
209  allocate( transx(0:iag) )
210  allocate( transy(0:jag) )
211  allocate( ctrxg( iag) )
212  allocate( ctryg( jag) )
213 
214  ! X-direction
215  ! calculate buffer grid size
216  buffx(:) = dx
217  transx(0) = dx
218  bufftotx = 0.0_rp
219  transtotx = 0.0_rp
220 
221  do i = ihalo+1, iag
222  if( abs(cbfxg(i) - 0.0_rp) < eps ) exit
223  buffx(i) = buffx(i-1) * bufffact
224  bufftotx = bufftotx + buffx(i)
225  enddo
226  ibuff = i - (ihalo+1)
227 
228  do i = 1, iag
229  if( transtotx >= copytopo_transition_dx ) exit
230  transx(i) = transx(i-1) * copytopo_transfact
231  transtotx = transtotx + transx(i)
232  enddo
233  itrans = i - 1
234  imain = iag - 2*ibuff - 2*itrans - 2*ihalo
235 
236  if ( imain < 1 ) then
237  write(*,*) 'xxx Not appropriate transition width for global domain(X).', copytopo_transition_dx
238  write(*,*) ' # of buffer region (one side)', ibuff
239  write(*,*) ' # of transion region (one side)', itrans
240  call prc_mpistop
241  endif
242 
243  ! calc transition factor (global domaim)
244  ctrxg(:) = 0.0_rp
245  do i = 1, ihalo+ibuff
246  ctrxg(i) = 1.0_rp
247  enddo
248 
249  if ( itrans > 0 ) then
250  copy_is = ihalo+ibuff+1
251  copy_ie = ihalo+ibuff+itrans
252  do i = copy_is, copy_ie
253  ctrxg(i) = (transtotx+bufftotx+fxg(ihalo )-cxg(i)) / transtotx
254  enddo
255  copy_is = ihalo+ibuff+itrans+imain+1
256  copy_ie = ihalo+ibuff+itrans+imain+itrans+ibuff
257  do i = copy_is, copy_ie
258  ctrxg(i) = (transtotx+bufftotx-fxg(iag-ihalo)+cxg(i)) / transtotx
259  enddo
260  endif
261 
262  copy_is = ihalo+ibuff+itrans+imain+itrans+ibuff+1
263  copy_ie = ihalo+ibuff+itrans+imain+itrans+ibuff+ihalo
264  do i = copy_is, copy_ie
265  ctrxg(i) = 1.0_rp
266  enddo
267  ctrxg(:) = max( min( ctrxg(:), 1.0_rp ), 0.0_rp )
268 
269  ! Y-direction
270  ! calculate buffer grid size
271  buffy(:) = dy
272  transy(0) = dy
273  bufftoty = 0.0_rp
274  transtoty = 0.0_rp
275 
276  do j = jhalo+1, jag
277  if( abs(cbfyg(j) - 0.0_rp) < eps ) exit
278  buffy(j) = buffy(j-1) * bufffact
279  bufftoty = bufftoty + buffy(j)
280  enddo
281  jbuff = j - (jhalo+1)
282 
283  do j = 1, jag
284  if( transtoty >= copytopo_transition_dy ) exit
285  transy(j) = transy(j-1) * copytopo_transfact
286  transtoty = transtoty + transy(j)
287  enddo
288  jtrans = j - 1
289  jmain = jag - 2*jbuff - 2*jtrans - 2*jhalo
290 
291  if ( jmain < 1 ) then
292  write(*,*) 'xxx Not appropriate transition width for global domain(Y).', copytopo_transition_dy
293  write(*,*) ' # of buffer region (one side)', jbuff
294  write(*,*) ' # of transion region (one side)', jtrans
295  call prc_mpistop
296  endif
297 
298  ! calc transition factor (global domaim)
299  ctryg(:) = 0.0_rp
300  do j = 1, jhalo+jbuff
301  ctryg(j) = 1.0_rp
302  enddo
303 
304  if ( jtrans > 0 ) then
305  copy_js = jhalo+jbuff+1
306  copy_je = jhalo+jbuff+jtrans
307  do j = copy_js, copy_je
308  ctryg(j) = (transtoty+bufftoty+fyg(jhalo )-cyg(j)) / transtoty
309  enddo
310  copy_js = jhalo+jbuff+jtrans+jmain+1
311  copy_je = jhalo+jbuff+jtrans+jmain+jtrans+jbuff
312  do j = copy_js, copy_je
313  ctryg(j) = (transtoty+bufftoty-fyg(jag-jhalo)+cyg(j)) / transtoty
314  enddo
315  endif
316 
317  copy_js = jhalo+jbuff+jtrans+jmain+jtrans+jbuff+1
318  copy_je = jhalo+jbuff+jtrans+jmain+jtrans+jbuff+jhalo
319  do j = copy_js, copy_je
320  ctryg(j) = 1.0_rp
321  enddo
322  ctryg(:) = max( min( ctryg(:), 1.0_rp ), 0.0_rp )
323 
324  ! horizontal coordinate (local domaim)
325  do i = 1, ia
326  ii = i + prc_2drank(prc_myrank,1) * imax
327  ctrx(i) = ctrxg(ii)
328  enddo
329  do j = 1, ja
330  jj = j + prc_2drank(prc_myrank,2) * jmax
331  ctry(j) = ctryg(jj)
332  enddo
333 
334  deallocate( transx )
335  deallocate( transy )
336 
337  return
338  end subroutine copytopo_transgrid
339 
340 
341  !-----------------------------------------------------------------------------
343  subroutine copytopo_setalpha
344  use scale_const, only: &
345  eps => const_eps
346  use scale_comm, only: &
347  comm_vars8, &
348  comm_wait
349 
350  real(RP) :: coef_x, alpha_x1
351  real(RP) :: coef_y, alpha_y1
352  real(RP) :: ee1
353 
354  integer :: i, j
355  !---------------------------------------------------------------------------
356 
357  ! check invalid fraction
358  copytopo_fracx = max( min( copytopo_fracx, 1.0_rp ), eps )
359  copytopo_fracy = max( min( copytopo_fracy, 1.0_rp ), eps )
360 
361  if ( copytopo_taux <= 0.0_rp ) then ! invalid tau
362  coef_x = 0.0_rp
363  else
364  coef_x = 1.0_rp / copytopo_taux
365  endif
366 
367  if ( copytopo_tauy <= 0.0_rp ) then ! invalid tau
368  coef_y = 0.0_rp
369  else
370  coef_y = 1.0_rp / copytopo_tauy
371  endif
372 
373  do j = 1, ja
374  do i = 1, ia
375  ee1 = ctrx(i)
376  if ( ee1 <= 1.0_rp - copytopo_fracx ) then
377  ee1 = 0.0_rp
378  else
379  ee1 = ( ee1 - 1.0_rp + copytopo_fracx ) / copytopo_fracx
380  endif
381 
382  if ( copytopo_linear_h ) then
383  alpha_x1 = coef_x * ee1
384  else
385  alpha_x1 = coef_x * ee1 * exp( -(1.0_rp-ee1) * copytopo_exp_h )
386  end if
387 
388  ee1 = ctry(j)
389  if ( ee1 <= 1.0_rp - copytopo_fracy ) then
390  ee1 = 0.0_rp
391  else
392  ee1 = ( ee1 - 1.0_rp + copytopo_fracy ) / copytopo_fracy
393  endif
394 
395  if ( copytopo_linear_h ) then
396  alpha_y1 = coef_y * ee1
397  else
398  alpha_y1 = coef_y * ee1 * exp( -(1.0_rp-ee1) * copytopo_exp_h )
399  end if
400 
401  copytopo_alpha(i,j) = max( alpha_x1, alpha_y1 )
402  enddo
403  enddo
404 
405  call comm_vars8( copytopo_alpha(:,:), 1 )
406  call comm_wait ( copytopo_alpha(:,:), 1 )
407 
408  return
409  end subroutine copytopo_setalpha
410 
411 
412  !-----------------------------------------------------------------------------
414  subroutine copytopo_input_data( &
415  topo_pd ) ! (out)
416  use scale_const, only: &
417  d2r => const_d2r
418  use scale_comm, only: &
419  comm_vars8, &
420  comm_wait
421  use scale_grid_nest, only: &
422  parent_imax, &
423  parent_jmax, &
424  nest_tile_num_x, &
425  nest_tile_num_y, &
426  nest_tile_id, &
428  use scale_interpolation_nest, only: &
432  use scale_grid_real, only: &
433  lat => real_lat, &
434  lon => real_lon
435  implicit none
436 
437  real(RP), intent(out) :: topo_pd(:,:)
438 
439  real(RP) :: dummy (1,1,1)
440  real(RP), allocatable :: read2D(:,:)
441  real(RP), allocatable :: lon_org (:,:)
442  real(RP), allocatable :: lat_org (:,:)
443  real(RP), allocatable :: topo_org(:,:)
444  real(RP), allocatable :: hfact(:,:,:)
445  integer, allocatable :: igrd (:,:,:)
446  integer, allocatable :: jgrd (:,:,:)
447 
448  integer :: IALL, JALL ! number of grids for whole domain
449  integer :: PTI, PTJ ! number of grids for a tile
450  integer :: tilei, tilej
451  integer :: rank
452  integer :: i
453  integer :: cxs, cxe, cys, cye ! for child domain
454  integer :: pxs, pxe, pys, pye ! for parent domain
455  !---------------------------------------------------------------------------
456 
457  pti = parent_imax(handle)
458  ptj = parent_jmax(handle)
459  iall = pti * nest_tile_num_x
460  jall = ptj * nest_tile_num_y
461 
462  allocate( hfact( ia, ja, itp_nh ) )
463  allocate( igrd( ia, ja, itp_nh ) )
464  allocate( jgrd( ia, ja, itp_nh ) )
465  allocate( lon_org( iall, jall ) )
466  allocate( lat_org( iall, jall ) )
467  allocate( topo_org( iall, jall ) )
468 
469  do i = 1, size( nest_tile_id(:) )
470  ! read data from split files
471  rank = nest_tile_id(i)
472 
473  call nest_domain_shape ( tilei, tilej, & ! [out]
474  cxs, cxe, & ! [out]
475  cys, cye, & ! [out]
476  pxs, pxe, & ! [out]
477  pys, pye, & ! [out]
478  i ) ! [in ]
479  allocate( read2d( tilei,tilej ) )
480 
481  call fileread( read2d(:,:), copytopo_in_basename, "lon", 1, rank )
482  lon_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
483  call fileread( read2d(:,:), copytopo_in_basename, "lat", 1, rank )
484  lat_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye) * d2r
485  call fileread( read2d(:,:), copytopo_in_basename, "TOPO", 1, rank )
486  topo_org(cxs:cxe,cys:cye) = read2d(pxs:pxe,pys:pye)
487 
488  deallocate( read2d )
489  end do
490 
491  call intrpnest_domain_compatibility( lon_org(:,:), lat_org(:,:), dummy(:,:,:), &
492  lon(:,:), lat(:,:), dummy(:,:,:), &
493  skip_z=.true. )
494  !call check_domain_compatibility( lon_org(:,:), lat_org(:,:), dummy(:,:,:), &
495  ! LON(:,:), LAT(:,:), dummy(:,:,:), &
496  ! skip_z=.true. )
497 
498  call intrpnest_interp_fact_latlon( hfact(:,:,:), & ! (out)
499  igrd(:,:,:), & ! (out)
500  jgrd(:,:,:), & ! (out)
501  lat(:,:), & ! (in)
502  lon(:,:), & ! (in)
503  ia, ja, & ! (in)
504  lat_org(:,:), & ! (in)
505  lon_org(:,:), & ! (in)
506  iall, jall ) ! (in)
507 
508  call intrpnest_interp_2d( topo_pd(:,:), &
509  topo_org(:,:), &
510  hfact(:,:,:), &
511  igrd(:,:,:), &
512  jgrd(:,:,:), &
513  ia, ja )
514 
515  call comm_vars8( topo_pd(:,:), 1 )
516  call comm_wait ( topo_pd(:,:), 1 )
517 
518  return
519  end subroutine copytopo_input_data
520 
521  !-----------------------------------------------------------------------------
523  subroutine copytopo_mix_data( &
524  topo_cd, & ! (inout)
525  topo_pd ) ! (in)
526  implicit none
527  real(RP), intent(inout) :: topo_cd(:,:) ! topography of child domain (mine)
528  real(RP), intent(in) :: topo_pd(:,:) ! topography of parent domain
529 
530  real(RP) :: frac
531  integer :: i, j
532  !---------------------------------------------------------------------------
533 
534  if ( copytopo_entire_region ) then
535  topo_cd(:,:) = topo_pd(:,:)
536  else
537  do j = 1, ja
538  do i = 1, ia
539  frac = copytopo_alpha(i,j)
540  topo_cd(i,j) = topo_cd(i,j) * ( 1.0_rp - frac ) &
541  + topo_pd(i,j) * frac
542  end do
543  end do
544  endif
545 
546  return
547  end subroutine copytopo_mix_data
548 
549 end module mod_copytopo
integer, public imax
of computational cells: x
integer, public prc_num_x
x length of 2D processor topology
real(rp), public buffer_dx
thickness of buffer region [m]: x
module GTOOL_FILE
Definition: gtool_file.f90:17
subroutine, public intrpnest_domain_compatibility(lon_org, lat_org, lev_org, lon_loc, lat_loc, lev_loc, skip_x, skip_y, skip_z)
real(rp), public dy
length in the main region [m]: y
real(rp), dimension(:), allocatable, public grid_cyg
center coordinate [m]: y, global
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:), allocatable, public grid_cbfyg
center buffer factor [0-1]: y, global
integer, parameter, public i_dem50m
integer, dimension(2), public parent_jmax
parent max number in y-direction
module GRID (nesting system)
real(rp), dimension(:), allocatable, public grid_cxg
center coordinate [m]: x, global
real(rp), public dx
length in the main region [m]: x
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:), allocatable, public grid_fxg
face coordinate [m]: x, global
real(rp), public const_radius
radius of the planet [m]
Definition: scale_const.F90:46
module STDIO
Definition: scale_stdio.F90:12
real(rp), public const_d2r
degree to radian
Definition: scale_const.F90:35
procedure(intrpnest_intfc_interp_2d), pointer, public intrpnest_interp_2d
integer, dimension(:), allocatable, public nest_tile_id
parent tile real id
integer, public nest_interp_level
horizontal interpolation level
subroutine, public intrpnest_interp_fact_latlon(hfact, igrd, jgrd, mylat, mylon, myIA, myJA, inlat, inlon, inIA, inJA)
module Copy topography
integer, parameter, public i_gtopo30
integer, public prc_num_y
y length of 2D processor topology
module grid index
module TRACER
subroutine, public nest_domain_shape(tilei, tilej, cxs, cxe, cys, cye, pxs, pxe, pys, pye, iloc)
Return shape of ParentDomain at the specified rank (for offline)
integer, public ia
of x whole cells (local, with HALO)
module GRID (real space)
integer, parameter, public i_gmted2010
integer, parameter, public i_ignore
integer, public jhalo
of halo cells: y
module COMMUNICATION
Definition: scale_comm.F90:23
module PROCESS
integer, public nest_tile_num_y
parent tile number in y-direction
module CONSTANT
Definition: scale_const.F90:14
integer, public prc_myrank
process num in local communicator
module GRID (cartesian)
module INTERPOLATION (nesting system)
module RM PROCESS
subroutine, public copytopo(topo_cd)
Setup and Main.
real(rp), public buffer_dy
thickness of buffer region [m]: y
module profiler
Definition: scale_prof.F90:10
integer, dimension(2), public parent_imax
parent max number in x-direction
real(rp), public const_eps
small number
Definition: scale_const.F90:36
logical, public io_lnml
output log or not? (for namelist, this process)
Definition: scale_stdio.F90:60
real(rp), dimension(:,:), allocatable, public real_lon
longitude [rad,0-2pi]
integer, public nest_tile_num_x
parent tile number in x-direction
integer, public cnvtopo_type
module PRECISION
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
module FILE I/O HEADER
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:,:), allocatable, public real_lat
latitude [rad,-pi,pi]
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
real(rp), public bufffact
default strech factor for dx/dy/dz of buffer region
integer, public jmax
of computational cells: y
real(rp), dimension(:), allocatable, public grid_fyg
face coordinate [m]: y, global
real(rp), dimension(:), allocatable, public grid_cbfxg
center buffer factor [0-1]: x, global
integer, public ihalo
of halo cells: x
integer, public ja
of y whole cells (local, with HALO)