SCALE-RM
Functions/Subroutines | Variables
scale_grid Module Reference

module GRID (cartesian) More...

Functions/Subroutines

subroutine, public grid_setup
 Setup. More...
 
subroutine, public grid_allocate
 Allocate arrays. More...
 
subroutine, public grid_generate
 Generate horizontal&vertical grid. More...
 

Variables

real(rp), public dz = 500.0_RP
 length in the main region [m]: z More...
 
real(rp), public dx = 500.0_RP
 length in the main region [m]: x More...
 
real(rp), public dy = 500.0_RP
 length in the main region [m]: y More...
 
real(rp), public buffer_dz = 0.0_RP
 thickness of buffer region [m]: z More...
 
real(rp), public buffer_dx = 0.0_RP
 thickness of buffer region [m]: x More...
 
real(rp), public buffer_dy = 0.0_RP
 thickness of buffer region [m]: y More...
 
real(rp), public bufffact = 1.0_RP
 strech factor for dx/dy/dz of buffer region More...
 
real(rp), public grid_domain_center_x
 center position of global domain [m]: x More...
 
real(rp), public grid_domain_center_y
 center position of global domain [m]: y More...
 
real(rp), dimension(:), allocatable, public grid_cz
 center coordinate [m]: z, local=global More...
 
real(rp), dimension(:), allocatable, public grid_cx
 center coordinate [m]: x, local More...
 
real(rp), dimension(:), allocatable, public grid_cy
 center coordinate [m]: y, local More...
 
real(rp), dimension(:), allocatable, public grid_cdz
 z-length of control volume [m] More...
 
real(rp), dimension(:), allocatable, public grid_cdx
 x-length of control volume [m] More...
 
real(rp), dimension(:), allocatable, public grid_cdy
 y-length of control volume [m] More...
 
real(rp), dimension(:), allocatable, public grid_rcdz
 reciprocal of center-dz More...
 
real(rp), dimension(:), allocatable, public grid_rcdx
 reciprocal of center-dx More...
 
real(rp), dimension(:), allocatable, public grid_rcdy
 reciprocal of center-dy More...
 
real(rp), dimension(:), allocatable, public grid_fz
 face coordinate [m]: z, local=global More...
 
real(rp), dimension(:), allocatable, public grid_fx
 face coordinate [m]: x, local More...
 
real(rp), dimension(:), allocatable, public grid_fy
 face coordinate [m]: y, local More...
 
real(rp), dimension(:), allocatable, public grid_fdz
 z-length of grid(k+1) to grid(k) [m] More...
 
real(rp), dimension(:), allocatable, public grid_fdx
 x-length of grid(i+1) to grid(i) [m] More...
 
real(rp), dimension(:), allocatable, public grid_fdy
 y-length of grid(j+1) to grid(j) [m] More...
 
real(rp), dimension(:), allocatable, public grid_rfdz
 reciprocal of face-dz More...
 
real(rp), dimension(:), allocatable, public grid_rfdx
 reciprocal of face-dx More...
 
real(rp), dimension(:), allocatable, public grid_rfdy
 reciprocal of face-dy More...
 
real(rp), dimension(:), allocatable, public grid_cbfz
 center buffer factor [0-1]: z More...
 
real(rp), dimension(:), allocatable, public grid_cbfx
 center buffer factor [0-1]: x More...
 
real(rp), dimension(:), allocatable, public grid_cbfy
 center buffer factor [0-1]: y More...
 
real(rp), dimension(:), allocatable, public grid_fbfz
 face buffer factor [0-1]: z More...
 
real(rp), dimension(:), allocatable, public grid_fbfx
 face buffer factor [0-1]: x More...
 
real(rp), dimension(:), allocatable, public grid_fbfy
 face buffer factor [0-1]: y More...
 
real(rp), dimension(:), allocatable, public grid_fxg
 face coordinate [m]: x, global More...
 
real(rp), dimension(:), allocatable, public grid_fyg
 face coordinate [m]: y, global More...
 
real(rp), dimension(:), allocatable, public grid_cxg
 center coordinate [m]: x, global More...
 
real(rp), dimension(:), allocatable, public grid_cyg
 center coordinate [m]: y, global More...
 
real(rp), dimension(:), allocatable, public grid_fbfxg
 face buffer factor [0-1]: x, global More...
 
real(rp), dimension(:), allocatable, public grid_fbfyg
 face buffer factor [0-1]: y, global More...
 
real(rp), dimension(:), allocatable, public grid_cbfxg
 center buffer factor [0-1]: x, global More...
 
real(rp), dimension(:), allocatable, public grid_cbfyg
 center buffer factor [0-1]: y, global More...
 

Detailed Description

module GRID (cartesian)

Description
Grid module for cartesian coordinate
Author
Team SCALE
History
  • 2011-11-11 (H.Yashiro) [new]
  • 2012-03-23 (H.Yashiro) [mod] Explicit index parameter inclusion
  • 2012-06-25 (Y.Sato) [mod] change for unisotropic grid
  • 2012-07-05 (S.Nishizawa) [mod] divided setup into some subroutines

Function/Subroutine Documentation

◆ grid_setup()

subroutine, public scale_grid::grid_setup ( )

Setup.

Definition at line 114 of file scale_grid_cartesian.F90.

References buffer_dx, buffer_dy, buffer_dz, bufffact, dx, dy, dz, grid_allocate(), grid_cx, grid_cy, grid_fx, grid_fy, grid_generate(), scale_grid_index::ia, scale_grid_index::ie, scale_stdio::io_fid_conf, scale_stdio::io_fid_log, scale_stdio::io_l, scale_stdio::io_lnml, scale_grid_index::is, scale_grid_index::ja, scale_grid_index::je, scale_grid_index::js, and scale_process::prc_mpistop().

Referenced by mod_rm_driver::scalerm(), and mod_rm_prep::scalerm_prep().

114  use scale_process, only: &
116  implicit none
117 
118  namelist / param_grid / &
119  grid_in_basename, &
120  grid_out_basename, &
121  grid_offset_x, &
122  grid_offset_y, &
123  dx, &
124  dy, &
125  dz, &
126  buffer_dz, &
127  buffer_dx, &
128  buffer_dy, &
129  bufffact, &
130  fz, &
131  debug
132 
133  integer :: ierr
134  !---------------------------------------------------------------------------
135 
136  if( io_l ) write(io_fid_log,*)
137  if( io_l ) write(io_fid_log,*) '++++++ Module[GRID] / Categ[ATMOS-RM GRID] / Origin[SCALElib]'
138 
139  call grid_allocate
140 
141  !--- read namelist
142  rewind(io_fid_conf)
143  read(io_fid_conf,nml=param_grid,iostat=ierr)
144  if( ierr < 0 ) then !--- missing
145  if( io_l ) write(io_fid_log,*) '*** Not found namelist. Default used.'
146  elseif( ierr > 0 ) then !--- fatal error
147  write(*,*) 'xxx Not appropriate names in namelist PARAM_GRID. Check!'
148  call prc_mpistop
149  endif
150  if( io_lnml ) write(io_fid_log,nml=param_grid)
151 
152  if( io_l ) write(io_fid_log,*)
153  if( io_l ) write(io_fid_log,*) '*** Atmosphere grid information ***'
154  if( io_l ) write(io_fid_log,'(1x,A,3(F7.1))') '*** delta Z, X, Y [m] :', dz, dx, dy
155 
156  if ( grid_in_basename /= '' ) then
157  call grid_read
158  else
159  if( io_l ) write(io_fid_log,*) '*** Not found input grid file. Grid position is calculated.'
160 
161  call grid_generate
162  endif
163 
164  if( io_l ) write(io_fid_log,*)
165  if( io_l ) write(io_fid_log,*) '*** Domain size [km] (local) :'
166  if( io_l ) write(io_fid_log,'(1x,6(A,F8.3))') ' X:', &
167  grid_fx(0) *1.e-3_rp, ' -HALO- ', grid_fx(is-1)*1.e-3_rp, ' | ', &
168  grid_cx(is)*1.e-3_rp, ' - ', grid_cx(ie) *1.e-3_rp, ' | ', &
169  grid_fx(ie)*1.e-3_rp, ' -HALO- ', grid_fx(ia) *1.e-3_rp
170  if( io_l ) write(io_fid_log,'(1x,6(A,F8.3))') ' Y:', &
171  grid_fy(0) *1.e-3_rp, ' -HALO- ', grid_fy(js-1)*1.e-3_rp, ' | ', &
172  grid_cy(js)*1.e-3_rp, ' - ', grid_cy(je) *1.e-3_rp, ' | ', &
173  grid_fy(je)*1.e-3_rp, ' -HALO- ', grid_fy(ja) *1.e-3_rp
174 
175  return
integer, public is
start point of inner domain: x, local
integer, public je
end point of inner domain: y, local
subroutine, public prc_mpistop
Abort MPI.
real(rp), dimension(:), allocatable, public grid_fx
face coordinate [m]: x, local
integer, public ia
of x whole cells (local, with HALO)
integer, public js
start point of inner domain: y, local
module PROCESS
real(rp), dimension(:), allocatable, public grid_cx
center coordinate [m]: x, local
integer, public ie
end point of inner domain: x, local
integer, public io_fid_conf
Config file ID.
Definition: scale_stdio.F90:55
real(rp), dimension(:), allocatable, public grid_cy
center coordinate [m]: y, local
integer, public ja
of y whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local
Here is the call graph for this function:
Here is the caller graph for this function:

◆ grid_allocate()

subroutine, public scale_grid::grid_allocate ( )

Allocate arrays.

Definition at line 181 of file scale_grid_cartesian.F90.

References grid_cbfx, grid_cbfxg, grid_cbfy, grid_cbfyg, grid_cbfz, grid_cdx, grid_cdy, grid_cdz, grid_cx, grid_cxg, grid_cy, grid_cyg, grid_cz, grid_fbfx, grid_fbfxg, grid_fbfy, grid_fbfyg, grid_fbfz, grid_fdx, grid_fdy, grid_fdz, grid_fx, grid_fxg, grid_fy, grid_fyg, grid_fz, grid_rcdx, grid_rcdy, grid_rcdz, grid_rfdx, grid_rfdy, grid_rfdz, scale_grid_index::ia, scale_grid_index::ihalo, scale_grid_index::imax, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::ja, scale_grid_index::jhalo, scale_grid_index::jmax, scale_grid_index::ka, scale_process::prc_myrank, scale_rm_process::prc_num_x, and scale_rm_process::prc_num_y.

Referenced by grid_setup().

181  use scale_rm_process, only: &
182  prc_num_x, &
183  prc_num_y
184  implicit none
185 
186  integer :: iag, jag
187  !---------------------------------------------------------------------------
188 
189  ! working
190  fz(:) = -1.0_rp
191 
192  ! local domain
193  allocate( grid_cz(ka) )
194  allocate( grid_cx(ia) )
195  allocate( grid_cy(ja) )
196  allocate( grid_cdz(ka) )
197  allocate( grid_cdx(ia) )
198  allocate( grid_cdy(ja) )
199  allocate( grid_rcdz(ka) )
200  allocate( grid_rcdx(ia) )
201  allocate( grid_rcdy(ja) )
202 
203  allocate( grid_fz(0:ka) )
204  allocate( grid_fx(0:ia) )
205  allocate( grid_fy(0:ja) )
206  allocate( grid_fdz(ka-1) )
207  allocate( grid_fdx(ia-1) )
208  allocate( grid_fdy(ja-1) )
209  allocate( grid_rfdz(ka-1) )
210  allocate( grid_rfdx(ia-1) )
211  allocate( grid_rfdy(ja-1) )
212 
213  allocate( grid_cbfz(ka) )
214  allocate( grid_cbfx(ia) )
215  allocate( grid_cbfy(ja) )
216  allocate( grid_fbfz(ka) )
217  allocate( grid_fbfx(ia) )
218  allocate( grid_fbfy(ja) )
219 
220  ! array size (global domain)
221  iag = ihalo + imax*prc_num_x + ihalo
222  jag = jhalo + jmax*prc_num_y + jhalo
223 
224  ! global domain
225  allocate( grid_fxg(0:iag) )
226  allocate( grid_fyg(0:jag) )
227  allocate( grid_cxg( iag) )
228  allocate( grid_cyg( jag) )
229  allocate( grid_cbfxg( iag) )
230  allocate( grid_cbfyg( jag) )
231  allocate( grid_fbfxg( iag) )
232  allocate( grid_fbfyg( jag) )
233 
234  return
integer, public imax
of computational cells: x
integer, public prc_num_x
x length of 2D processor topology
real(rp), dimension(:), allocatable, public grid_cyg
center coordinate [m]: y, global
real(rp), dimension(:), allocatable, public grid_cbfyg
center buffer factor [0-1]: y, global
real(rp), dimension(:), allocatable, public grid_cxg
center coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fdy
y-length of grid(j+1) to grid(j) [m]
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fxg
face coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fbfy
face buffer factor [0-1]: y
real(rp), dimension(:), allocatable, public grid_fx
face coordinate [m]: x, local
integer, public prc_num_y
y length of 2D processor topology
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fbfx
face buffer factor [0-1]: x
real(rp), dimension(:), allocatable, public grid_fdz
z-length of grid(k+1) to grid(k) [m]
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fbfxg
face buffer factor [0-1]: x, global
real(rp), dimension(:), allocatable, public grid_fz
face coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fbfz
face buffer factor [0-1]: z
integer, public jhalo
of halo cells: y
real(rp), dimension(:), allocatable, public grid_cbfx
center buffer factor [0-1]: x
real(rp), dimension(:), allocatable, public grid_fbfyg
face buffer factor [0-1]: y, global
real(rp), dimension(:), allocatable, public grid_cbfz
center buffer factor [0-1]: z
real(rp), dimension(:), allocatable, public grid_cx
center coordinate [m]: x, local
module RM PROCESS
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_fdx
x-length of grid(i+1) to grid(i) [m]
real(rp), dimension(:), allocatable, public grid_cdy
y-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_cbfy
center buffer factor [0-1]: y
real(rp), dimension(:), allocatable, public grid_cdx
x-length of control volume [m]
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_cy
center coordinate [m]: y, local
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)
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local
Here is the caller graph for this function:

◆ grid_generate()

subroutine, public scale_grid::grid_generate ( )

Generate horizontal&vertical grid.

Definition at line 293 of file scale_grid_cartesian.F90.

References buffer_dx, buffer_dy, buffer_dz, bufffact, dx, dy, dz, grid_cbfx, grid_cbfxg, grid_cbfy, grid_cbfyg, grid_cbfz, grid_cdx, grid_cdy, grid_cdz, grid_cx, grid_cxg, grid_cy, grid_cyg, grid_cz, grid_domain_center_x, grid_domain_center_y, grid_fbfx, grid_fbfxg, grid_fbfy, grid_fbfyg, grid_fbfz, grid_fdx, grid_fdy, grid_fdz, grid_fx, grid_fxg, grid_fy, grid_fyg, grid_fz, grid_rcdx, grid_rcdy, grid_rcdz, grid_rfdx, grid_rfdy, grid_rfdz, scale_grid_index::ia, scale_grid_index::ihalo, scale_grid_index::imax, scale_stdio::io_fid_log, scale_stdio::io_l, scale_grid_index::ja, scale_grid_index::jhalo, scale_grid_index::jmax, scale_grid_index::ka, scale_grid_index::ke, scale_grid_index::kmax, scale_grid_index::ks, scale_rm_process::prc_2drank, scale_process::prc_mpistop(), scale_process::prc_myrank, scale_rm_process::prc_num_x, and scale_rm_process::prc_num_y.

Referenced by grid_setup().

293  use scale_process, only: &
294  prc_mpistop, &
295  prc_myrank
296  use scale_rm_process, only: &
297  prc_2drank, &
298  prc_num_x, &
299  prc_num_y
300  implicit none
301 
302  integer :: iag ! # of x whole cells (global, with HALO)
303  integer :: jag ! # of y whole cells (global, with HALO)
304 
305  real(RP), allocatable :: buffz(:), buffx(:), buffy(:)
306  real(RP) :: bufftotz, bufftotx, bufftoty
307 
308  integer :: kbuff, ibuff, jbuff
309  integer :: kmain, imain, jmain
310 
311  integer :: k, i, j, ii, jj
312  !---------------------------------------------------------------------------
313 
314  !##### coordinate in global domain #####
315 
316  ! array size (global domain)
317  iag = ihalo + imax*prc_num_x + ihalo
318  jag = jhalo + jmax*prc_num_y + jhalo
319 
320  allocate( buffx(0:iag) )
321  allocate( buffy(0:jag) )
322 
323  ! X-direction
324  ! calculate buffer grid size
325  buffx(0) = dx
326  bufftotx = 0.0_rp
327 
328  do i = 1, iag
329  if( bufftotx >= buffer_dx ) exit
330 
331  buffx(i) = buffx(i-1) * bufffact
332  bufftotx = bufftotx + buffx(i)
333  enddo
334  ibuff = i - 1
335  imain = iag - 2*ibuff - 2*ihalo
336 
337  if ( imain < 1 ) then
338  write(*,*) 'xxx Buffer size (', bufftotx*2.0_rp, ') must be smaller than global domain size (X). Use smaller BUFFER_DX!'
339  call prc_mpistop
340  endif
341 
342  ! horizontal coordinate (global domaim)
343  grid_fxg(ihalo) = grid_offset_x
344  do i = ihalo-1, 0, -1
345  grid_fxg(i) = grid_fxg(i+1) - buffx(ibuff)
346  enddo
347 
348  do i = 1, ihalo
349  grid_cxg(i) = 0.5_rp * ( grid_fxg(i)+grid_fxg(i-1) )
350  enddo
351 
352  if ( ibuff > 0 ) then
353  do i = ihalo+1, ihalo+ibuff
354  grid_fxg(i) = grid_fxg(i-1) + buffx(ibuff+ihalo+1-i)
355  grid_cxg(i) = 0.5_rp * ( grid_fxg(i)+grid_fxg(i-1) )
356  enddo
357  endif
358 
359  do i = ihalo+ibuff+1, ihalo+ibuff+imain
360  grid_fxg(i) = grid_fxg(i-1) + dx
361  grid_cxg(i) = 0.5_rp * ( grid_fxg(i)+grid_fxg(i-1) )
362  enddo
363 
364  if ( ibuff > 0 ) then
365  do i = ihalo+ibuff+imain+1, ihalo+ibuff+imain+ibuff
366  grid_fxg(i) = grid_fxg(i-1) + buffx(i-ihalo-ibuff-imain)
367  grid_cxg(i) = 0.5_rp * ( grid_fxg(i)+grid_fxg(i-1) )
368  enddo
369  endif
370 
371  do i = ihalo+ibuff+imain+ibuff+1, ihalo+ibuff+imain+ibuff+ihalo
372  grid_fxg(i) = grid_fxg(i-1) + buffx(ibuff)
373  grid_cxg(i) = 0.5_rp * ( grid_fxg(i)+grid_fxg(i-1) )
374  enddo
375 
376  ! calc buffer factor (global domaim)
377  grid_cbfxg(:) = 0.0_rp
378  grid_fbfxg(:) = 0.0_rp
379  do i = 1, ihalo
380  grid_cbfxg(i) = 1.0_rp
381  grid_fbfxg(i) = 1.0_rp
382  enddo
383 
384  if ( ibuff > 0 ) then
385  do i = ihalo+1, ihalo+ibuff
386  grid_cbfxg(i) = (bufftotx+grid_fxg(ihalo)-grid_cxg(i)) / bufftotx
387  grid_fbfxg(i) = (bufftotx+grid_fxg(ihalo)-grid_fxg(i)) / bufftotx
388  enddo
389 
390  do i = ihalo+ibuff+imain+1, ihalo+ibuff+imain+ibuff
391  grid_cbfxg(i) = (bufftotx-grid_fxg(iag-ihalo)+grid_cxg(i)) / bufftotx
392  grid_fbfxg(i) = (bufftotx-grid_fxg(iag-ihalo)+grid_fxg(i)) / bufftotx
393  enddo
394  endif
395 
396  do i = ihalo+ibuff+imain+ibuff+1, ihalo+ibuff+imain+ibuff+ihalo
397  grid_cbfxg(i) = 1.0_rp
398  grid_fbfxg(i) = 1.0_rp
399  enddo
400 
401  grid_cbfxg(:) = max( min( grid_cbfxg(:), 1.0_rp ), 0.0_rp )
402  grid_fbfxg(:) = max( min( grid_fbfxg(:), 1.0_rp ), 0.0_rp )
403 
404  ! Y-direction
405  ! calculate buffer grid size
406  buffy(0) = dy
407  bufftoty = 0.0_rp
408 
409  do j = 1, jag
410  if( bufftoty >= buffer_dy ) exit
411 
412  buffy(j) = buffy(j-1) * bufffact
413  bufftoty = bufftoty + buffy(j)
414  enddo
415  jbuff = j - 1
416  jmain = jag - 2*jbuff - 2*jhalo
417 
418  if ( jmain < 1 ) then
419  write(*,*) 'xxx Buffer size (', bufftoty*2.0_rp, ') must be smaller than global domain size (Y). Use smaller BUFFER_DY!'
420  call prc_mpistop
421  endif
422 
423  ! horizontal coordinate (global domaim)
424  grid_fyg(jhalo) = grid_offset_y
425  do j = jhalo-1, 0, -1
426  grid_fyg(j) = grid_fyg(j+1) - buffy(jbuff)
427  enddo
428 
429  do j = 1, jhalo
430  grid_cyg(j) = 0.5_rp * ( grid_fyg(j)+grid_fyg(j-1) )
431  enddo
432 
433  if ( jbuff > 0 ) then
434  do j = jhalo+1, jhalo+jbuff
435  grid_fyg(j) = grid_fyg(j-1) + buffy(jbuff+jhalo+1-j)
436  grid_cyg(j) = 0.5_rp * ( grid_fyg(j)+grid_fyg(j-1) )
437  enddo
438  endif
439 
440  do j = jhalo+jbuff+1, jhalo+jbuff+jmain
441  grid_fyg(j) = grid_fyg(j-1) + dy
442  grid_cyg(j) = 0.5_rp * ( grid_fyg(j)+grid_fyg(j-1) )
443  enddo
444 
445  if ( jbuff > 0 ) then
446  do j = jhalo+jbuff+jmain+1, jhalo+jbuff+jmain+jbuff
447  grid_fyg(j) = grid_fyg(j-1) + buffy(j-jhalo-jbuff-jmain)
448  grid_cyg(j) = 0.5_rp * ( grid_fyg(j)+grid_fyg(j-1) )
449  enddo
450  endif
451 
452  do j = jhalo+jbuff+jmain+jbuff+1, jhalo+jbuff+jmain+jbuff+jhalo
453  grid_fyg(j) = grid_fyg(j-1) + buffy(jbuff)
454  grid_cyg(j) = 0.5_rp * ( grid_fyg(j)+grid_fyg(j-1) )
455  enddo
456 
457  ! calc buffer factor (global domaim)
458  grid_cbfyg(:) = 0.0_rp
459  grid_fbfyg(:) = 0.0_rp
460  do j = 1, jhalo
461  grid_cbfyg(j) = 1.0_rp
462  grid_fbfyg(j) = 1.0_rp
463  enddo
464 
465  if ( jbuff > 0 ) then
466  do j = jhalo+1, jhalo+jbuff
467  grid_cbfyg(j) = (bufftoty+grid_fyg(jhalo)-grid_cyg(j)) / bufftoty
468  grid_fbfyg(j) = (bufftoty+grid_fyg(jhalo)-grid_fyg(j)) / bufftoty
469  enddo
470 
471  do j = jhalo+jbuff+jmain+1, jhalo+jbuff+jmain+jbuff
472  grid_cbfyg(j) = (bufftoty-grid_fyg(jag-jhalo)+grid_cyg(j)) / bufftoty
473  grid_fbfyg(j) = (bufftoty-grid_fyg(jag-jhalo)+grid_fyg(j)) / bufftoty
474  enddo
475  endif
476 
477  do j = jhalo+jbuff+jmain+jbuff+1, jhalo+jbuff+jmain+jbuff+jhalo
478  grid_cbfyg(j) = 1.0_rp
479  grid_fbfyg(j) = 1.0_rp
480  enddo
481  grid_cbfyg(:) = max( min( grid_cbfyg(:), 1.0_rp ), 0.0_rp )
482  grid_fbfyg(:) = max( min( grid_fbfyg(:), 1.0_rp ), 0.0_rp )
483 
484  deallocate( buffx )
485  deallocate( buffy )
486 
487  !##### coordinate in local domain #####
488 
489  allocate( buffz(0:ka) )
490 
491  if ( minval(fz(1:kmax)) > 0.0_rp ) then ! input from namelist
492  if( io_l ) write(io_fid_log,*) '*** Z coordinate is given from NAMELIST.'
493 
494  if ( kmax < 2 ) then
495  write(*,*) 'xxx If you use FZ, KMAX must be larger than 1. Check!', kmax
496  call prc_mpistop
497  endif
498 
499  ! Z-direction
500  ! calculate buffer grid size
501  bufftotz = 0.0_rp
502 
503  do k = kmax, 2, -1
504  if( bufftotz >= buffer_dz ) exit
505 
506  bufftotz = bufftotz + ( fz(k) - fz(k-1) )
507  enddo
508  kbuff = kmax - k
509  kmain = k
510 
511  if ( kmain == 1 ) then
512  write(*,*) 'xxx Buffer size (', bufftotz, ') must be smaller than domain size (z). Use smaller BUFFER_DZ!'
513  call prc_mpistop
514  endif
515 
516  ! vartical coordinate (local=global domaim)
517  grid_fz(ks-1) = 0.0_rp
518 
519  dz = fz(1)
520  do k = ks-2, 0, -1
521  grid_fz(k) = grid_fz(k+1) - dz
522  enddo
523 
524  do k = ks, ke
525  grid_fz(k) = fz(k-ks+1)
526  enddo
527 
528  dz = fz(kmax) - fz(kmax-1)
529  do k = ke+1, ka
530  grid_fz(k) = grid_fz(k-1) + dz
531  enddo
532 
533  do k = 1, ka
534  grid_cz(k) = 0.5_rp * ( grid_fz(k)+grid_fz(k-1) )
535  enddo
536 
537  else ! calc using DZ
538 
539  ! Z-direction
540  ! calculate buffer grid size
541  buffz(0) = dz
542  bufftotz = 0.0_rp
543 
544  do k = 1, ka
545  if( bufftotz >= buffer_dz ) exit
546 
547  buffz(k) = buffz(k-1) * bufffact
548  bufftotz = bufftotz + buffz(k)
549  enddo
550  kbuff = k - 1
551  kmain = ke - ks + 1 - kbuff
552 
553  if ( kmain < 1 ) then
554  write(*,*) 'xxx Buffer size (', bufftotz, ') must be smaller than domain size (z). Use smaller BUFFER_DZ!'
555  call prc_mpistop
556  endif
557 
558  ! vartical coordinate (local=global domaim)
559  grid_fz(ks-1) = 0.0_rp
560  do k = ks-2, 0, -1
561  grid_fz(k) = grid_fz(k+1) - dz
562  enddo
563 
564  do k = 1, ks-1
565  grid_cz(k) = 0.5_rp * ( grid_fz(k)+grid_fz(k-1) )
566  enddo
567 
568  do k = ks, ks+kmain-1
569  grid_fz(k) = grid_fz(k-1) + dz
570  grid_cz(k) = 0.5_rp * ( grid_fz(k)+grid_fz(k-1) )
571  enddo
572 
573  if ( kbuff > 0 ) then
574  do k = ks+kmain, ke
575  grid_fz(k) = grid_fz(k-1) + buffz(k-ks-kmain+1)
576  grid_cz(k) = 0.5_rp * ( grid_fz(k)+grid_fz(k-1) )
577  enddo
578  endif
579 
580  do k = ke+1, ka
581  grid_fz(k) = grid_fz(k-1) + buffz(kbuff)
582  grid_cz(k) = 0.5_rp * ( grid_fz(k)+grid_fz(k-1) )
583  enddo
584 
585  endif
586 
587  ! calc buffer factor (global domaim)
588  grid_cbfz(:) = 0.0_rp
589  grid_fbfz(:) = 0.0_rp
590  if ( kbuff > 0 ) then
591  do k = ks+kmain, ke
592  grid_cbfz(k) = (bufftotz-grid_fz(ke)+grid_cz(k)) / bufftotz
593  grid_fbfz(k) = (bufftotz-grid_fz(ke)+grid_fz(k)) / bufftotz
594  enddo
595  endif
596 
597  do k = ke+1, ka
598  grid_cbfz(k) = 1.0_rp
599  grid_fbfz(k) = 1.0_rp
600  enddo
601  grid_cbfz(:) = max( min( grid_cbfz(:), 1.0_rp ), 0.0_rp )
602  grid_fbfz(:) = max( min( grid_fbfz(:), 1.0_rp ), 0.0_rp )
603 
604  deallocate( buffz )
605 
606  ! vartical coordinate (local domaim)
607  do k = 1, ka
608  grid_cdz(k) = grid_fz(k) - grid_fz(k-1)
609  grid_rcdz(k) = 1.0_rp / grid_cdz(k)
610  enddo
611 
612  do k = 1, ka-1
613  grid_fdz(k) = grid_cz(k+1)-grid_cz(k)
614  grid_rfdz(k) = 1.0_rp / grid_fdz(k)
615  enddo
616 
617  ! X-direction
618  ! horizontal coordinate (local domaim)
619  do i = 0, ia
620  ii = i + prc_2drank(prc_myrank,1) * imax
621 
622  grid_fx(i) = grid_fxg(ii)
623  enddo
624 
625  do i = 1, ia
626  ii = i + prc_2drank(prc_myrank,1) * imax
627 
628  grid_cx(i) = grid_cxg(ii)
629  grid_cbfx(i) = grid_cbfxg(ii)
630  grid_fbfx(i) = grid_fbfxg(ii)
631 
632  grid_cdx(i) = grid_fx(i) - grid_fx(i-1)
633  grid_rcdx(i) = 1.0_rp / grid_cdx(i)
634  enddo
635 
636  do i = 1, ia-1
637  grid_fdx(i) = grid_cx(i+1)-grid_cx(i)
638  grid_rfdx(i) = 1.0_rp / grid_fdx(i)
639  enddo
640 
641  ! Y-direction
642  ! horizontal coordinate (local domaim)
643  do j = 0, ja
644  jj = j + prc_2drank(prc_myrank,2) * jmax
645 
646  grid_fy(j) = grid_fyg(jj)
647  enddo
648 
649  do j = 1, ja
650  jj = j + prc_2drank(prc_myrank,2) * jmax
651 
652  grid_cy(j) = grid_cyg(jj)
653  grid_cbfy(j) = grid_cbfyg(jj)
654  grid_fbfy(j) = grid_fbfyg(jj)
655 
656  grid_cdy(j) = grid_fy(j) - grid_fy(j-1)
657  grid_rcdy(j) = 1.0_rp / grid_cdy(j)
658  enddo
659 
660  do j = 1, ja-1
661  grid_fdy(j) = grid_cy(j+1)-grid_cy(j)
662  grid_rfdy(j) = 1.0_rp / grid_fdy(j)
663  enddo
664 
665  grid_domain_center_x = 0.5_rp * ( grid_fxg(0) + grid_fxg(iag) )
666  grid_domain_center_y = 0.5_rp * ( grid_fyg(0) + grid_fyg(jag) )
667 
668  ! report
669  if( io_l ) write(io_fid_log,*)
670  if( io_l ) write(io_fid_log,*) '*** Main/buffer Grid (global) :'
671  if( io_l ) write(io_fid_log,'(1x,2(A,I6))') ' Z: buffer = ', kbuff,' x 1, main = ',kmain
672  if( io_l ) write(io_fid_log,'(1x,2(A,I6))') ' X: buffer = ', ibuff,' x 2, main = ',imain
673  if( io_l ) write(io_fid_log,'(1x,2(A,I6))') ' Y: buffer = ', jbuff,' x 2, main = ',jmain
674  if( io_l ) write(io_fid_log,*)
675  if( io_l ) write(io_fid_log,*) '*** Domain size [km] (global) :'
676  if( io_l ) write(io_fid_log,'(1x,7(A,F8.3))') ' Z:', &
677  grid_fz(0) *1.e-3_rp, ' -HALO- ', &
678  grid_fz(ks-1) *1.e-3_rp, ' | ', &
679  grid_cz(ks) *1.e-3_rp, ' - ', &
680  grid_cz(ke-kbuff)*1.e-3_rp, ' | ', &
681  grid_fz(ke-kbuff)*1.e-3_rp, ' -buffer- ', &
682  grid_fz(ke) *1.e-3_rp, ' -HALO- ', &
683  grid_fz(ka) *1.e-3_rp
684  if( io_l ) write(io_fid_log,'(1x,8(A,F8.3))') ' X:', &
685  grid_fxg(0) *1.e-3_rp, ' -HALO- ', &
686  grid_fxg(ihalo) *1.e-3_rp, ' -buffer- ', &
687  grid_fxg(ihalo+ibuff) *1.e-3_rp, ' | ', &
688  grid_cxg(ihalo+ibuff+1) *1.e-3_rp, ' - ', &
689  grid_cxg(iag-ihalo-ibuff)*1.e-3_rp, ' | ', &
690  grid_fxg(iag-ihalo-ibuff)*1.e-3_rp, ' -buffer- ', &
691  grid_fxg(iag-ihalo) *1.e-3_rp, ' -HALO- ', &
692  grid_fxg(iag) *1.e-3_rp
693  if( io_l ) write(io_fid_log,'(1x,8(A,F8.3))') ' Y:', &
694  grid_fyg(0) *1.e-3_rp, ' -HALO- ', &
695  grid_fyg(jhalo) *1.e-3_rp, ' -buffer- ', &
696  grid_fyg(jhalo+jbuff) *1.e-3_rp, ' | ', &
697  grid_cyg(jhalo+jbuff+1) *1.e-3_rp, ' - ', &
698  grid_cyg(jag-jhalo-jbuff)*1.e-3_rp, ' | ', &
699  grid_fyg(jag-jhalo-jbuff)*1.e-3_rp, ' -buffer- ', &
700  grid_fyg(jag-jhalo) *1.e-3_rp, ' -HALO- ', &
701  grid_fyg(jag) *1.e-3_rp
702  if( io_l ) write(io_fid_log,*)
703  if( io_l ) write(io_fid_log,*) '*** Center Position of Grid (global) :'
704  if( io_l ) write(io_fid_log,'(1x,A,F12.3)') ' X: ', grid_domain_center_x
705  if( io_l ) write(io_fid_log,'(1x,A,F12.3)') ' Y: ', grid_domain_center_y
706 
707  if ( debug ) then
708  if( io_l ) write(io_fid_log,*)
709  if( io_l ) write(io_fid_log,'(1x,A)') &
710  '|============= Vertical Coordinate =============|'
711  if( io_l ) write(io_fid_log,'(1x,A)') &
712  '| k z zh dz buffer k |'
713  if( io_l ) write(io_fid_log,'(1x,A)') &
714  '| [m] [m] [m] factor |'
715 
716  do k = ka, ke+1, -1
717  if( io_l ) write(io_fid_log,'(1x,A,F9.2,A,F9.2,I5,A)') &
718  '| ',grid_fz(k),' ', grid_fbfz(k),k,' |'
719  if( io_l ) write(io_fid_log,'(1x,A,I5,F9.2,A,2F9.2,A)') &
720  '|',k,grid_cz(k),' ',grid_cdz(k), grid_cbfz(k),' |'
721  enddo
722 
723  k = ke
724  if( io_l ) write(io_fid_log,'(1x,A,F9.2,A,F9.2,I5,A)') &
725  '| ',grid_fz(k),' ', grid_fbfz(k),k,' | KE = TOA'
726  if( io_l ) write(io_fid_log,'(1x,A,I5,F9.2,A,2F9.2,A)') &
727  '|',k,grid_cz(k),' ',grid_cdz(k), grid_cbfz(k),' |'
728 
729  do k = ke-1, ks, -1
730  if( io_l ) write(io_fid_log,'(1x,A,F9.2,A,F9.2,I5,A)') &
731  '| ',grid_fz(k),' ', grid_fbfz(k),k,' |'
732  if( io_l ) write(io_fid_log,'(1x,A,I5,F9.2,A,2F9.2,A)') &
733  '|',k,grid_cz(k),' ',grid_cdz(k), grid_cbfz(k),' |'
734  enddo
735 
736  k = ks-1
737  if( io_l ) write(io_fid_log,'(1x,A,F9.2,A,F9.2,I5,A)') &
738  '| ',grid_fz(k),' ', grid_fbfz(k),k,' | KS-1 = surface'
739  if( io_l ) write(io_fid_log,'(1x,A,I5,F9.2,A,2F9.2,A)') &
740  '|',k,grid_cz(k),' ',grid_cdz(k), grid_cbfz(k),' |'
741 
742  do k = ks-2, 1, -1
743  if( io_l ) write(io_fid_log,'(1x,A,F9.2,A,F9.2,I5,A)') &
744  '| ',grid_fz(k),' ', grid_fbfz(k),k,' |'
745  if( io_l ) write(io_fid_log,'(1x,A,I5,F9.2,A,2F9.2,A)') &
746  '|',k,grid_cz(k),' ',grid_cdz(k), grid_cbfz(k),' |'
747  enddo
748 
749  k = 0
750  if( io_l ) write(io_fid_log,'(1x,A,F9.2,A,F9.2,I5,A)') &
751  '| ',grid_fz(k),' ', grid_fbfz(k),k,' |'
752 
753  if( io_l ) write(io_fid_log,'(1x,A)') &
754  '|===============================================|'
755 
756 ! if( IO_L ) write(IO_FID_LOG,*)
757 ! if( IO_L ) write(IO_FID_LOG,*) ' ', 0, GRID_FZ(0)
758 ! do k = 1, KA-1
759 ! if( IO_L ) write(IO_FID_LOG,*) k, GRID_CZ(k), GRID_CBFZ(k), GRID_CDZ(k)
760 ! if( IO_L ) write(IO_FID_LOG,*) ' ', k, GRID_FZ(k), GRID_FBFZ(k), GRID_FDZ(k)
761 ! enddo
762 ! k = KA
763 ! if( IO_L ) write(IO_FID_LOG,*) k, GRID_CZ(k), GRID_CBFZ(k), GRID_CDZ(k)
764 ! if( IO_L ) write(IO_FID_LOG,*) ' ', k, GRID_FZ(k), GRID_FBFZ(k)
765 !
766 ! if( IO_L ) write(IO_FID_LOG,*)
767 ! if( IO_L ) write(IO_FID_LOG,*) ' ', 0, GRID_FX(0)
768 ! do i = 1, IA-1
769 ! if( IO_L ) write(IO_FID_LOG,*) i, GRID_CX(i), GRID_CBFX(i), GRID_CDX(i)
770 ! if( IO_L ) write(IO_FID_LOG,*) ' ', i, GRID_FX(i), GRID_FBFX(i), GRID_FDX(i)
771 ! enddo
772 ! i = IA
773 ! if( IO_L ) write(IO_FID_LOG,*) i, GRID_CX(i), GRID_CBFX(i), GRID_CDX(i)
774 ! if( IO_L ) write(IO_FID_LOG,*) ' ', i, GRID_FX(i), GRID_FBFX(i)
775 !
776 ! if( IO_L ) write(IO_FID_LOG,*)
777 ! if( IO_L ) write(IO_FID_LOG,*) ' ', 0, GRID_FY(0)
778 ! do j = 1, JA-1
779 ! if( IO_L ) write(IO_FID_LOG,*) j, GRID_CY(j), GRID_CBFY(j), GRID_CDY(j)
780 ! if( IO_L ) write(IO_FID_LOG,*) ' ', j, GRID_FY(j), GRID_FBFY(j), GRID_FDY(j)
781 ! enddo
782 ! j = JA
783 ! if( IO_L ) write(IO_FID_LOG,*) j, GRID_CY(j), GRID_CBFY(j), GRID_CDY(j)
784 ! if( IO_L ) write(IO_FID_LOG,*) ' ', j, GRID_FY(j), GRID_FBFY(j)
785  endif
786 
787  return
integer, public imax
of computational cells: x
integer, public prc_num_x
x length of 2D processor topology
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
real(rp), dimension(:), allocatable, public grid_cxg
center coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fdy
y-length of grid(j+1) to grid(j) [m]
real(rp), dimension(:), allocatable, public grid_cz
center coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fxg
face coordinate [m]: x, global
real(rp), dimension(:), allocatable, public grid_fbfy
face buffer factor [0-1]: y
integer, public ke
end point of inner domain: z, local
real(rp), dimension(:), allocatable, public grid_fx
face coordinate [m]: x, local
integer, public prc_num_y
y length of 2D processor topology
integer, public ia
of x whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fbfx
face buffer factor [0-1]: x
real(rp), dimension(:), allocatable, public grid_fdz
z-length of grid(k+1) to grid(k) [m]
integer, public ka
of z whole cells (local, with HALO)
real(rp), dimension(:), allocatable, public grid_fbfxg
face buffer factor [0-1]: x, global
integer, public kmax
of computational cells: z
real(rp), dimension(:), allocatable, public grid_fz
face coordinate [m]: z, local=global
real(rp), dimension(:), allocatable, public grid_fbfz
face buffer factor [0-1]: z
integer, public jhalo
of halo cells: y
real(rp), dimension(:), allocatable, public grid_cbfx
center buffer factor [0-1]: x
module PROCESS
real(rp), dimension(:), allocatable, public grid_fbfyg
face buffer factor [0-1]: y, global
integer, public ks
start point of inner domain: z, local
real(rp), dimension(:), allocatable, public grid_cbfz
center buffer factor [0-1]: z
integer, public prc_myrank
process num in local communicator
real(rp), dimension(:), allocatable, public grid_cx
center coordinate [m]: x, local
module RM PROCESS
integer, dimension(:,:), allocatable, public prc_2drank
node index in 2D topology
real(rp), dimension(:), allocatable, public grid_cdz
z-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_fdx
x-length of grid(i+1) to grid(i) [m]
real(rp), dimension(:), allocatable, public grid_cdy
y-length of control volume [m]
real(rp), dimension(:), allocatable, public grid_cbfy
center buffer factor [0-1]: y
real(rp), dimension(:), allocatable, public grid_cdx
x-length of control volume [m]
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_cy
center coordinate [m]: y, local
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)
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local
Here is the call graph for this function:
Here is the caller graph for this function:

Variable Documentation

◆ dz

real(rp), public scale_grid::dz = 500.0_RP

length in the main region [m]: z

Definition at line 41 of file scale_grid_cartesian.F90.

Referenced by scale_atmos_phy_mp_sn14::atmos_phy_mp_sn14(), grid_generate(), grid_setup(), and scale_grid_real::real_calc_areavol().

41  real(RP), public :: dz = 500.0_rp

◆ dx

real(rp), public scale_grid::dx = 500.0_RP

length in the main region [m]: x

Definition at line 42 of file scale_grid_cartesian.F90.

Referenced by mod_cnvtopo::cnvtopo(), mod_cnvtopo::cnvtopo_setup(), mod_copytopo::copytopo(), grid_generate(), grid_setup(), and scale_grid_real::real_calc_areavol().

42  real(RP), public :: dx = 500.0_rp

◆ dy

real(rp), public scale_grid::dy = 500.0_RP

length in the main region [m]: y

Definition at line 43 of file scale_grid_cartesian.F90.

Referenced by mod_cnvtopo::cnvtopo(), mod_cnvtopo::cnvtopo_setup(), mod_copytopo::copytopo(), grid_generate(), grid_setup(), and scale_grid_real::real_calc_areavol().

43  real(RP), public :: dy = 500.0_rp

◆ buffer_dz

real(rp), public scale_grid::buffer_dz = 0.0_RP

thickness of buffer region [m]: z

Definition at line 45 of file scale_grid_cartesian.F90.

Referenced by grid_generate(), and grid_setup().

45  real(RP), public :: buffer_dz = 0.0_rp

◆ buffer_dx

real(rp), public scale_grid::buffer_dx = 0.0_RP

thickness of buffer region [m]: x

Definition at line 46 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), grid_generate(), and grid_setup().

46  real(RP), public :: buffer_dx = 0.0_rp

◆ buffer_dy

real(rp), public scale_grid::buffer_dy = 0.0_RP

thickness of buffer region [m]: y

Definition at line 47 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), grid_generate(), and grid_setup().

47  real(RP), public :: buffer_dy = 0.0_rp

◆ bufffact

real(rp), public scale_grid::bufffact = 1.0_RP

strech factor for dx/dy/dz of buffer region

Definition at line 48 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), grid_generate(), and grid_setup().

48  real(RP), public :: bufffact = 1.0_rp

◆ grid_domain_center_x

real(rp), public scale_grid::grid_domain_center_x

center position of global domain [m]: x

Definition at line 50 of file scale_grid_cartesian.F90.

Referenced by grid_generate(), scale_grid_real::real_setup(), and scale_grid_real::real_update_z().

50  real(RP), public :: grid_domain_center_x

◆ grid_domain_center_y

real(rp), public scale_grid::grid_domain_center_y

center position of global domain [m]: y

Definition at line 51 of file scale_grid_cartesian.F90.

Referenced by grid_generate(), scale_grid_real::real_setup(), and scale_grid_real::real_update_z().

51  real(RP), public :: grid_domain_center_y

◆ grid_cz

real(rp), dimension (:), allocatable, public scale_grid::grid_cz

◆ grid_cx

real(rp), dimension (:), allocatable, public scale_grid::grid_cx

◆ grid_cy

real(rp), dimension (:), allocatable, public scale_grid::grid_cy

center coordinate [m]: y, local

Definition at line 55 of file scale_grid_cartesian.F90.

Referenced by scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), grid_setup(), scale_gridtrans::gtrans_rotcoef(), scale_history::hist_switch(), mod_mkinit::mkinit(), mod_mkinit::read_sounding(), scale_grid_real::real_update_z(), mod_mkinit::rect_setup(), and scale_topography::topo_fillhalo().

55  real(RP), public, allocatable :: grid_cy (:)
real(rp), dimension(:), allocatable, public grid_cy
center coordinate [m]: y, local

◆ grid_cdz

real(rp), dimension (:), allocatable, public scale_grid::grid_cdz

◆ grid_cdx

real(rp), dimension (:), allocatable, public scale_grid::grid_cdx

x-length of control volume [m]

Definition at line 57 of file scale_grid_cartesian.F90.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver(), mod_atmos_dyn_driver::atmos_dyn_driver_setup(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver_setup(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), scale_history::hist_switch(), and scale_grid_real::real_update_z().

57  real(RP), public, allocatable :: grid_cdx (:)
real(rp), dimension(:), allocatable, public grid_cdx
x-length of control volume [m]

◆ grid_cdy

real(rp), dimension (:), allocatable, public scale_grid::grid_cdy

y-length of control volume [m]

Definition at line 58 of file scale_grid_cartesian.F90.

Referenced by mod_atmos_dyn_driver::atmos_dyn_driver(), mod_atmos_dyn_driver::atmos_dyn_driver_setup(), mod_atmos_phy_tb_driver::atmos_phy_tb_driver_setup(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), scale_history::hist_switch(), and scale_grid_real::real_update_z().

58  real(RP), public, allocatable :: grid_cdy (:)
real(rp), dimension(:), allocatable, public grid_cdy
y-length of control volume [m]

◆ grid_rcdz

real(rp), dimension(:), allocatable, public scale_grid::grid_rcdz

◆ grid_rcdx

real(rp), dimension(:), allocatable, public scale_grid::grid_rcdx

◆ grid_rcdy

real(rp), dimension(:), allocatable, public scale_grid::grid_rcdy

◆ grid_fz

real(rp), dimension (:), allocatable, public scale_grid::grid_fz

◆ grid_fx

real(rp), dimension (:), allocatable, public scale_grid::grid_fx

face coordinate [m]: x, local

Definition at line 64 of file scale_grid_cartesian.F90.

Referenced by scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), grid_setup(), scale_gridtrans::gtrans_rotcoef(), scale_history::hist_switch(), mod_mkinit::read_sounding(), and scale_grid_real::real_update_z().

64  real(RP), public, allocatable :: grid_fx (:)
real(rp), dimension(:), allocatable, public grid_fx
face coordinate [m]: x, local

◆ grid_fy

real(rp), dimension (:), allocatable, public scale_grid::grid_fy

face coordinate [m]: y, local

Definition at line 65 of file scale_grid_cartesian.F90.

Referenced by scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), grid_setup(), scale_gridtrans::gtrans_rotcoef(), scale_history::hist_switch(), mod_mkinit::read_sounding(), and scale_grid_real::real_update_z().

65  real(RP), public, allocatable :: grid_fy (:)
real(rp), dimension(:), allocatable, public grid_fy
face coordinate [m]: y, local

◆ grid_fdz

real(rp), dimension (:), allocatable, public scale_grid::grid_fdz

◆ grid_fdx

real(rp), dimension (:), allocatable, public scale_grid::grid_fdx

◆ grid_fdy

real(rp), dimension (:), allocatable, public scale_grid::grid_fdy

◆ grid_rfdz

real(rp), dimension(:), allocatable, public scale_grid::grid_rfdz

◆ grid_rfdx

real(rp), dimension(:), allocatable, public scale_grid::grid_rfdx

◆ grid_rfdy

real(rp), dimension(:), allocatable, public scale_grid::grid_rfdy

◆ grid_cbfz

real(rp), dimension(:), allocatable, public scale_grid::grid_cbfz

center buffer factor [0-1]: z

Definition at line 73 of file scale_grid_cartesian.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), scale_history::hist_switch(), and scale_land_grid::land_grid_setup().

73  real(RP), public, allocatable :: grid_cbfz(:)
real(rp), dimension(:), allocatable, public grid_cbfz
center buffer factor [0-1]: z

◆ grid_cbfx

real(rp), dimension(:), allocatable, public scale_grid::grid_cbfx

center buffer factor [0-1]: x

Definition at line 74 of file scale_grid_cartesian.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume(), mod_cnvtopo::cnvtopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), scale_history::hist_switch(), and scale_land_grid::land_grid_setup().

74  real(RP), public, allocatable :: grid_cbfx(:)
real(rp), dimension(:), allocatable, public grid_cbfx
center buffer factor [0-1]: x

◆ grid_cbfy

real(rp), dimension(:), allocatable, public scale_grid::grid_cbfy

center buffer factor [0-1]: y

Definition at line 75 of file scale_grid_cartesian.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume(), mod_cnvtopo::cnvtopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), scale_history::hist_switch(), and scale_land_grid::land_grid_setup().

75  real(RP), public, allocatable :: grid_cbfy(:)
real(rp), dimension(:), allocatable, public grid_cbfy
center buffer factor [0-1]: y

◆ grid_fbfz

real(rp), dimension(:), allocatable, public scale_grid::grid_fbfz

face buffer factor [0-1]: z

Definition at line 76 of file scale_grid_cartesian.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

76  real(RP), public, allocatable :: grid_fbfz(:)
real(rp), dimension(:), allocatable, public grid_fbfz
face buffer factor [0-1]: z

◆ grid_fbfx

real(rp), dimension(:), allocatable, public scale_grid::grid_fbfx

face buffer factor [0-1]: x

Definition at line 77 of file scale_grid_cartesian.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

77  real(RP), public, allocatable :: grid_fbfx(:)
real(rp), dimension(:), allocatable, public grid_fbfx
face buffer factor [0-1]: x

◆ grid_fbfy

real(rp), dimension(:), allocatable, public scale_grid::grid_fbfy

face buffer factor [0-1]: y

Definition at line 78 of file scale_grid_cartesian.F90.

Referenced by scale_atmos_boundary::atmos_boundary_resume(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

78  real(RP), public, allocatable :: grid_fbfy(:)
real(rp), dimension(:), allocatable, public grid_fbfy
face buffer factor [0-1]: y

◆ grid_fxg

real(rp), dimension (:), allocatable, public scale_grid::grid_fxg

face coordinate [m]: x, global

Definition at line 80 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

80  real(RP), public, allocatable :: grid_fxg (:)
real(rp), dimension(:), allocatable, public grid_fxg
face coordinate [m]: x, global

◆ grid_fyg

real(rp), dimension (:), allocatable, public scale_grid::grid_fyg

face coordinate [m]: y, global

Definition at line 81 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

81  real(RP), public, allocatable :: grid_fyg (:)
real(rp), dimension(:), allocatable, public grid_fyg
face coordinate [m]: y, global

◆ grid_cxg

real(rp), dimension (:), allocatable, public scale_grid::grid_cxg

center coordinate [m]: x, global

Definition at line 82 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), scale_history::hist_switch(), and mod_mkinit::interporation_fact().

82  real(RP), public, allocatable :: grid_cxg (:)
real(rp), dimension(:), allocatable, public grid_cxg
center coordinate [m]: x, global

◆ grid_cyg

real(rp), dimension (:), allocatable, public scale_grid::grid_cyg

center coordinate [m]: y, global

Definition at line 83 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

83  real(RP), public, allocatable :: grid_cyg (:)
real(rp), dimension(:), allocatable, public grid_cyg
center coordinate [m]: y, global

◆ grid_fbfxg

real(rp), dimension(:), allocatable, public scale_grid::grid_fbfxg

face buffer factor [0-1]: x, global

Definition at line 84 of file scale_grid_cartesian.F90.

Referenced by scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

84  real(RP), public, allocatable :: grid_fbfxg(:)
real(rp), dimension(:), allocatable, public grid_fbfxg
face buffer factor [0-1]: x, global

◆ grid_fbfyg

real(rp), dimension(:), allocatable, public scale_grid::grid_fbfyg

face buffer factor [0-1]: y, global

Definition at line 85 of file scale_grid_cartesian.F90.

Referenced by scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

85  real(RP), public, allocatable :: grid_fbfyg(:)
real(rp), dimension(:), allocatable, public grid_fbfyg
face buffer factor [0-1]: y, global

◆ grid_cbfxg

real(rp), dimension(:), allocatable, public scale_grid::grid_cbfxg

center buffer factor [0-1]: x, global

Definition at line 86 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

86  real(RP), public, allocatable :: grid_cbfxg(:)
real(rp), dimension(:), allocatable, public grid_cbfxg
center buffer factor [0-1]: x, global

◆ grid_cbfyg

real(rp), dimension(:), allocatable, public scale_grid::grid_cbfyg

center buffer factor [0-1]: y, global

Definition at line 87 of file scale_grid_cartesian.F90.

Referenced by mod_copytopo::copytopo(), scale_fileio::fileio_set_axes(), grid_allocate(), grid_generate(), and scale_history::hist_switch().

87  real(RP), public, allocatable :: grid_cbfyg(:)
real(rp), dimension(:), allocatable, public grid_cbfyg
center buffer factor [0-1]: y, global