SCALE-RM
scale_atmos_grid_cartesC_metric.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 !-------------------------------------------------------------------------------
11 #include "scalelib.h"
13  !-----------------------------------------------------------------------------
14  !
15  !++ used modules
16  !
17  use scale_precision
18  use scale_io
19  use scale_prof
21  !-----------------------------------------------------------------------------
22  implicit none
23  private
24  !-----------------------------------------------------------------------------
25  !
26  !++ Public procedure
27  !
29 
30  !-----------------------------------------------------------------------------
31  !
32  !++ Public parameters & variables
33  !
34  real(RP), public, allocatable :: atmos_grid_cartesc_metric_mapf (:,:,:,:)
35  real(RP), public, allocatable :: atmos_grid_cartesc_metric_rotc (:,:,:)
36 
37  real(RP), public, allocatable :: atmos_grid_cartesc_metric_gsqrt(:,:,:,:)
38  real(RP), public, allocatable :: atmos_grid_cartesc_metric_j13g (:,:,:,:)
39  real(RP), public, allocatable :: atmos_grid_cartesc_metric_j23g (:,:,:,:)
40  real(RP), public :: atmos_grid_cartesc_metric_j33g
41 
42  real(RP), public, allocatable :: atmos_grid_cartesc_metric_limyz(:,:,:,:)
43  real(RP), public, allocatable :: atmos_grid_cartesc_metric_limxz(:,:,:,:)
44  real(RP), public, allocatable :: atmos_grid_cartesc_metric_limxy(:,:,:,:)
45 
46  !-----------------------------------------------------------------------------
47  !
48  !++ Private procedure
49  !
50  private :: atmos_grid_cartesc_metric_mapfactor
51  private :: atmos_grid_cartesc_metric_terrainfollowing
52  private :: atmos_grid_cartesc_metric_thin_wall
53  private :: atmos_grid_cartesc_metric_step_mountain
54  private :: atmos_grid_cartesc_metric_write
55 
56  !-----------------------------------------------------------------------------
57  !
58  !++ Private parameters & variables
59  !
60  character(len=H_LONG), private :: atmos_grid_cartesc_metric_out_basename = ''
61  character(len=H_MID), private :: atmos_grid_cartesc_metric_out_title = 'SCALE-RM GEOMETRICS'
62  character(len=H_SHORT), private :: atmos_grid_cartesc_metric_out_dtype = 'DEFAULT'
63 
64  character(len=H_SHORT), private :: atmos_grid_cartesc_metric_topo_type = 'TERRAINFOLLOWING'
65  integer, private :: atmos_grid_cartesc_metric_thinwall_xdiv = 50
66  integer, private :: atmos_grid_cartesc_metric_thinwall_ydiv = 50
67  !-----------------------------------------------------------------------------
68 contains
69  !-----------------------------------------------------------------------------
72  use scale_prc, only: &
73  prc_abort
74  implicit none
75 
76  namelist / param_atmos_grid_cartesc_metric / &
77  atmos_grid_cartesc_metric_out_basename, &
78  atmos_grid_cartesc_metric_out_dtype, &
79  atmos_grid_cartesc_metric_topo_type, &
80  atmos_grid_cartesc_metric_thinwall_xdiv, &
81  atmos_grid_cartesc_metric_thinwall_ydiv
82 
83  integer :: ierr
84  !---------------------------------------------------------------------------
85 
86  log_newline
87  log_info("ATMOS_GRID_CARTESC_METRIC_setup",*) 'Setup'
88 
89  !--- read namelist
90  rewind(io_fid_conf)
91  read(io_fid_conf,nml=param_atmos_grid_cartesc_metric,iostat=ierr)
92  if( ierr < 0 ) then !--- missing
93  log_info("ATMOS_GRID_CARTESC_METRIC_setup",*) 'Not found namelist. Default used.'
94  elseif( ierr > 0 ) then !--- fatal error
95  log_error("ATMOS_GRID_CARTESC_METRIC_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_GRID_CARTESC_METRIC. Check!'
96  call prc_abort
97  endif
98  log_nml(param_atmos_grid_cartesc_metric)
99 
100  allocate( atmos_grid_cartesc_metric_mapf(ia,ja,2,4) )
101 
102  allocate( atmos_grid_cartesc_metric_rotc(ia,ja,2) )
103 
104  allocate( atmos_grid_cartesc_metric_gsqrt(ka,ia,ja,7) )
105  allocate( atmos_grid_cartesc_metric_j13g(ka,ia,ja,7) )
106  allocate( atmos_grid_cartesc_metric_j23g(ka,ia,ja,7) )
107 
108  atmos_grid_cartesc_metric_gsqrt(:,:,:,:) = 1.0_rp
109  atmos_grid_cartesc_metric_j13g(:,:,:,:) = 0.0_rp
110  atmos_grid_cartesc_metric_j23g(:,:,:,:) = 0.0_rp
112 
113  allocate( atmos_grid_cartesc_metric_limyz(ka,ia,ja,7) )
114  allocate( atmos_grid_cartesc_metric_limxz(ka,ia,ja,7) )
115  allocate( atmos_grid_cartesc_metric_limxy(ka,ia,ja,7) )
116  atmos_grid_cartesc_metric_limyz(:,:,:,:) = 1.0_rp
117  atmos_grid_cartesc_metric_limxz(:,:,:,:) = 1.0_rp
118  atmos_grid_cartesc_metric_limxy(:,:,:,:) = 1.0_rp
119 
120  ! calc metrics for orthogonal curvelinear coordinate
121  call atmos_grid_cartesc_metric_mapfactor
122 
123  ! calc coeficient for rotaion of velocity vector
125 
126  ! calc metrics for terrain-following,step-mountain,thin-wall coordinate
127  log_newline
128  log_info("ATMOS_GRID_CARTESC_METRIC_setup",*) 'Terrain coordinate type : ', trim(atmos_grid_cartesc_metric_topo_type)
129  select case(atmos_grid_cartesc_metric_topo_type)
130  case('TERRAINFOLLOWING')
131  log_info_cont(*) '=> Terrain-following method'
132  call atmos_grid_cartesc_metric_terrainfollowing
133  case('STEPMOUNTAIN')
134  log_info_cont(*) '=> Step-mountain method'
135  call atmos_grid_cartesc_metric_thin_wall
136  call atmos_grid_cartesc_metric_step_mountain
137  case('THINWALL')
138  log_info_cont(*) '=> Thin-wall approximation method'
139  call atmos_grid_cartesc_metric_thin_wall
140  case default
141  log_error("ATMOS_GRID_CARTESC_METRIC_setup",*) 'Unsupported ATMOS_GRID_CARTESC_METRIC_TOPO_type. STOP'
142  call prc_abort
143  end select
144 
145  ! output metrics (for debug)
146  call atmos_grid_cartesc_metric_write
147 
148  return
149  end subroutine atmos_grid_cartesc_metric_setup
150 
151  !-----------------------------------------------------------------------------
153  subroutine atmos_grid_cartesc_metric_mapfactor
154  use scale_mapprojection, only: &
156  use scale_atmos_grid_cartesc_real, only: &
162  implicit none
163  !---------------------------------------------------------------------------
164 
165  call mapprojection_mapfactor( ia, 1, ia, ja, 1, ja, &
167  call mapprojection_mapfactor( ia, 1, ia, ja, 1, ja, &
169  call mapprojection_mapfactor( ia, 1, ia, ja, 1, ja, &
171  call mapprojection_mapfactor( ia, 1, ia, ja, 1, ja, &
173 
175 
176  return
177  end subroutine atmos_grid_cartesc_metric_mapfactor
178 
179  !-----------------------------------------------------------------------------
182  use scale_mapprojection, only: &
184  use scale_atmos_grid_cartesc_real, only: &
187  implicit none
188  !---------------------------------------------------------------------------
189 
190  call mapprojection_rotcoef( ia, 1, ia, ja, 1, ja, &
191  atmos_grid_cartesc_real_lon(:,:), & ! [IN]
192  atmos_grid_cartesc_real_lat(:,:), & ! [IN]
193  atmos_grid_cartesc_metric_rotc(:,:,:) ) ! [OUT]
194 
195  return
196  end subroutine atmos_grid_cartesc_metric_rotcoef
197 
198  !-----------------------------------------------------------------------------
200  subroutine atmos_grid_cartesc_metric_terrainfollowing
201  use scale_atmos_grid_cartesc, only: &
208  use scale_atmos_grid_cartesc_real, only: &
217  use scale_comm_cartesc, only: &
218  comm_vars8, &
219  comm_wait
220  implicit none
221 
222  integer :: k, i, j
223  !---------------------------------------------------------------------------
224 
225  ! G^1/2
226  do j = js, je
227  do i = is, ie
228  ! at (x,y,z)
229  do k = 1, ka
231  enddo
232 
233  ! at (x,y,w)
234  do k = 1, ka-1
236  enddo
238 
239  ! at (u,y,w)
240  do k = 1, ka-1
242  enddo
244 
245  ! at (x,v,w)
246  do k = 1, ka-1
248  enddo
250 
251  ! at (u,y,z)
252  do k = 1, ka
254  enddo
255 
256  ! at (x,v,z)
257  do k = 1, ka
259  enddo
260 
261  ! at (u,v,z)
262  do k = 1, ka
264  enddo
265  enddo
266  enddo
267 
268  call comm_vars8( atmos_grid_cartesc_metric_gsqrt(:,:,:,1), 1 )
269  call comm_vars8( atmos_grid_cartesc_metric_gsqrt(:,:,:,2), 2 )
270  call comm_vars8( atmos_grid_cartesc_metric_gsqrt(:,:,:,3), 3 )
271  call comm_vars8( atmos_grid_cartesc_metric_gsqrt(:,:,:,4), 4 )
272  call comm_vars8( atmos_grid_cartesc_metric_gsqrt(:,:,:,5), 5 )
273  call comm_vars8( atmos_grid_cartesc_metric_gsqrt(:,:,:,6), 6 )
274  call comm_vars8( atmos_grid_cartesc_metric_gsqrt(:,:,:,7), 7 )
275  call comm_wait ( atmos_grid_cartesc_metric_gsqrt(:,:,:,1), 1 )
276  call comm_wait ( atmos_grid_cartesc_metric_gsqrt(:,:,:,2), 2 )
277  call comm_wait ( atmos_grid_cartesc_metric_gsqrt(:,:,:,3), 3 )
278  call comm_wait ( atmos_grid_cartesc_metric_gsqrt(:,:,:,4), 4 )
279  call comm_wait ( atmos_grid_cartesc_metric_gsqrt(:,:,:,5), 5 )
280  call comm_wait ( atmos_grid_cartesc_metric_gsqrt(:,:,:,6), 6 )
281  call comm_wait ( atmos_grid_cartesc_metric_gsqrt(:,:,:,7), 7 )
282 
283  ! Jacobian * G^1/2
284  do j = js, je
285  do i = is, ie
286  do k = 1, ka
294  enddo
295  enddo
296  enddo
297 
298  do j = js, je
299  do i = is, ie
300  do k = 1, ka
308  enddo
309  enddo
310  enddo
311 
312  atmos_grid_cartesc_metric_j33g = 1.0_rp ! - 1 / G^1/2 * G^1/2
313 
314  call comm_vars8( atmos_grid_cartesc_metric_j13g(:,:,:,i_xyz), 8 )
315  call comm_vars8( atmos_grid_cartesc_metric_j13g(:,:,:,i_xyw), 9 )
316  call comm_vars8( atmos_grid_cartesc_metric_j13g(:,:,:,i_uyw), 10 )
317  call comm_vars8( atmos_grid_cartesc_metric_j13g(:,:,:,i_xvw), 11 )
318  call comm_vars8( atmos_grid_cartesc_metric_j13g(:,:,:,i_uyz), 12 )
319  call comm_vars8( atmos_grid_cartesc_metric_j13g(:,:,:,i_xvz), 13 )
320  call comm_vars8( atmos_grid_cartesc_metric_j13g(:,:,:,i_uvz), 14 )
321  call comm_wait ( atmos_grid_cartesc_metric_j13g(:,:,:,i_xyz), 8 )
322  call comm_wait ( atmos_grid_cartesc_metric_j13g(:,:,:,i_xyw), 9 )
323  call comm_wait ( atmos_grid_cartesc_metric_j13g(:,:,:,i_uyw), 10 )
324  call comm_wait ( atmos_grid_cartesc_metric_j13g(:,:,:,i_xvw), 11 )
325  call comm_wait ( atmos_grid_cartesc_metric_j13g(:,:,:,i_uyz), 12 )
326  call comm_wait ( atmos_grid_cartesc_metric_j13g(:,:,:,i_xvz), 13 )
327  call comm_wait ( atmos_grid_cartesc_metric_j13g(:,:,:,i_uvz), 14 )
328 
329  call comm_vars8( atmos_grid_cartesc_metric_j23g(:,:,:,i_xyz), 15 )
330  call comm_vars8( atmos_grid_cartesc_metric_j23g(:,:,:,i_xyw), 16 )
331  call comm_vars8( atmos_grid_cartesc_metric_j23g(:,:,:,i_uyw), 17 )
332  call comm_vars8( atmos_grid_cartesc_metric_j23g(:,:,:,i_xvw), 18 )
333  call comm_vars8( atmos_grid_cartesc_metric_j23g(:,:,:,i_uyz), 19 )
334  call comm_vars8( atmos_grid_cartesc_metric_j23g(:,:,:,i_xvz), 20 )
335  call comm_vars8( atmos_grid_cartesc_metric_j23g(:,:,:,i_uvz), 21 )
336  call comm_wait ( atmos_grid_cartesc_metric_j23g(:,:,:,i_xyz), 15 )
337  call comm_wait ( atmos_grid_cartesc_metric_j23g(:,:,:,i_xyw), 16 )
338  call comm_wait ( atmos_grid_cartesc_metric_j23g(:,:,:,i_uyw), 17 )
339  call comm_wait ( atmos_grid_cartesc_metric_j23g(:,:,:,i_xvw), 18 )
340  call comm_wait ( atmos_grid_cartesc_metric_j23g(:,:,:,i_uyz), 19 )
341  call comm_wait ( atmos_grid_cartesc_metric_j23g(:,:,:,i_xvz), 20 )
342  call comm_wait ( atmos_grid_cartesc_metric_j23g(:,:,:,i_uvz), 21 )
343 
344  return
345  end subroutine atmos_grid_cartesc_metric_terrainfollowing
346 
347  !-----------------------------------------------------------------------------
348  subroutine atmos_grid_cartesc_metric_thin_wall
349  use scale_prc, only: &
350  prc_abort
351  use scale_atmos_grid_cartesc, only: &
358  use scale_topography, only : &
359  topo_zsfc
360  use scale_comm_cartesc, only : &
361  comm_vars8, &
362  comm_wait
363  implicit none
364 
365  real(RP) :: TOPO_ZsfcALL(2*ia,2*ja)
366  real(RP) :: TOPO_ZsfcXY (ia,ja)
367  real(RP) :: TOPO_ZsfcUY (ia,ja)
368  real(RP) :: TOPO_ZsfcXV (ia,ja)
369  real(RP) :: TOPO_ZsfcUV (ia,ja)
370 
371  real(RP) :: ATMOS_GRID_CARTESC_METRIC_QLIM (2*ka,2*ia,2*ja,3)
372  real(RP) :: QDZ(2*ka)
373  real(RP) :: QDX(2*ia)
374  real(RP) :: QDY(2*ja)
375  real(RP) :: AQAF (3)
376  real(RP) :: XSLOPE, YSLOPE
377  real(RP) :: Ztop
378  real(RP) :: DX_piece, DY_piece
379  real(RP) :: DX, DY, DZ
380 
381  integer :: I_QLIMtoLIM(3,7)
382 
383  integer :: iii, jjj, n
384  integer :: k, i, j, kk, ii, jj
385  !---------------------------------------------------------------------------
386 
387  ! calc absolute height at staggered position
388  ! at (x,y)
389  do j = 1, ja
390  do i = 1, ia
391  topo_zsfcxy(i,j) = topo_zsfc(i,j)
392  enddo
393  enddo
394  ! at (u,y)
395  do j = 1, ja
396  do i = 1, ia-1
397  topo_zsfcuy(i,j) = 0.5_rp * ( topo_zsfc(i,j) + topo_zsfc(i+1,j) )
398  enddo
399  enddo
400  ! at (x,v)
401  do j = 1, ja-1
402  do i = 1, ia
403  topo_zsfcxv(i,j) = 0.5_rp * ( topo_zsfc(i,j) + topo_zsfc(i,j+1) )
404  enddo
405  enddo
406  ! at (u,v)
407  do j = 1, ja-1
408  do i = 1, ia-1
409  topo_zsfcuv(i,j) = 0.25_rp * ( topo_zsfc(i ,j ) + topo_zsfc(i ,j+1) &
410  + topo_zsfc(i ,j+1) + topo_zsfc(i+1,j+1) )
411  enddo
412  enddo
413 
414  ! reset topography
415  topo_zsfc(:,:) = 0.d0
416 
417  call comm_vars8( topo_zsfcxy(:,:), 1 )
418  call comm_vars8( topo_zsfcuy(:,:), 2 )
419  call comm_vars8( topo_zsfcxv(:,:), 3 )
420  call comm_vars8( topo_zsfcuv(:,:), 4 )
421  call comm_wait ( topo_zsfcxy(:,:), 1 )
422  call comm_wait ( topo_zsfcuy(:,:), 2 )
423  call comm_wait ( topo_zsfcxv(:,:), 3 )
424  call comm_wait ( topo_zsfcuv(:,:), 4 )
425 
426  ! all height
427  do j = 1, ja
428  do i = 1, ia
429  ii = (i-1) * 2 + 1
430  jj = (j-1) * 2 + 1
431 
432  topo_zsfcall(ii ,jj ) = topo_zsfcxy(i,j)
433  topo_zsfcall(ii+1,jj ) = topo_zsfcuy(i,j)
434  topo_zsfcall(ii ,jj+1) = topo_zsfcxv(i,j)
435  topo_zsfcall(ii+1,jj+1) = topo_zsfcuv(i,j)
436  enddo
437  enddo
438 
439  ! length of control volume of quarter cell
440  do k = 1, ka
441  kk = (k-1) * 2 + 1
442 
443  qdz(kk ) = atmos_grid_cartesc_cz(k) - atmos_grid_cartesc_fz(k-1)
444  qdz(kk+1) = atmos_grid_cartesc_fz(k) - atmos_grid_cartesc_cz(k )
445  enddo
446 
447  do i = 1, ia-1
448  ii = (i-1) * 2 + 1
449 
450  qdx(ii ) = atmos_grid_cartesc_fx(i ) - atmos_grid_cartesc_cx(i)
451  qdx(ii+1) = atmos_grid_cartesc_cx(i+1) - atmos_grid_cartesc_fx(i)
452  enddo
453 
454  do j = 1, ja-1
455  jj = (j-1) * 2 + 1
456 
457  qdy(jj ) = atmos_grid_cartesc_fy(j ) - atmos_grid_cartesc_cy(j)
458  qdy(jj+1) = atmos_grid_cartesc_cy(j+1) - atmos_grid_cartesc_fy(j)
459  enddo
460 
461  ! quarter flux limiter
462  do jj = 1, 2*(ja-1)
463  do ii = 1, 2*(ia-1)
464  do kk = ks, 2*ke
465  dx_piece = qdx(ii) / real(atmos_grid_cartesc_metric_thinwall_xdiv,kind=rp)
466  dy_piece = qdy(jj) / real(atmos_grid_cartesc_metric_thinwall_ydiv,kind=rp)
467 
468  aqaf(1:3) = 0.0_rp
469  ztop = sum(qdz(ks:kk))
470 
471  !--- y-z face ---
472  yslope = ( topo_zsfcall(ii,jj+1) - topo_zsfcall(ii,jj) ) / qdy(jj)
473 
474  do jjj = 1, atmos_grid_cartesc_metric_thinwall_ydiv
475  dy = ( real(jjj,kind=RP) - 0.5_rp ) * dy_piece
476  dz = ztop - topo_zsfcall(ii,jj) - yslope * dy
477 
478  if ( dz > 0.0_rp ) then
479  if ( dz < qdz(kk) ) then
480  aqaf(i_fyz) = aqaf(i_fyz) + dz * dy_piece
481  else
482  aqaf(i_fyz) = aqaf(i_fyz) + qdz(kk) * dy_piece
483  endif
484  endif
485  enddo
486 
487  !--- x-z face ---
488  xslope = ( topo_zsfcall(ii+1,jj) - topo_zsfcall(ii,jj) ) / qdx(ii)
489 
490  do iii = 1, atmos_grid_cartesc_metric_thinwall_xdiv
491  dx = ( real(iii,kind=RP) - 0.5_RP ) * dx_piece
492  dz = ztop - topo_zsfcall(ii,jj) + xslope * dx
493 
494  if ( dz > 0.0_rp ) then
495  if ( dz < qdz(kk) ) then
496  aqaf(i_fxz) = aqaf(i_fxz) + dz * dx_piece
497  else
498  aqaf(i_fxz) = aqaf(i_fxz) + qdz(kk) * dx_piece
499  endif
500  endif
501  enddo
502 
503  !--- x-y face ---
504  do jjj = 1, atmos_grid_cartesc_metric_thinwall_ydiv
505  do iii = 1, atmos_grid_cartesc_metric_thinwall_xdiv
506  dx = ( real(iii,kind=RP) - 0.5_RP ) * dx_piece
507  dy = ( real(jjj,kind=RP) - 0.5_RP ) * dy_piece
508  dz = ztop - topo_zsfcall(ii,jj) - xslope * dx - yslope * dy
509 
510  if ( dz > 0.0_rp ) then
511  aqaf(i_fxy) = aqaf(i_fxy) + dx_piece * dy_piece
512  endif
513  enddo
514  enddo
515 
516  atmos_grid_cartesc_metric_qlim(kk,ii,jj,i_fyz) = aqaf(i_fyz) / ( qdy(jj) * qdz(kk) )
517  atmos_grid_cartesc_metric_qlim(kk,ii,jj,i_fxz) = aqaf(i_fxz) / ( qdx(ii) * qdz(kk) )
518  atmos_grid_cartesc_metric_qlim(kk,ii,jj,i_fxy) = aqaf(i_fxy) / ( qdy(jj) * qdx(ii) )
519  enddo
520  enddo
521  enddo
522 
523  ! index i,j,k
524  i_qlimtolim(1:3,i_xyz) = (/ 1, 1, 1 /)
525  i_qlimtolim(1:3,i_xyw) = (/ 1, 1, 0 /)
526  i_qlimtolim(1:3,i_uyw) = (/ 0, 1, 0 /)
527  i_qlimtolim(1:3,i_xvw) = (/ 1, 0, 0 /)
528  i_qlimtolim(1:3,i_uyz) = (/ 0, 1, 1 /)
529  i_qlimtolim(1:3,i_xvz) = (/ 1, 0, 1 /)
530  i_qlimtolim(1:3,i_uvz) = (/ 0, 0, 1 /)
531 
532  do n = 1, 7
533  do j = js, je
534  do i = is, ie
535  do k = ks, ke
536  ii = (i-1) * 2 + 1 - i_qlimtolim(1,n)
537  jj = (j-1) * 2 + 1 - i_qlimtolim(2,n)
538  kk = (k-1) * 2 + 1 - i_qlimtolim(3,n)
539 
540  atmos_grid_cartesc_metric_limyz(k,i,j,n) = 0.25_rp * ( atmos_grid_cartesc_metric_qlim(kk ,ii,jj,i_fyz) + atmos_grid_cartesc_metric_qlim(kk ,ii ,jj+1,i_fyz) &
541  + atmos_grid_cartesc_metric_qlim(kk+1,ii,jj,i_fyz) + atmos_grid_cartesc_metric_qlim(kk+1,ii ,jj+1,i_fyz) )
542  atmos_grid_cartesc_metric_limxz(k,i,j,n) = 0.25_rp * ( atmos_grid_cartesc_metric_qlim(kk ,ii,jj,i_fxz) + atmos_grid_cartesc_metric_qlim(kk ,ii+1,jj ,i_fxz) &
543  + atmos_grid_cartesc_metric_qlim(kk+1,ii,jj,i_fxz) + atmos_grid_cartesc_metric_qlim(kk+1,ii+1,jj ,i_fxz) )
544  atmos_grid_cartesc_metric_limxy(k,i,j,n) = 0.25_rp * ( atmos_grid_cartesc_metric_qlim(kk ,ii,jj,i_fxy) + atmos_grid_cartesc_metric_qlim(kk ,ii+1,jj+1,i_fxy) &
545  + atmos_grid_cartesc_metric_qlim(kk ,ii,jj,i_fxy) + atmos_grid_cartesc_metric_qlim(kk ,ii+1,jj ,i_fxy) )
546 
547  if ( atmos_grid_cartesc_metric_limyz(k,i,j,n) > 1.d0 ) then
548  log_error("ATMOS_GRID_CARTESC_METRIC_thin_wall",*) 'Facter miss! Check!'
549  log_error_cont(*) k,i,j,n,atmos_grid_cartesc_metric_limyz(k,i,j,n)
550  call prc_abort
551  endif
552  enddo
553  enddo
554  enddo
555  enddo
556 
557  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyz), 1 )
558  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyw), 2 )
559  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyw), 3 )
560  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvw), 4 )
561  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyz), 5 )
562  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvz), 6 )
563  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_uvz), 7 )
564  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyz), 1 )
565  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyw), 2 )
566  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyw), 3 )
567  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvw), 4 )
568  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyz), 5 )
569  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvz), 6 )
570  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_uvz), 7 )
571 
572  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyz), 8 )
573  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyw), 9 )
574  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyw), 10 )
575  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvw), 11 )
576  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyz), 12 )
577  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvz), 13 )
578  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_uvz), 14 )
579  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyz), 8 )
580  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyw), 9 )
581  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyw), 10 )
582  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvw), 11 )
583  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyz), 12 )
584  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvz), 13 )
585  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_uvz), 14 )
586 
587  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyz), 15 )
588  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyw), 16 )
589  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyw), 17 )
590  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvw), 18 )
591  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyz), 19 )
592  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvz), 20 )
593  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_uvz), 21 )
594  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyz), 15 )
595  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyw), 16 )
596  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyw), 17 )
597  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvw), 18 )
598  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyz), 19 )
599  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvz), 20 )
600  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_uvz), 21 )
601 
602  return
603  end subroutine atmos_grid_cartesc_metric_thin_wall
604 
605  !-----------------------------------------------------------------------------
606  subroutine atmos_grid_cartesc_metric_step_mountain
607  use scale_comm_cartesc, only: &
608  comm_vars8, &
609  comm_wait
610  implicit none
611 
612  integer :: i,j,k,n
613  !---------------------------------------------------------------------------
614 
615  do n = 1, 7
616  do j = js, je
617  do i = is, ie
618  do k = ks, ke
619 
620  if ( atmos_grid_cartesc_metric_limyz(k,i,j,n) > 0.0_rp ) then
621  atmos_grid_cartesc_metric_limyz(k,i,j,n) = 1.0_rp
622  endif
623 
624  if ( atmos_grid_cartesc_metric_limxz(k,i,j,n) > 0.0_rp ) then
625  atmos_grid_cartesc_metric_limxz(k,i,j,n) = 1.0_rp
626  endif
627 
628  if ( atmos_grid_cartesc_metric_limxy(k,i,j,n) > 0.0_rp ) then
629  atmos_grid_cartesc_metric_limxy(k,i,j,n) = 1.0_rp
630  endif
631  enddo
632  enddo
633  enddo
634  enddo
635 
636  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyz), 1 )
637  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyw), 2 )
638  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyw), 3 )
639  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvw), 4 )
640  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyz), 5 )
641  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvz), 6 )
642  call comm_vars8( atmos_grid_cartesc_metric_limyz(:,:,:,i_uvz), 7 )
643  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyz), 1 )
644  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xyw), 2 )
645  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyw), 3 )
646  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvw), 4 )
647  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_uyz), 5 )
648  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_xvz), 6 )
649  call comm_wait ( atmos_grid_cartesc_metric_limyz(:,:,:,i_uvz), 7 )
650 
651  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyz), 8 )
652  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyw), 9 )
653  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyw), 10 )
654  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvw), 11 )
655  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyz), 12 )
656  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvz), 13 )
657  call comm_vars8( atmos_grid_cartesc_metric_limxz(:,:,:,i_uvz), 14 )
658  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyz), 8 )
659  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xyw), 9 )
660  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyw), 10 )
661  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvw), 11 )
662  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_uyz), 12 )
663  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_xvz), 13 )
664  call comm_wait ( atmos_grid_cartesc_metric_limxz(:,:,:,i_uvz), 14 )
665 
666  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyz), 15 )
667  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyw), 16 )
668  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyw), 17 )
669  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvw), 18 )
670  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyz), 19 )
671  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvz), 20 )
672  call comm_vars8( atmos_grid_cartesc_metric_limxy(:,:,:,i_uvz), 21 )
673  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyz), 15 )
674  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xyw), 16 )
675  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyw), 17 )
676  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvw), 18 )
677  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_uyz), 19 )
678  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_xvz), 20 )
679  call comm_wait ( atmos_grid_cartesc_metric_limxy(:,:,:,i_uvz), 21 )
680 
681  return
682  end subroutine atmos_grid_cartesc_metric_step_mountain
683 
684  !-----------------------------------------------------------------------------
686  subroutine atmos_grid_cartesc_metric_write
687  use scale_const, only: &
689  use scale_vector, only: &
691  use scale_file_cartesc, only: &
692  file_cartesc_write
693  use scale_atmos_grid_cartesc_real, only: &
698  use scale_mapprojection, only: &
699  mapprojection_lonlat2xy
700  implicit none
701 
702  real(RP) :: check_X_XY(ia,ja)
703  real(RP) :: check_Y_XY(ia,ja)
704  real(RP) :: distance (ia,ja)
705 
706  integer :: i, j
707  !---------------------------------------------------------------------------
708 
709  if ( atmos_grid_cartesc_metric_out_basename /= '' ) then
710 
711  log_newline
712  log_info("ATMOS_GRID_CARTESC_METRIC_write",*) 'Output metrics file '
713 
714  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,1,i_xy), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
715  'MAPF_X_XY', 'Map factor x-dir at XY', 'NIL', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
716  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,2,i_xy), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
717  'MAPF_Y_XY', 'Map factor y-dir at XY', 'NIL', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
718  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,1,i_uy), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
719  'MAPF_X_UY', 'Map factor x-dir at UY', 'NIL', 'UY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
720  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,2,i_uy), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
721  'MAPF_Y_UY', 'Map factor y-dir at UY', 'NIL', 'UY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
722  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,1,i_xv), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
723  'MAPF_X_XV', 'Map factor x-dir at XV', 'NIL', 'XV', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
724  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,2,i_xv), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
725  'MAPF_Y_XV', 'Map factor y-dir at XV', 'NIL', 'XV', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
726  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,1,i_uv), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
727  'MAPF_X_UV', 'Map factor x-dir at UV', 'NIL', 'UV', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
728  call file_cartesc_write( atmos_grid_cartesc_metric_mapf(:,:,2,i_uv), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
729  'MAPF_Y_UV', 'Map factor y-dir at UV', 'NIL', 'UV', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
730 
731  call file_cartesc_write( atmos_grid_cartesc_metric_rotc(:,:,1), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
732  'ROTC_COS', 'Rotation factor (cos)', 'NIL', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
733  call file_cartesc_write( atmos_grid_cartesc_metric_rotc(:,:,2), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
734  'ROTC_SIN', 'Rotation factor (sin)', 'NIL', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
735 
736  call file_cartesc_write( atmos_grid_cartesc_metric_rotc(:,:,1), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
737  'ROTC_COS', 'Rotation factor (cos)', 'NIL', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
738 
739  do j = 1, ja
740  do i = 1, ia
741  call mapprojection_lonlat2xy( atmos_grid_cartesc_real_lon(i,j), atmos_grid_cartesc_real_lat(i,j), check_x_xy(i,j), check_y_xy(i,j) )
742  enddo
743  enddo
744 
745  call file_cartesc_write( check_x_xy(:,:), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
746  'X_XY', 'x at XY for check', 'NIL', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
747  call file_cartesc_write( check_y_xy(:,:), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
748  'Y_XY', 'y at XY for check', 'NIL', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
749 
750  do j = 1, ja
751  do i = 1, ia
752  call vectr_distance( const_radius, & ! [IN]
755  atmos_grid_cartesc_real_lon(i,j), & ! [IN]
756  atmos_grid_cartesc_real_lat(i,j), & ! [IN]
757  distance(i,j) ) ! [OUT]
758  enddo
759  enddo
760 
761  call file_cartesc_write( distance(:,:), atmos_grid_cartesc_metric_out_basename, atmos_grid_cartesc_metric_out_title, & ! [IN]
762  'distance', 'distance from basepoint', 'm', 'XY', atmos_grid_cartesc_metric_out_dtype ) ! [IN]
763 
764  endif
765 
766  return
767  end subroutine atmos_grid_cartesc_metric_write
768 
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_gsqrt
transformation metrics from Z to Xi, {G}^1/2
real(rp), public atmos_grid_cartesc_real_basepoint_lon
position of base point in real world [rad,0-2pi]
integer, public ia
of whole cells: x, local, with HALO
subroutine, public mapprojection_rotcoef(IA, IS, IE, JA, JS, JE, lon, lat, rotc)
u(lat,lon) = cos u(x,y) - sin v(x,y) v(lat,lon) = sin u(x,y) + cos v(x,y)
module Atmosphere Grid CartesianC metirc
module vector
subroutine, public atmos_grid_cartesc_metric_setup
Setup.
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_limxz
flux limiter x-z face
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
real(rp), public const_radius
radius of the planet [m]
Definition: scale_const.F90:44
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fzuy
geopotential height [m] (wuy)
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fzuv
geopotential height [m] (wuv)
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fz
geopotential height [m] (wxy)
integer, public ja
of whole cells: y, local, with HALO
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:55
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fz
face coordinate [m]: z, local
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lon
longitude [rad,0-2pi]
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_czxv
geopotential height [m] (zxv)
module COMMUNICATION
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdy
reciprocal of center-dy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cy
center coordinate [m]: y, local
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdz
reciprocal of face-dz
integer, public is
start point of inner domain: x, local
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_czuy
geopotential height [m] (zuy)
integer, public ie
end point of inner domain: x, local
real(rp), public atmos_grid_cartesc_metric_j33g
(3,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdz
reciprocal of center-dz
subroutine, public vectr_distance(r, lon1, lat1, lon2, lat2, dist)
Get horizontal distance on the sphere.
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuy
latitude at staggered point (uy) [rad,-pi,pi]
module atmosphere / grid / cartesC index
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_limxy
flux limiter x-y face
module PROCESS
Definition: scale_prc.F90:11
real(rp), public atmos_grid_cartesc_real_basepoint_lat
position of base point in real world [rad,-pi,pi]
integer, public je
end point of inner domain: y, local
subroutine atmos_grid_cartesc_metric_rotcoef
Calculate rotation coeffient.
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cx
center coordinate [m]: x, local
module atmosphere / grid / cartesC
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_metric_rotc
rotation coefficient
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_limyz
flux limiter y-z face
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_czuv
geopotential height [m] (zuv)
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:338
subroutine, public mapprojection_mapfactor(IA, IS, IE, JA, JS, JE, lat, m1, m2)
(x,y) -> (lon,lat)
module CONSTANT
Definition: scale_const.F90:11
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fy
face coordinate [m]: y, local
integer, public js
start point of inner domain: y, local
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_mapf
map factor
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j23g
(2,3) element of Jacobian matrix * {G}^1/2
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_cz
geopotential height [m] (zxy)
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j13g
(1,3) element of Jacobian matrix * {G}^1/2
module profiler
Definition: scale_prof.F90:11
module Atmosphere GRID CartesC Real(real space)
module Map projection
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdx
reciprocal of center-dx
module PRECISION
module file / cartesianC
real(rp), dimension(:,:), allocatable, public topo_zsfc
absolute ground height [m]
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fx
face coordinate [m]: x, local
integer, public ka
of whole cells: z, local, with HALO
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_fzxv
geopotential height [m] (wxv)
module TOPOGRAPHY
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latuv
latitude at staggered point (uv) [rad,-pi,pi]
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_latxv
latitude at staggered point (xv) [rad,-pi,pi]
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdy
reciprocal of face-dy
module STDIO
Definition: scale_io.F90:10
integer, parameter, public rp
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cz
center coordinate [m]: z, local
subroutine, public atmos_grid_cartesc_real_calc_areavol(MAPF)
Calc control area/volume.
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdx
reciprocal of face-dx