SCALE-RM
scale_gridtrans.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
14 !-------------------------------------------------------------------------------
16  !-----------------------------------------------------------------------------
17  !
18  !++ used modules
19  !
20  use scale_precision
21  use scale_stdio
22  use scale_prof
24  !-----------------------------------------------------------------------------
25  implicit none
26  private
27  !-----------------------------------------------------------------------------
28  !
29  !++ Public procedure
30  !
31  public :: gtrans_setup
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  real(RP), public, allocatable :: gtrans_mapf (:,:,:,:)
38  real(RP), public, allocatable :: gtrans_rotc (:,:,:)
39 
40  real(RP), public, allocatable :: gtrans_gsqrt(:,:,:,:)
41  real(RP), public, allocatable :: gtrans_j13g (:,:,:,:)
42  real(RP), public, allocatable :: gtrans_j23g (:,:,:,:)
43  real(RP), public :: gtrans_j33g
44 
45  real(RP), public, allocatable :: gtrans_limyz(:,:,:,:)
46  real(RP), public, allocatable :: gtrans_limxz(:,:,:,:)
47  real(RP), public, allocatable :: gtrans_limxy(:,:,:,:)
48 
49  integer, public :: i_xyz = 1 ! at (x,y,z)
50  integer, public :: i_xyw = 2 ! at (x,y,w)
51  integer, public :: i_uyw = 3 ! at (u,y,w)
52  integer, public :: i_xvw = 4 ! at (x,v,w)
53  integer, public :: i_uyz = 5 ! at (u,y,z)
54  integer, public :: i_xvz = 6 ! at (x,v,z)
55  integer, public :: i_uvz = 7 ! at (u,v,z)
56 
57  integer, public :: i_xy = 1 ! at (x,y)
58  integer, public :: i_uy = 2 ! at (u,y)
59  integer, public :: i_xv = 3 ! at (x,v)
60  integer, public :: i_uv = 4 ! at (u,v)
61 
62  integer, public :: i_fyz = 1 ! y-z face limiting x-flux
63  integer, public :: i_fxz = 2 ! x-z face limiting y-flux
64  integer, public :: i_fxy = 3 ! x-y face limiting z-flux
65 
66  !-----------------------------------------------------------------------------
67  !
68  !++ Private procedure
69  !
70  private :: gtrans_mapfactor
71  private :: gtrans_terrainfollowing
72  private :: gtrans_thin_wall
73  private :: gtrans_step_mountain
74  private :: gtrans_write
75 
76  !-----------------------------------------------------------------------------
77  !
78  !++ Private parameters & variables
79  !
80  character(len=H_LONG), private :: gtrans_out_basename = ''
81  character(len=H_MID), private :: gtrans_out_title = 'SCALE-RM GEOMETRICS'
82  character(len=H_MID), private :: gtrans_out_dtype = 'DEFAULT'
83 
84  character(len=H_MID), private :: gtrans_topo_type = 'TERRAINFOLLOWING'
85  integer, private :: gtrans_thinwall_xdiv = 50
86  integer, private :: gtrans_thinwall_ydiv = 50
87  logical, private :: debug = .false.
88  !-----------------------------------------------------------------------------
89 contains
90  !-----------------------------------------------------------------------------
92  subroutine gtrans_setup
93  use scale_process, only: &
95  implicit none
96 
97  namelist / param_gtrans / &
98  gtrans_out_basename, &
99  gtrans_out_dtype, &
100  gtrans_topo_type, &
101  gtrans_thinwall_xdiv, &
102  gtrans_thinwall_ydiv, &
103  debug
104 
105  integer :: ierr
106  !---------------------------------------------------------------------------
107 
108  if( io_l ) write(io_fid_log,*)
109  if( io_l ) write(io_fid_log,*) '++++++ Module[GRIDTRANS] / Categ[ATMOS-RM GRID] / Origin[SCALElib]'
110 
111  !--- read namelist
112  rewind(io_fid_conf)
113  read(io_fid_conf,nml=param_gtrans,iostat=ierr)
114  if( ierr < 0 ) then !--- missing
115  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
116  elseif( ierr > 0 ) then !--- fatal error
117  write(*,*) 'xxx Not appropriate names in namelist PARAM_GTRANS. Check!'
118  call prc_mpistop
119  endif
120  if( io_lnml ) write(io_fid_log,nml=param_gtrans)
121 
122  allocate( gtrans_mapf(ia,ja,2,4) )
123 
124  allocate( gtrans_rotc(ia,ja,2) )
125 
126  allocate( gtrans_gsqrt(ka,ia,ja,7) )
127  allocate( gtrans_j13g(ka,ia,ja,7) )
128  allocate( gtrans_j23g(ka,ia,ja,7) )
129 
130  gtrans_gsqrt(:,:,:,:) = 1.0_rp
131  gtrans_j13g(:,:,:,:) = 0.0_rp
132  gtrans_j23g(:,:,:,:) = 0.0_rp
133  gtrans_j33g = 1.0_rp
134 
135  allocate( gtrans_limyz(ka,ia,ja,7) )
136  allocate( gtrans_limxz(ka,ia,ja,7) )
137  allocate( gtrans_limxy(ka,ia,ja,7) )
138  gtrans_limyz(:,:,:,:) = 1.0_rp
139  gtrans_limxz(:,:,:,:) = 1.0_rp
140  gtrans_limxy(:,:,:,:) = 1.0_rp
141 
142  ! calc metrics for orthogonal curvelinear coordinate
143  call gtrans_mapfactor
144 
145  ! calc coeficient for rotaion of velocity vector
146  call gtrans_rotcoef
147 
148  ! calc metrics for terrain-following,step-mountain,thin-wall coordinate
149  select case(gtrans_topo_type)
150  case ('TERRAINFOLLOWING')
151  if( io_l ) write(io_fid_log,*) '=> Use terrain-following coordinate'
152  call gtrans_terrainfollowing
153  case ('STEPMOUNTAIN')
154  if( io_l ) write(io_fid_log,*) '=> Use step mountain method'
155  call gtrans_thin_wall
156  call gtrans_step_mountain
157  case ('THINWALL')
158  if( io_l ) write(io_fid_log,*) '=> Use thin-wall approximation'
159  call gtrans_thin_wall
160  case default
161  write(*,*) 'xxx Not appropriate name for GTRANS_TOPO_TYPE : ', trim(gtrans_topo_type)
162  call prc_mpistop
163  end select
164 
165  ! output metrics (for debug)
166  call gtrans_write
167 
168  return
169  end subroutine gtrans_setup
170 
171  !-----------------------------------------------------------------------------
173  subroutine gtrans_mapfactor
174  use scale_mapproj, only: &
176  use scale_grid_real, only: &
178  real_lat, &
179  real_latx, &
180  real_laty, &
181  real_latxy
182  implicit none
183  !---------------------------------------------------------------------------
184 
185  call mprj_mapfactor( real_lat , gtrans_mapf(:,:,1,i_xy), gtrans_mapf(:,:,2,i_xy))
186  call mprj_mapfactor( real_latx , gtrans_mapf(:,:,1,i_uy), gtrans_mapf(:,:,2,i_uy))
187  call mprj_mapfactor( real_laty , gtrans_mapf(:,:,1,i_xv), gtrans_mapf(:,:,2,i_xv))
189 
190  call real_calc_areavol( gtrans_mapf(:,:,:,i_xy) )
191 
192  return
193  end subroutine gtrans_mapfactor
194 
195  !-----------------------------------------------------------------------------
197  subroutine gtrans_rotcoef
198  use scale_mapproj, only: &
199  mprj_rotcoef
200  use scale_grid_real, only: &
201  real_lon, &
202  real_lat
203  implicit none
204  !---------------------------------------------------------------------------
205 
206  call mprj_rotcoef( gtrans_rotc(:,:,:), & ! [OUT]
207  real_lon(:,:), & ! [IN]
208  real_lat(:,:) ) ! [IN]
209 
210  return
211  end subroutine gtrans_rotcoef
212 
213  !-----------------------------------------------------------------------------
215  subroutine gtrans_terrainfollowing
216  use scale_grid, only: &
217  grid_rcdz, &
218  grid_rcdx, &
219  grid_rcdy, &
220  grid_rfdz, &
221  grid_rfdx, &
222  grid_rfdy
223  use scale_grid_real, only: &
224  real_cz, &
225  real_fz
226  use scale_comm, only: &
227  comm_vars8, &
228  comm_wait
229  implicit none
230 
231  real(RP) :: REAL_CZ_U ( ka,ia,ja)
232  real(RP) :: REAL_CZ_V ( ka,ia,ja)
233  real(RP) :: REAL_CZ_UV( ka,ia,ja)
234  real(RP) :: REAL_FZ_U (0:ka,ia,ja)
235  real(RP) :: REAL_FZ_V (0:ka,ia,ja)
236  real(RP) :: REAL_FZ_UV(0:ka,ia,ja)
237 
238  integer :: k, i, j
239  !---------------------------------------------------------------------------
240 
241  ! calc Z-coordinate height at staggered position
242  do j = 1, ja
243  do i = 1, ia-1
244  do k = 1, ka
245  real_cz_u(k,i,j) = 0.5_rp * ( real_cz(k,i+1,j) + real_cz(k,i,j) )
246  enddo
247  enddo
248  enddo
249 
250  do j = 1, ja
251  do i = 1, ia-1
252  do k = 0, ka
253  real_fz_u(k,i,j) = 0.5_rp * ( real_fz(k,i+1,j) + real_fz(k,i,j) )
254  enddo
255  enddo
256  enddo
257 
258  do j = 1, ja-1
259  do i = 1, ia
260  do k = 1, ka
261  real_cz_v(k,i,j) = 0.5_rp * ( real_cz(k,i,j+1) + real_cz(k,i,j) )
262  enddo
263  enddo
264  enddo
265 
266  do j = 1, ja-1
267  do i = 1, ia
268  do k = 0, ka
269  real_fz_v(k,i,j) = 0.5_rp * ( real_fz(k,i,j+1) + real_fz(k,i,j) )
270  enddo
271  enddo
272  enddo
273 
274  do j = 1, ja-1
275  do i = 1, ia-1
276  do k = 1, ka
277  real_cz_uv(k,i,j) = 0.25_rp * ( real_cz(k,i+1,j+1) + real_cz(k,i+1,j) &
278  + real_cz(k,i ,j+1) + real_cz(k,i ,j) )
279  enddo
280  enddo
281  enddo
282 
283  do j = 1, ja-1
284  do i = 1, ia-1
285  do k = 0, ka
286  real_fz_uv(k,i,j) = 0.25_rp * ( real_fz(k,i+1,j+1) + real_fz(k,i+1,j) &
287  + real_fz(k,i ,j+1) + real_fz(k,i ,j) )
288  enddo
289  enddo
290  enddo
291 
292  ! G^1/2
293  do j = js, je
294  do i = is, ie
295  ! at (x,y,z)
296  do k = 1, ka
297  gtrans_gsqrt(k,i,j,i_xyz) = ( real_fz(k,i,j) - real_fz(k-1,i,j) ) * grid_rcdz(k)
298  enddo
299 
300  ! at (x,y,w)
301  do k = 1, ka-1
302  gtrans_gsqrt(k,i,j,i_xyw) = ( real_cz(k+1,i,j) - real_cz(k,i,j) ) * grid_rfdz(k)
303  enddo
304  gtrans_gsqrt(ka,i,j,i_xyw) = gtrans_gsqrt(ka-1,i,j,i_xyw)
305 
306  ! at (u,y,w)
307  do k = 1, ka-1
308  gtrans_gsqrt(k,i,j,i_uyw) = ( real_cz_u(k+1,i,j) - real_cz_u(k,i,j) ) * grid_rfdz(k)
309  enddo
310  gtrans_gsqrt(ka,i,j,i_uyw) = gtrans_gsqrt(ka-1,i,j,i_uyw)
311 
312  ! at (x,v,w)
313  do k = 1, ka-1
314  gtrans_gsqrt(k,i,j,i_xvw) = ( real_cz_v(k+1,i,j) - real_cz_v(k,i,j) ) * grid_rfdz(k)
315  enddo
316  gtrans_gsqrt(ka,i,j,i_xvw) = gtrans_gsqrt(ka-1,i,j,i_xvw)
317 
318  ! at (u,y,z)
319  do k = 1, ka
320  gtrans_gsqrt(k,i,j,i_uyz) = ( real_fz_u(k,i,j) - real_fz_u(k-1,i,j) ) * grid_rcdz(k)
321  enddo
322 
323  ! at (x,v,z)
324  do k = 1, ka
325  gtrans_gsqrt(k,i,j,i_xvz) = ( real_fz_v(k,i,j) - real_fz_v(k-1,i,j) ) * grid_rcdz(k)
326  enddo
327 
328  ! at (u,v,z)
329  do k = 1, ka
330  gtrans_gsqrt(k,i,j,i_uvz) = ( real_fz_uv(k,i,j) - real_fz_uv(k-1,i,j) ) * grid_rcdz(k)
331  enddo
332  enddo
333  enddo
334 
335  call comm_vars8( gtrans_gsqrt(:,:,:,1), 1 )
336  call comm_vars8( gtrans_gsqrt(:,:,:,2), 2 )
337  call comm_vars8( gtrans_gsqrt(:,:,:,3), 3 )
338  call comm_vars8( gtrans_gsqrt(:,:,:,4), 4 )
339  call comm_vars8( gtrans_gsqrt(:,:,:,5), 5 )
340  call comm_vars8( gtrans_gsqrt(:,:,:,6), 6 )
341  call comm_vars8( gtrans_gsqrt(:,:,:,7), 7 )
342  call comm_wait ( gtrans_gsqrt(:,:,:,1), 1 )
343  call comm_wait ( gtrans_gsqrt(:,:,:,2), 2 )
344  call comm_wait ( gtrans_gsqrt(:,:,:,3), 3 )
345  call comm_wait ( gtrans_gsqrt(:,:,:,4), 4 )
346  call comm_wait ( gtrans_gsqrt(:,:,:,5), 5 )
347  call comm_wait ( gtrans_gsqrt(:,:,:,6), 6 )
348  call comm_wait ( gtrans_gsqrt(:,:,:,7), 7 )
349 
350  ! Jacobian * G^1/2
351  do j = js, je
352  do i = is, ie
353  do k = 1, ka
354  gtrans_j13g(k,i,j,i_xyz) = -( real_cz_u(k,i ,j) - real_cz_u(k,i-1,j) ) * grid_rcdx(i)
355  gtrans_j13g(k,i,j,i_xyw) = -( real_fz_u(k,i ,j) - real_fz_u(k,i-1,j) ) * grid_rcdx(i)
356  gtrans_j13g(k,i,j,i_uyw) = -( real_fz(k,i+1,j) - real_fz(k,i ,j) ) * grid_rfdx(i)
357  gtrans_j13g(k,i,j,i_xvw) = -( real_fz_uv(k,i ,j) - real_fz_uv(k,i-1,j) ) * grid_rcdx(i)
358  gtrans_j13g(k,i,j,i_uyz) = -( real_cz(k,i ,j) - real_cz(k,i-1,j) ) * grid_rfdx(i)
359  gtrans_j13g(k,i,j,i_xvz) = -( real_cz_uv(k,i ,j) - real_cz_uv(k,i-1,j) ) * grid_rcdx(i)
360  gtrans_j13g(k,i,j,i_uvz) = -( real_cz_v(k,i ,j) - real_cz_v(k,i-1,j) ) * grid_rfdx(i)
361  enddo
362  enddo
363  enddo
364 
365  do j = js, je
366  do i = is, ie
367  do k = 1, ka
368  gtrans_j23g(k,i,j,i_xyz) = -( real_cz_v(k,i,j ) - real_cz_v(k,i,j-1) ) * grid_rcdy(j)
369  gtrans_j23g(k,i,j,i_xyw) = -( real_fz_v(k,i,j ) - real_fz_v(k,i,j-1) ) * grid_rcdy(j)
370  gtrans_j23g(k,i,j,i_uyw) = -( real_fz(k,i,j+1) - real_fz(k,i,j ) ) * grid_rfdy(j)
371  gtrans_j23g(k,i,j,i_xvw) = -( real_fz_uv(k,i,j ) - real_fz_uv(k,i,j-1) ) * grid_rcdy(j)
372  gtrans_j23g(k,i,j,i_uyz) = -( real_cz_uv(k,i,j ) - real_cz_uv(k,i,j-1) ) * grid_rcdy(j)
373  gtrans_j23g(k,i,j,i_xvz) = -( real_cz(k,i,j ) - real_cz(k,i,j-1) ) * grid_rfdy(j)
374  gtrans_j23g(k,i,j,i_uvz) = -( real_cz_u(k,i,j ) - real_cz_u(k,i,j-1) ) * grid_rfdy(j)
375  enddo
376  enddo
377  enddo
378 
379  gtrans_j33g = 1.0_rp ! - 1 / G^1/2 * G^1/2
380 
381  call comm_vars8( gtrans_j13g(:,:,:,i_xyz), 8 )
382  call comm_vars8( gtrans_j13g(:,:,:,i_xyw), 9 )
383  call comm_vars8( gtrans_j13g(:,:,:,i_uyw), 10 )
384  call comm_vars8( gtrans_j13g(:,:,:,i_xvw), 11 )
385  call comm_vars8( gtrans_j13g(:,:,:,i_uyz), 12 )
386  call comm_vars8( gtrans_j13g(:,:,:,i_xvz), 13 )
387  call comm_vars8( gtrans_j13g(:,:,:,i_uvz), 14 )
388  call comm_wait ( gtrans_j13g(:,:,:,i_xyz), 8 )
389  call comm_wait ( gtrans_j13g(:,:,:,i_xyw), 9 )
390  call comm_wait ( gtrans_j13g(:,:,:,i_uyw), 10 )
391  call comm_wait ( gtrans_j13g(:,:,:,i_xvw), 11 )
392  call comm_wait ( gtrans_j13g(:,:,:,i_uyz), 12 )
393  call comm_wait ( gtrans_j13g(:,:,:,i_xvz), 13 )
394  call comm_wait ( gtrans_j13g(:,:,:,i_uvz), 14 )
395 
396  call comm_vars8( gtrans_j23g(:,:,:,i_xyz), 15 )
397  call comm_vars8( gtrans_j23g(:,:,:,i_xyw), 16 )
398  call comm_vars8( gtrans_j23g(:,:,:,i_uyw), 17 )
399  call comm_vars8( gtrans_j23g(:,:,:,i_xvw), 18 )
400  call comm_vars8( gtrans_j23g(:,:,:,i_uyz), 19 )
401  call comm_vars8( gtrans_j23g(:,:,:,i_xvz), 20 )
402  call comm_vars8( gtrans_j23g(:,:,:,i_uvz), 21 )
403  call comm_wait ( gtrans_j23g(:,:,:,i_xyz), 15 )
404  call comm_wait ( gtrans_j23g(:,:,:,i_xyw), 16 )
405  call comm_wait ( gtrans_j23g(:,:,:,i_uyw), 17 )
406  call comm_wait ( gtrans_j23g(:,:,:,i_xvw), 18 )
407  call comm_wait ( gtrans_j23g(:,:,:,i_uyz), 19 )
408  call comm_wait ( gtrans_j23g(:,:,:,i_xvz), 20 )
409  call comm_wait ( gtrans_j23g(:,:,:,i_uvz), 21 )
410 
411  return
412  end subroutine gtrans_terrainfollowing
413 
414  !-----------------------------------------------------------------------------
415  subroutine gtrans_thin_wall
416  use scale_process, only: &
418  use scale_grid, only: &
419  grid_cz, &
420  grid_cx, &
421  grid_cy, &
422  grid_fz, &
423  grid_fx, &
424  grid_fy
425  use scale_topography, only : &
426  topo_zsfc
427  use scale_comm, only : &
428  comm_vars8, &
429  comm_wait
430  implicit none
431 
432  real(RP) :: TOPO_ZsfcALL(2*ia,2*ja)
433  real(RP) :: TOPO_ZsfcXY (ia,ja)
434  real(RP) :: TOPO_ZsfcUY (ia,ja)
435  real(RP) :: TOPO_ZsfcXV (ia,ja)
436  real(RP) :: TOPO_ZsfcUV (ia,ja)
437 
438  real(RP) :: GTRANS_QLIM (2*ka,2*ia,2*ja,3)
439  real(RP) :: QDZ(2*ka)
440  real(RP) :: QDX(2*ia)
441  real(RP) :: QDY(2*ja)
442  real(RP) :: AQAF (3)
443  real(RP) :: XSLOPE, YSLOPE
444  real(RP) :: Ztop
445  real(RP) :: DX_piece, DY_piece
446  real(RP) :: DX, DY, DZ
447 
448  integer :: I_QLIMtoLIM(3,7)
449 
450  integer :: iii, jjj, n
451  integer :: k, i, j, kk, ii, jj
452  !---------------------------------------------------------------------------
453 
454  ! calc absolute height at staggered position
455  ! at (x,y)
456  do j = 1, ja
457  do i = 1, ia
458  topo_zsfcxy(i,j) = topo_zsfc(i,j)
459  enddo
460  enddo
461  ! at (u,y)
462  do j = 1, ja
463  do i = 1, ia-1
464  topo_zsfcuy(i,j) = 0.5_rp * ( topo_zsfc(i,j) + topo_zsfc(i+1,j) )
465  enddo
466  enddo
467  ! at (x,v)
468  do j = 1, ja-1
469  do i = 1, ia
470  topo_zsfcxv(i,j) = 0.5_rp * ( topo_zsfc(i,j) + topo_zsfc(i,j+1) )
471  enddo
472  enddo
473  ! at (u,v)
474  do j = 1, ja-1
475  do i = 1, ia-1
476  topo_zsfcuv(i,j) = 0.25_rp * ( topo_zsfc(i ,j ) + topo_zsfc(i ,j+1) &
477  + topo_zsfc(i ,j+1) + topo_zsfc(i+1,j+1) )
478  enddo
479  enddo
480 
481  ! reset topography
482  topo_zsfc(:,:) = 0.d0
483 
484  call comm_vars8( topo_zsfcxy(:,:), 1 )
485  call comm_vars8( topo_zsfcuy(:,:), 2 )
486  call comm_vars8( topo_zsfcxv(:,:), 3 )
487  call comm_vars8( topo_zsfcuv(:,:), 4 )
488  call comm_wait ( topo_zsfcxy(:,:), 1 )
489  call comm_wait ( topo_zsfcuy(:,:), 2 )
490  call comm_wait ( topo_zsfcxv(:,:), 3 )
491  call comm_wait ( topo_zsfcuv(:,:), 4 )
492 
493  ! all height
494  do j = 1, ja
495  do i = 1, ia
496  ii = (i-1) * 2 + 1
497  jj = (j-1) * 2 + 1
498 
499  topo_zsfcall(ii ,jj ) = topo_zsfcxy(i,j)
500  topo_zsfcall(ii+1,jj ) = topo_zsfcuy(i,j)
501  topo_zsfcall(ii ,jj+1) = topo_zsfcxv(i,j)
502  topo_zsfcall(ii+1,jj+1) = topo_zsfcuv(i,j)
503  enddo
504  enddo
505 
506  ! length of control volume of quarter cell
507  do k = 1, ka
508  kk = (k-1) * 2 + 1
509 
510  qdz(kk ) = grid_cz(k) - grid_fz(k-1)
511  qdz(kk+1) = grid_fz(k) - grid_cz(k )
512  enddo
513 
514  do i = 1, ia-1
515  ii = (i-1) * 2 + 1
516 
517  qdx(ii ) = grid_fx(i ) - grid_cx(i)
518  qdx(ii+1) = grid_cx(i+1) - grid_fx(i)
519  enddo
520 
521  do j = 1, ja-1
522  jj = (j-1) * 2 + 1
523 
524  qdy(jj ) = grid_fy(j ) - grid_cy(j)
525  qdy(jj+1) = grid_cy(j+1) - grid_fy(j)
526  enddo
527 
528  ! quarter flux limiter
529  do jj = 1, 2*(ja-1)
530  do ii = 1, 2*(ia-1)
531  do kk = ks, 2*ke
532  dx_piece = qdx(ii) / real(gtrans_thinwall_xdiv,kind=rp)
533  dy_piece = qdy(jj) / real(gtrans_thinwall_ydiv,kind=rp)
534 
535  aqaf(1:3) = 0.0_rp
536  ztop = sum(qdz(ks:kk))
537 
538  !--- y-z face ---
539  yslope = ( topo_zsfcall(ii,jj+1) - topo_zsfcall(ii,jj) ) / qdy(jj)
540 
541  do jjj = 1, gtrans_thinwall_ydiv
542  dy = ( real(jjj,kind=RP) - 0.5_rp ) * dy_piece
543  dz = ztop - topo_zsfcall(ii,jj) - yslope * dy
544 
545  if ( dz > 0.0_rp ) then
546  if ( dz < qdz(kk) ) then
547  aqaf(i_fyz) = aqaf(i_fyz) + dz * dy_piece
548  else
549  aqaf(i_fyz) = aqaf(i_fyz) + qdz(kk) * dy_piece
550  endif
551  endif
552  enddo
553 
554  !--- x-z face ---
555  xslope = ( topo_zsfcall(ii+1,jj) - topo_zsfcall(ii,jj) ) / qdx(ii)
556 
557  do iii = 1, gtrans_thinwall_xdiv
558  dx = ( real(iii,kind=RP) - 0.5_RP ) * dx_piece
559  dz = ztop - topo_zsfcall(ii,jj) + xslope * dx
560 
561  if ( dz > 0.0_rp ) then
562  if ( dz < qdz(kk) ) then
563  aqaf(i_fxz) = aqaf(i_fxz) + dz * dx_piece
564  else
565  aqaf(i_fxz) = aqaf(i_fxz) + qdz(kk) * dx_piece
566  endif
567  endif
568  enddo
569 
570  !--- x-y face ---
571  do jjj = 1, gtrans_thinwall_ydiv
572  do iii = 1, gtrans_thinwall_xdiv
573  dx = ( real(iii,kind=RP) - 0.5_RP ) * dx_piece
574  dy = ( real(jjj,kind=RP) - 0.5_RP ) * dy_piece
575  dz = ztop - topo_zsfcall(ii,jj) - xslope * dx - yslope * dy
576 
577  if ( dz > 0.0_rp ) then
578  aqaf(i_fxy) = aqaf(i_fxy) + dx_piece * dy_piece
579  endif
580  enddo
581  enddo
582 
583  gtrans_qlim(kk,ii,jj,i_fyz) = aqaf(i_fyz) / ( qdy(jj) * qdz(kk) )
584  gtrans_qlim(kk,ii,jj,i_fxz) = aqaf(i_fxz) / ( qdx(ii) * qdz(kk) )
585  gtrans_qlim(kk,ii,jj,i_fxy) = aqaf(i_fxy) / ( qdy(jj) * qdx(ii) )
586  enddo
587  enddo
588  enddo
589 
590  ! index i,j,k
591  i_qlimtolim(1:3,i_xyz) = (/ 1, 1, 1 /)
592  i_qlimtolim(1:3,i_xyw) = (/ 1, 1, 0 /)
593  i_qlimtolim(1:3,i_uyw) = (/ 0, 1, 0 /)
594  i_qlimtolim(1:3,i_xvw) = (/ 1, 0, 0 /)
595  i_qlimtolim(1:3,i_uyz) = (/ 0, 1, 1 /)
596  i_qlimtolim(1:3,i_xvz) = (/ 1, 0, 1 /)
597  i_qlimtolim(1:3,i_uvz) = (/ 0, 0, 1 /)
598 
599  do n = 1, 7
600  do j = js, je
601  do i = is, ie
602  do k = ks, ke
603  ii = (i-1) * 2 + 1 - i_qlimtolim(1,n)
604  jj = (j-1) * 2 + 1 - i_qlimtolim(2,n)
605  kk = (k-1) * 2 + 1 - i_qlimtolim(3,n)
606 
607  gtrans_limyz(k,i,j,n) = 0.25_rp * ( gtrans_qlim(kk ,ii,jj,i_fyz) + gtrans_qlim(kk ,ii ,jj+1,i_fyz) &
608  + gtrans_qlim(kk+1,ii,jj,i_fyz) + gtrans_qlim(kk+1,ii ,jj+1,i_fyz) )
609  gtrans_limxz(k,i,j,n) = 0.25_rp * ( gtrans_qlim(kk ,ii,jj,i_fxz) + gtrans_qlim(kk ,ii+1,jj ,i_fxz) &
610  + gtrans_qlim(kk+1,ii,jj,i_fxz) + gtrans_qlim(kk+1,ii+1,jj ,i_fxz) )
611  gtrans_limxy(k,i,j,n) = 0.25_rp * ( gtrans_qlim(kk ,ii,jj,i_fxy) + gtrans_qlim(kk ,ii+1,jj+1,i_fxy) &
612  + gtrans_qlim(kk ,ii,jj,i_fxy) + gtrans_qlim(kk ,ii+1,jj ,i_fxy) )
613 
614  if ( gtrans_limyz(k,i,j,n) > 1.d0 ) then
615  write(*,*) 'xxx Facter miss! Check!'
616  write(*,*) k,i,j,n,gtrans_limyz(k,i,j,n)
617  call prc_mpistop
618  endif
619  enddo
620  enddo
621  enddo
622  enddo
623 
624  call comm_vars8( gtrans_limyz(:,:,:,i_xyz), 1 )
625  call comm_vars8( gtrans_limyz(:,:,:,i_xyw), 2 )
626  call comm_vars8( gtrans_limyz(:,:,:,i_uyw), 3 )
627  call comm_vars8( gtrans_limyz(:,:,:,i_xvw), 4 )
628  call comm_vars8( gtrans_limyz(:,:,:,i_uyz), 5 )
629  call comm_vars8( gtrans_limyz(:,:,:,i_xvz), 6 )
630  call comm_vars8( gtrans_limyz(:,:,:,i_uvz), 7 )
631  call comm_wait ( gtrans_limyz(:,:,:,i_xyz), 1 )
632  call comm_wait ( gtrans_limyz(:,:,:,i_xyw), 2 )
633  call comm_wait ( gtrans_limyz(:,:,:,i_uyw), 3 )
634  call comm_wait ( gtrans_limyz(:,:,:,i_xvw), 4 )
635  call comm_wait ( gtrans_limyz(:,:,:,i_uyz), 5 )
636  call comm_wait ( gtrans_limyz(:,:,:,i_xvz), 6 )
637  call comm_wait ( gtrans_limyz(:,:,:,i_uvz), 7 )
638 
639  call comm_vars8( gtrans_limxz(:,:,:,i_xyz), 8 )
640  call comm_vars8( gtrans_limxz(:,:,:,i_xyw), 9 )
641  call comm_vars8( gtrans_limxz(:,:,:,i_uyw), 10 )
642  call comm_vars8( gtrans_limxz(:,:,:,i_xvw), 11 )
643  call comm_vars8( gtrans_limxz(:,:,:,i_uyz), 12 )
644  call comm_vars8( gtrans_limxz(:,:,:,i_xvz), 13 )
645  call comm_vars8( gtrans_limxz(:,:,:,i_uvz), 14 )
646  call comm_wait ( gtrans_limxz(:,:,:,i_xyz), 8 )
647  call comm_wait ( gtrans_limxz(:,:,:,i_xyw), 9 )
648  call comm_wait ( gtrans_limxz(:,:,:,i_uyw), 10 )
649  call comm_wait ( gtrans_limxz(:,:,:,i_xvw), 11 )
650  call comm_wait ( gtrans_limxz(:,:,:,i_uyz), 12 )
651  call comm_wait ( gtrans_limxz(:,:,:,i_xvz), 13 )
652  call comm_wait ( gtrans_limxz(:,:,:,i_uvz), 14 )
653 
654  call comm_vars8( gtrans_limxy(:,:,:,i_xyz), 15 )
655  call comm_vars8( gtrans_limxy(:,:,:,i_xyw), 16 )
656  call comm_vars8( gtrans_limxy(:,:,:,i_uyw), 17 )
657  call comm_vars8( gtrans_limxy(:,:,:,i_xvw), 18 )
658  call comm_vars8( gtrans_limxy(:,:,:,i_uyz), 19 )
659  call comm_vars8( gtrans_limxy(:,:,:,i_xvz), 20 )
660  call comm_vars8( gtrans_limxy(:,:,:,i_uvz), 21 )
661  call comm_wait ( gtrans_limxy(:,:,:,i_xyz), 15 )
662  call comm_wait ( gtrans_limxy(:,:,:,i_xyw), 16 )
663  call comm_wait ( gtrans_limxy(:,:,:,i_uyw), 17 )
664  call comm_wait ( gtrans_limxy(:,:,:,i_xvw), 18 )
665  call comm_wait ( gtrans_limxy(:,:,:,i_uyz), 19 )
666  call comm_wait ( gtrans_limxy(:,:,:,i_xvz), 20 )
667  call comm_wait ( gtrans_limxy(:,:,:,i_uvz), 21 )
668 
669  return
670  end subroutine gtrans_thin_wall
671 
672  !-----------------------------------------------------------------------------
673  subroutine gtrans_step_mountain
674  use scale_comm, only: &
675  comm_vars8, &
676  comm_wait
677  implicit none
678 
679  integer :: i,j,k,n
680  !---------------------------------------------------------------------------
681 
682  do n = 1, 7
683  do j = js, je
684  do i = is, ie
685  do k = ks, ke
686 
687  if ( gtrans_limyz(k,i,j,n) > 0.0_rp ) then
688  gtrans_limyz(k,i,j,n) = 1.0_rp
689  endif
690 
691  if ( gtrans_limxz(k,i,j,n) > 0.0_rp ) then
692  gtrans_limxz(k,i,j,n) = 1.0_rp
693  endif
694 
695  if ( gtrans_limxy(k,i,j,n) > 0.0_rp ) then
696  gtrans_limxy(k,i,j,n) = 1.0_rp
697  endif
698  enddo
699  enddo
700  enddo
701  enddo
702 
703  call comm_vars8( gtrans_limyz(:,:,:,i_xyz), 1 )
704  call comm_vars8( gtrans_limyz(:,:,:,i_xyw), 2 )
705  call comm_vars8( gtrans_limyz(:,:,:,i_uyw), 3 )
706  call comm_vars8( gtrans_limyz(:,:,:,i_xvw), 4 )
707  call comm_vars8( gtrans_limyz(:,:,:,i_uyz), 5 )
708  call comm_vars8( gtrans_limyz(:,:,:,i_xvz), 6 )
709  call comm_vars8( gtrans_limyz(:,:,:,i_uvz), 7 )
710  call comm_wait ( gtrans_limyz(:,:,:,i_xyz), 1 )
711  call comm_wait ( gtrans_limyz(:,:,:,i_xyw), 2 )
712  call comm_wait ( gtrans_limyz(:,:,:,i_uyw), 3 )
713  call comm_wait ( gtrans_limyz(:,:,:,i_xvw), 4 )
714  call comm_wait ( gtrans_limyz(:,:,:,i_uyz), 5 )
715  call comm_wait ( gtrans_limyz(:,:,:,i_xvz), 6 )
716  call comm_wait ( gtrans_limyz(:,:,:,i_uvz), 7 )
717 
718  call comm_vars8( gtrans_limxz(:,:,:,i_xyz), 8 )
719  call comm_vars8( gtrans_limxz(:,:,:,i_xyw), 9 )
720  call comm_vars8( gtrans_limxz(:,:,:,i_uyw), 10 )
721  call comm_vars8( gtrans_limxz(:,:,:,i_xvw), 11 )
722  call comm_vars8( gtrans_limxz(:,:,:,i_uyz), 12 )
723  call comm_vars8( gtrans_limxz(:,:,:,i_xvz), 13 )
724  call comm_vars8( gtrans_limxz(:,:,:,i_uvz), 14 )
725  call comm_wait ( gtrans_limxz(:,:,:,i_xyz), 8 )
726  call comm_wait ( gtrans_limxz(:,:,:,i_xyw), 9 )
727  call comm_wait ( gtrans_limxz(:,:,:,i_uyw), 10 )
728  call comm_wait ( gtrans_limxz(:,:,:,i_xvw), 11 )
729  call comm_wait ( gtrans_limxz(:,:,:,i_uyz), 12 )
730  call comm_wait ( gtrans_limxz(:,:,:,i_xvz), 13 )
731  call comm_wait ( gtrans_limxz(:,:,:,i_uvz), 14 )
732 
733  call comm_vars8( gtrans_limxy(:,:,:,i_xyz), 15 )
734  call comm_vars8( gtrans_limxy(:,:,:,i_xyw), 16 )
735  call comm_vars8( gtrans_limxy(:,:,:,i_uyw), 17 )
736  call comm_vars8( gtrans_limxy(:,:,:,i_xvw), 18 )
737  call comm_vars8( gtrans_limxy(:,:,:,i_uyz), 19 )
738  call comm_vars8( gtrans_limxy(:,:,:,i_xvz), 20 )
739  call comm_vars8( gtrans_limxy(:,:,:,i_uvz), 21 )
740  call comm_wait ( gtrans_limxy(:,:,:,i_xyz), 15 )
741  call comm_wait ( gtrans_limxy(:,:,:,i_xyw), 16 )
742  call comm_wait ( gtrans_limxy(:,:,:,i_uyw), 17 )
743  call comm_wait ( gtrans_limxy(:,:,:,i_xvw), 18 )
744  call comm_wait ( gtrans_limxy(:,:,:,i_uyz), 19 )
745  call comm_wait ( gtrans_limxy(:,:,:,i_xvz), 20 )
746  call comm_wait ( gtrans_limxy(:,:,:,i_uvz), 21 )
747 
748  return
749  end subroutine gtrans_step_mountain
750 
751  !-----------------------------------------------------------------------------
753  subroutine gtrans_write
754  use scale_fileio, only: &
755  fileio_write
756  implicit none
757  !---------------------------------------------------------------------------
758 
759  if ( gtrans_out_basename /= '' ) then
760 
761  if( io_l ) write(io_fid_log,*)
762  if( io_l ) write(io_fid_log,*) '*** Output metrics file ***'
763 
764  call fileio_write( gtrans_mapf(:,:,1,i_xy), gtrans_out_basename, gtrans_out_title, & ! [IN]
765  'MAPF_X_XY', 'Map factor x-dir at XY', 'NIL', 'XY', gtrans_out_dtype ) ! [IN]
766  call fileio_write( gtrans_mapf(:,:,2,i_xy), gtrans_out_basename, gtrans_out_title, & ! [IN]
767  'MAPF_Y_XY', 'Map factor y-dir at XY', 'NIL', 'XY', gtrans_out_dtype ) ! [IN]
768  call fileio_write( gtrans_mapf(:,:,1,i_uy), gtrans_out_basename, gtrans_out_title, & ! [IN]
769  'MAPF_X_UY', 'Map factor x-dir at UY', 'NIL', 'UY', gtrans_out_dtype ) ! [IN]
770  call fileio_write( gtrans_mapf(:,:,2,i_uy), gtrans_out_basename, gtrans_out_title, & ! [IN]
771  'MAPF_Y_UY', 'Map factor y-dir at UY', 'NIL', 'UY', gtrans_out_dtype ) ! [IN]
772  call fileio_write( gtrans_mapf(:,:,1,i_xv), gtrans_out_basename, gtrans_out_title, & ! [IN]
773  'MAPF_X_XV', 'Map factor x-dir at XV', 'NIL', 'XY', gtrans_out_dtype ) ! [IN]
774  call fileio_write( gtrans_mapf(:,:,2,i_xv), gtrans_out_basename, gtrans_out_title, & ! [IN]
775  'MAPF_Y_XV', 'Map factor y-dir at XV', 'NIL', 'XY', gtrans_out_dtype ) ! [IN]
776  call fileio_write( gtrans_mapf(:,:,1,i_uv), gtrans_out_basename, gtrans_out_title, & ! [IN]
777  'MAPF_X_UV', 'Map factor x-dir at UV', 'NIL', 'UY', gtrans_out_dtype ) ! [IN]
778  call fileio_write( gtrans_mapf(:,:,2,i_uv), gtrans_out_basename, gtrans_out_title, & ! [IN]
779  'MAPF_Y_UV', 'Map factor y-dir at UV', 'NIL', 'UY', gtrans_out_dtype ) ! [IN]
780 
781  call fileio_write( gtrans_rotc(:,:,1), gtrans_out_basename, gtrans_out_title, & ! [IN]
782  'ROTC_COS', 'Rotation factor (cos)', 'NIL', 'XY', gtrans_out_dtype ) ! [IN]
783  call fileio_write( gtrans_rotc(:,:,2), gtrans_out_basename, gtrans_out_title, & ! [IN]
784  'ROTC_SIN', 'Rotation factor (sin)', 'NIL', 'XY', gtrans_out_dtype ) ! [IN]
785 
786  call fileio_write( gtrans_rotc(:,:,1), gtrans_out_basename, gtrans_out_title, & ! [IN]
787  'ROTC_COS', 'Rotation factor (cos)', 'NIL', 'XY', gtrans_out_dtype ) ! [IN]
788 
789  endif
790 
791  return
792  end subroutine gtrans_write
793 
794 end module scale_gridtrans
integer, public is
start point of inner domain: x, local
integer, public i_xvz
integer, public je
end point of inner domain: y, local
real(rp), dimension(:), allocatable, public grid_rcdy
reciprocal of center-dy
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:,:,:,:), allocatable, public gtrans_limxy
flux limiter x-y face
real(rp), dimension(:,:,:,:), allocatable, public gtrans_j23g
(2,3) element of Jacobian matrix * {G}^1/2
logical, public io_l
output log or not? (this process)
Definition: scale_stdio.F90:59
real(rp), dimension(:), allocatable, public grid_rcdx
reciprocal of center-dx
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
module STDIO
Definition: scale_stdio.F90:12
integer, public ke
end point of inner domain: z, local
real(rp), public gtrans_j33g
(3,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:,:,:), allocatable, public real_fz
geopotential height [m] (cell face )
real(rp), dimension(:,:,:), allocatable, public real_cz
geopotential height [m] (cell center)
integer, public i_xy
integer, public i_fyz
real(rp), dimension(:), allocatable, public grid_rfdy
reciprocal of face-dy
real(rp), dimension(:), allocatable, public grid_rcdz
reciprocal of center-dz
real(rp), dimension(:,:,:,:), allocatable, public gtrans_limyz
flux limiter y-z face
real(rp), dimension(:,:,:,:), allocatable, public gtrans_limxz
flux limiter x-z face
module FILE I/O (netcdf)
integer, public i_xvw
real(rp), dimension(:), allocatable, public grid_fx
face coordinate [m]: x, local
module grid index
real(rp), dimension(:,:,:,:), allocatable, public gtrans_mapf
map factor
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:,:), allocatable, public real_latx
latitude at staggered point (uy) [rad,-pi,pi]
module GRIDTRANS
module GRID (real space)
integer, public ka
of z whole cells (local, with HALO)
integer, public i_uy
integer, public i_xyw
real(rp), dimension(:), allocatable, public grid_fz
face coordinate [m]: z, local=global
module COMMUNICATION
Definition: scale_comm.F90:23
integer, public js
start point of inner domain: y, local
module PROCESS
integer, public i_uyw
real(rp), dimension(:,:,:), allocatable, public gtrans_rotc
rotation coefficient
integer, public i_fxy
real(rp), dimension(:,:,:,:), allocatable, public gtrans_j13g
(1,3) element of Jacobian matrix * {G}^1/2
subroutine, public gtrans_setup
Setup.
subroutine, public mprj_mapfactor(lat, m1, m2)
(x,y) -> (lon,lat)
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:), allocatable, public grid_cx
center coordinate [m]: x, local
module GRID (cartesian)
integer, public i_xv
module profiler
Definition: scale_prof.F90:10
integer, public ie
end point of inner domain: x, local
real(rp), dimension(:,:,:,:), allocatable, public gtrans_gsqrt
transformation metrics from Z to Xi, {G}^1/2
subroutine, public real_calc_areavol(MAPF)
Calc control area/volume.
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 i_uyz
module PRECISION
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
real(rp), dimension(:), allocatable, public grid_rfdx
reciprocal of face-dx
module TOPOGRAPHY
real(rp), dimension(:), allocatable, public grid_rfdz
reciprocal of face-dz
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
integer, public i_uvz
integer, public i_xyz
real(rp), dimension(:,:), allocatable, public real_lat
latitude [rad,-pi,pi]
integer, public io_fid_log
Log file ID.
Definition: scale_stdio.F90:56
integer, public i_fxz
integer, public i_uv
integer, parameter, public rp
real(rp), dimension(:,:), allocatable, public real_latxy
latitude at staggered point (uv) [rad,-pi,pi]
real(rp), dimension(:), allocatable, public grid_cy
center coordinate [m]: y, local
subroutine gtrans_rotcoef
Calculate rotation coeffient.
module Map projection
real(rp), dimension(:,:), allocatable, public real_laty
latitude at staggered point (xv) [rad,-pi,pi]
integer, public ja
of y whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local