SCALE-RM
mod_atmos_dyn_driver.F90
Go to the documentation of this file.
1 !-------------------------------------------------------------------------------
10 #include "scalelib.h"
12  !-----------------------------------------------------------------------------
13  !
14  !++ used modules
15  !
16  use scale_precision
17  use scale_io
18  use scale_prof
20  use scale_index
21  use scale_tracer
22  !-----------------------------------------------------------------------------
23  implicit none
24  private
25  !-----------------------------------------------------------------------------
26  !
27  !++ Public procedure
28  !
29  public :: atmos_dyn_driver_setup
31  public :: atmos_dyn_driver
32 
33  !-----------------------------------------------------------------------------
34  !
35  !++ Public parameters & variables
36  !
37  character(len=H_SHORT), public :: atmos_dyn_tstep_large_type = 'FVM-HEVE'
38  character(len=H_SHORT), public :: atmos_dyn_tstep_tracer_type = 'FVM-HEVE'
39 
40  character(len=H_SHORT), public :: atmos_dyn_tinteg_large_type = 'EULER' ! Type of time integration
41  ! 'RK3'
42  character(len=H_SHORT), public :: atmos_dyn_tinteg_short_type = 'RK4'
43  ! 'RK3WS2002'
44  ! 'RK3'
45  ! 'RK7s6o' (only for FVM-HEVE)
46  character(len=H_SHORT), public :: atmos_dyn_tinteg_tracer_type = 'RK3WS2002'
47  ! 'EULER'
48 
49  character(len=H_SHORT), public :: atmos_dyn_fvm_flux_type = 'CD4' ! Type of advective flux scheme (FVM)
50  character(len=H_SHORT), public :: atmos_dyn_fvm_flux_tracer_type = 'UD3KOREN1993'
51  ! 'CD2'
52  ! 'UD3'
53  ! 'CD4'
54  ! 'UD5'
55  ! 'CD6'
56 
57  !-----------------------------------------------------------------------------
58  !
59  !++ Private procedure
60  !
61  !-----------------------------------------------------------------------------
62  !
63  !++ Private parameters & variables
64  !
65  ! Numerical filter
66  integer, private :: atmos_dyn_numerical_diff_laplacian_num
67  real(rp), private :: atmos_dyn_numerical_diff_coef ! nondimensional numerical diffusion
68  real(rp), private :: atmos_dyn_numerical_diff_coef_tracer ! nondimensional numerical diffusion for tracer
69  real(rp), private :: atmos_dyn_numerical_diff_sfc_fact
70  logical , private :: atmos_dyn_numerical_diff_use_refstate
71 
72  real(rp), private :: atmos_dyn_wdamp_tau ! maximum tau for Rayleigh damping of w [s]
73  real(rp), private :: atmos_dyn_wdamp_height ! height to start apply Rayleigh damping [m]
74  integer, private :: atmos_dyn_wdamp_layer ! layer number to start apply Rayleigh damping [num]
75 
76  ! Divergence damping
77  real(rp), private :: atmos_dyn_divdmp_coef ! Divergence dumping coef
78 
79  ! Flux-Corrected Transport limiter
80  logical, private :: atmos_dyn_flag_tracer_split_tend
81  logical, private :: atmos_dyn_flag_fct_momentum
82  logical, private :: atmos_dyn_flag_fct_t
83  logical, private :: atmos_dyn_flag_fct_tracer
84  logical, private :: atmos_dyn_flag_fct_along_stream
85 
86  !-----------------------------------------------------------------------------
87 contains
88  !-----------------------------------------------------------------------------
90  subroutine atmos_dyn_driver_setup
91  use scale_prc, only: &
92  prc_abort
93  use scale_atmos_grid_cartesc, only: &
94  domain_center_y => atmos_grid_cartesc_domain_center_y, &
95  cy => atmos_grid_cartesc_cy, &
96  fz => atmos_grid_cartesc_fz, &
97  cdz => atmos_grid_cartesc_cdz, &
98  cdx => atmos_grid_cartesc_cdx, &
99  cdy => atmos_grid_cartesc_cdy, &
100  fdz => atmos_grid_cartesc_fdz, &
101  fdx => atmos_grid_cartesc_fdx, &
103  use scale_atmos_grid_cartesc_real, only: &
104  real_lat => atmos_grid_cartesc_real_lat
105  use scale_time, only: &
107  use mod_atmos_admin, only: &
108  atmos_sw_dyn, &
110  use mod_atmos_vars, only: &
111  dens, &
112  momz, &
113  momx, &
114  momy, &
115  rhot, &
116  qtrc
117  use mod_atmos_dyn_vars, only: &
118  prog
119  use scale_atmos_dyn, only: &
121  implicit none
122 
123  namelist / param_atmos_dyn / &
129  atmos_dyn_numerical_diff_laplacian_num, &
130  atmos_dyn_numerical_diff_coef, &
131  atmos_dyn_numerical_diff_coef_tracer, &
132  atmos_dyn_numerical_diff_sfc_fact, &
133  atmos_dyn_numerical_diff_use_refstate, &
134  atmos_dyn_wdamp_tau, &
135  atmos_dyn_wdamp_height, &
136  atmos_dyn_wdamp_layer, &
137  atmos_dyn_divdmp_coef, &
138  atmos_dyn_flag_tracer_split_tend, &
139  atmos_dyn_flag_fct_momentum, &
140  atmos_dyn_flag_fct_t, &
141  atmos_dyn_flag_fct_tracer, &
142  atmos_dyn_flag_fct_along_stream
143 
144  real(rp) :: dt
145  integer :: ierr
146  !---------------------------------------------------------------------------
147 
148  log_newline
149  log_info("ATMOS_DYN_driver_setup",*) 'Setup'
150 
151  if ( atmos_sw_dyn ) then
152 
153  atmos_dyn_numerical_diff_laplacian_num = 2
154  atmos_dyn_numerical_diff_coef = 1.0e-4_rp
155  atmos_dyn_numerical_diff_coef_tracer = 0.0_rp
156  atmos_dyn_numerical_diff_sfc_fact = 1.0_rp
157  atmos_dyn_numerical_diff_use_refstate = .true.
158 
159  atmos_dyn_wdamp_tau = -1.0_rp
160  atmos_dyn_wdamp_height = -1.0_rp
161  atmos_dyn_wdamp_layer = -1
162 
163  atmos_dyn_divdmp_coef = 0.0_rp
164 
165  atmos_dyn_flag_tracer_split_tend = .false.
166  atmos_dyn_flag_fct_momentum = .false.
167  atmos_dyn_flag_fct_t = .false.
168  atmos_dyn_flag_fct_tracer = .false.
169  atmos_dyn_flag_fct_along_stream = .true.
170 
171  !--- read namelist
172  rewind(io_fid_conf)
173  read(io_fid_conf,nml=param_atmos_dyn,iostat=ierr)
174  if( ierr < 0 ) then !--- missing
175  log_info("ATMOS_DYN_driver_setup",*) 'Not found namelist. Default used.'
176  elseif( ierr > 0 ) then !--- fatal error
177  log_error("ATMOS_DYN_driver_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_DYN. Check!'
178  call prc_abort
179  endif
180  log_nml(param_atmos_dyn)
181 
182  dt = real(time_dtsec_atmos_dyn,kind=rp)
183 
184  if ( atmos_dyn_wdamp_layer > kmax ) then
185  log_error("ATMOS_DYN_driver_setup",*) 'ATMOS_DYN_wdamp_layer should be less than total number of vertical layer(KA). Check!'
186  call prc_abort
187  elseif( atmos_dyn_wdamp_layer > 0 ) then
188  atmos_dyn_wdamp_height = fz(atmos_dyn_wdamp_layer+ks-1)
189  endif
190 
191  if ( atmos_dyn_wdamp_tau < 0.0_rp ) then
192  atmos_dyn_wdamp_tau = dt * 10.0_rp
193  elseif ( atmos_dyn_wdamp_tau < dt ) then
194  log_error("ATMOS_DYN_driver_setup",*) 'ATMOS_DYN_wdamp_tau should be larger than TIME_DT_ATMOS_DYN. Check!'
195  call prc_abort
196  end if
197 
198  if ( atmos_sw_dyn ) then
199  log_newline
200  log_info("ATMOS_DYN_driver_setup",*) 'Scheme for Large time step : ', trim(atmos_dyn_tinteg_large_type)
201  log_info("ATMOS_DYN_driver_setup",*) 'Scheme for Short time step : ', trim(atmos_dyn_tinteg_short_type)
202  log_info("ATMOS_DYN_driver_setup",*) 'Scheme for Tracer advection : ', trim(atmos_dyn_tinteg_tracer_type)
203  endif
204 
210  atmos_dyn_type, & ! [IN] Note that ATMOS_DYN_TYPE corresponds to ATMOS_DYN_TSTEP_SHORT_TYPE.
211  atmos_dyn_fvm_flux_type, & ! [IN]
213  dens, momz, momx, momy, rhot, qtrc, & ! [IN]
214  prog, & ! [IN]
215  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
216  atmos_dyn_wdamp_tau, & ! [IN]
217  atmos_dyn_wdamp_height, & ! [IN]
218  fz, & ! [IN]
219  none = atmos_dyn_type=='NONE' ) ! [IN]
220  endif
221 
222  return
223  end subroutine atmos_dyn_driver_setup
224 
225  !-----------------------------------------------------------------------------
227  subroutine atmos_dyn_driver_finalize
228  use scale_atmos_dyn, only: &
230  use mod_atmos_admin, only: &
232 
233  if ( atmos_sw_dyn ) then
234  call atmos_dyn_finalize
235  end if
236 
237  return
238  end subroutine atmos_dyn_driver_finalize
239  !-----------------------------------------------------------------------------
241  subroutine atmos_dyn_driver( do_flag )
242  use scale_prc_cartesc, only: &
243  prc_twod
244  use scale_atmos_grid_cartesc, only: &
245  cdz => atmos_grid_cartesc_cdz, &
246  cdx => atmos_grid_cartesc_cdx, &
247  cdy => atmos_grid_cartesc_cdy, &
248  fdz => atmos_grid_cartesc_fdz, &
249  fdx => atmos_grid_cartesc_fdx, &
250  fdy => atmos_grid_cartesc_fdy, &
251  rcdz => atmos_grid_cartesc_rcdz, &
252  rcdx => atmos_grid_cartesc_rcdx, &
253  rcdy => atmos_grid_cartesc_rcdy, &
254  rfdz => atmos_grid_cartesc_rfdz, &
255  rfdx => atmos_grid_cartesc_rfdx, &
257  use scale_atmos_grid_cartesc_real, only: &
258  real_phi => atmos_grid_cartesc_real_phi
265  use scale_time, only: &
266  time_dtsec, &
268  use mod_atmos_admin, only: &
270  use scale_atmos_hydrometeor, only: &
271  i_qv
272  use mod_atmos_vars, only: &
274  dens, &
275  momz, &
276  momx, &
277  momy, &
278  rhot, &
279  qtrc, &
280  dens_av, &
281  momz_av, &
282  momx_av, &
283  momy_av, &
284  rhot_av, &
285  qtrc_av, &
286  dens_tp, &
287  rhou_tp, &
288  rhov_tp, &
289  rhot_tp, &
290  rhoq_tp, &
291  rhoh_p, &
292  momz_tp, &
293  momx_tp, &
294  momy_tp, &
295  cptot, &
296  exner
297  use mod_atmos_dyn_vars, only: &
298  prog
299  use scale_coriolis, only: &
300  coriolis_f
301  use scale_atmos_refstate, only: &
306  use mod_atmos_bnd_driver, only: &
321  bnd_qa, &
322  bnd_iq, &
324  use scale_atmos_dyn, only: &
325  atmos_dyn
326  use scale_comm_cartesc, only: &
327  comm_vars8, &
328  comm_wait
329  implicit none
330 
331  logical, intent(in) :: do_flag
332 
333  integer :: k, i, j, iq
334  !---------------------------------------------------------------------------
335 
336 #if defined DEBUG || defined QUICKDEBUG
337  rhot_tp( 1:ks-1,:,:) = 0.0_rp
338  rhot_tp(ke+1:ka, :,:) = 0.0_rp
339  momx_tp( 1:ks-1,:,:) = 0.0_rp
340  momx_tp(ke+1:ka, :,:) = 0.0_rp
341  momy_tp( 1:ks-1,:,:) = 0.0_rp
342  momy_tp(ke+1:ka, :,:) = 0.0_rp
343 #endif
344 
345  if ( do_flag ) then
346 
347  !$acc data &
348  !$acc copy(DENS,MOMZ,MOMX,MOMY,RHOT,QTRC) &
349  !$acc copyin(DENS_tp,RHOU_tp,RHOV_tp,RHOT_tp,RHOQ_tp,RHOH_p,MOMZ_tp,MOMX_tp,MOMY_tp, &
350  !$acc coriolis_f, &
351  !$acc atmos_boundary_dens,atmos_boundary_velz,atmos_boundary_velx,atmos_boundary_vely,atmos_boundary_pott,atmos_boundary_qtrc, &
352  !$acc atmos_boundary_alpha_dens,atmos_boundary_alpha_velz,atmos_boundary_alpha_velx,atmos_boundary_alpha_vely,atmos_boundary_alpha_pott,atmos_boundary_alpha_qtrc, &
353  !$acc atmos_boundary_mflux_offset_x,atmos_boundary_mflux_offset_y)
354 
355 
356  if ( .not. prc_twod ) &
357  call comm_vars8( rhou_tp, 1 )
358  call comm_vars8( rhov_tp, 2 )
359  call comm_vars8( momz_tp, 3 )
360  do iq = 1, qa
361  call comm_vars8( rhoq_tp(:,:,:,iq), 4+iq)
362  end do
363 
364  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
365  !$omp private(k,i,j) &
366  !$omp shared (KA,KS,KE,IS,IE,JS,JE, &
367  !$omp RHOT_tp,RHOH_p,CPtot,EXNER)
368  !$acc kernels
369  do j = js, je
370  do i = is, ie
371  do k = ks, ke
372  rhot_tp(k,i,j) = rhot_tp(k,i,j) &
373  + rhoh_p(k,i,j) / ( cptot(k,i,j) * exner(k,i,j) )
374  end do
375  end do
376  end do
377  !$acc end kernels
378  call comm_vars8( rhot_tp, 4 )
379 
380  if ( prc_twod ) then
381  !$omp parallel do default(none) OMP_SCHEDULE_ &
382  !$omp private(k,j) &
383  !$omp shared (KA,KS,KE,IS,JS,JE, &
384  !$omp MOMX_tp,RHOU_tp)
385  !$acc kernels
386  do j = js, je
387  do k = ks, ke
388  momx_tp(k,is,j) = momx_tp(k,is,j) + rhou_tp(k,is,j)
389  end do
390  end do
391  !$acc end kernels
392  else
393  call comm_wait ( rhou_tp, 1 )
394  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
395  !$omp private(k,i,j) &
396  !$omp shared (KA,KS,KE,IS,IE,JS,JE, &
397  !$omp MOMX_tp,RHOU_tp)
398  !$acc kernels
399  do j = js, je
400  do i = is, ie
401  do k = ks, ke
402  momx_tp(k,i,j) = momx_tp(k,i,j) &
403  + 0.5_rp * ( rhou_tp(k,i,j) + rhou_tp(k,i+1,j) )
404  end do
405  end do
406  end do
407  !$acc end kernels
408  call comm_vars8( momx_tp, 1 )
409  end if
410 
411  call comm_wait ( rhov_tp, 2 )
412  !$omp parallel do default(none) OMP_SCHEDULE_ collapse(2) &
413  !$omp private(k,i,j) &
414  !$omp shared (KA,KS,KE,IS,IE,JS,JE, &
415  !$omp MOMY_tp,RHOV_tp)
416  !$acc kernels
417  do j = js, je
418  do i = is, ie
419  do k = ks, ke
420  momy_tp(k,i,j) = momy_tp(k,i,j) &
421  + 0.5_rp * ( rhov_tp(k,i,j) + rhov_tp(k,i,j+1) )
422  end do
423  end do
424  end do
425  !$acc end kernels
426  call comm_vars8( momy_tp, 2 )
427 
428  call comm_wait ( momz_tp, 3 )
429  call comm_wait ( rhot_tp, 4, .false. )
430  do iq = 1, qa
431  call comm_wait ( rhoq_tp(:,:,:,iq), 4+iq, .false. )
432  end do
433  if ( .not. prc_twod ) &
434  call comm_wait ( momx_tp, 1 )
435  call comm_wait ( momy_tp, 2 )
436 
437 
438  call atmos_dyn( dens, momz, momx, momy, rhot, qtrc, & ! [INOUT]
439  prog, & ! [IN]
440  dens_av, momz_av, momx_av, momy_av, rhot_av, qtrc_av, & ! [INOUT]
442  coriolis_f, & ! [IN]
443  cdz, cdx, cdy, fdz, fdx, fdy, & ! [IN]
444  rcdz, rcdx, rcdy, rfdz, rfdx, rfdy, & ! [IN]
445  real_phi, & ! [IN]
446  gsqrt, j13g, j23g, j33g, mapf, & ! [IN]
448  atmos_refstate_dens, & ! [IN]
449  atmos_refstate_pott, & ! [IN]
450  atmos_refstate_qv, & ! [IN]
451  atmos_refstate_pres, & ! [IN]
452  atmos_dyn_numerical_diff_coef, & ! [IN]
453  atmos_dyn_numerical_diff_coef_tracer, & ! [IN]
454  atmos_dyn_numerical_diff_laplacian_num, & ! [IN]
455  atmos_dyn_numerical_diff_sfc_fact, & ! [IN]
456  atmos_dyn_numerical_diff_use_refstate, & ! [IN]
458  atmos_boundary_dens, & ! [IN]
459  atmos_boundary_velz, & ! [IN]
460  atmos_boundary_velx, & ! [IN]
461  atmos_boundary_vely, & ! [IN]
462  atmos_boundary_pott, & ! [IN]
463  atmos_boundary_qtrc, & ! [IN]
464  atmos_boundary_alpha_dens, & ! [IN]
465  atmos_boundary_alpha_velz, & ! [IN]
466  atmos_boundary_alpha_velx, & ! [IN]
467  atmos_boundary_alpha_vely, & ! [IN]
468  atmos_boundary_alpha_pott, & ! [IN]
469  atmos_boundary_alpha_qtrc, & ! [IN]
472  atmos_dyn_divdmp_coef, & ! [IN]
473  atmos_dyn_flag_tracer_split_tend, & ! [IN]
474  atmos_dyn_flag_fct_momentum, & ! [IN]
475  atmos_dyn_flag_fct_t, & ! [IN]
476  atmos_dyn_flag_fct_tracer, & ! [IN]
477  atmos_dyn_flag_fct_along_stream, & ! [IN]
478  atmos_use_average, & ! [IN]
479  i_qv, & ! [IN]
480  time_dtsec, & ! [IN]
481  time_dtsec_atmos_dyn ) ! [IN]
482 
483  call atmos_vars_check
484 
485  !$acc end data
486  endif
487 
488  return
489  end subroutine atmos_dyn_driver
490 
491 end module mod_atmos_dyn_driver
mod_atmos_vars::momz_av
real(rp), dimension(:,:,:), pointer, public momz_av
Definition: mod_atmos_vars.F90:91
scale_atmos_grid_cartesc_index::ke
integer, public ke
end point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:52
mod_atmos_bnd_driver::atmos_boundary_alpha_pott
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_pott
Definition: mod_atmos_bnd_driver.F90:56
scale_prc::prc_abort
subroutine, public prc_abort
Abort Process.
Definition: scale_prc.F90:350
mod_atmos_dyn_vars
module Atmosphere / Dynamics
Definition: mod_atmos_dyn_vars.F90:12
mod_atmos_dyn_driver::atmos_dyn_driver_finalize
subroutine, public atmos_dyn_driver_finalize
finalize
Definition: mod_atmos_dyn_driver.F90:228
mod_atmos_dyn_driver::atmos_dyn_tstep_large_type
character(len=h_short), public atmos_dyn_tstep_large_type
Definition: mod_atmos_dyn_driver.F90:37
scale_tracer::qa
integer, public qa
Definition: scale_tracer.F90:35
mod_atmos_vars::momx_av
real(rp), dimension(:,:,:), pointer, public momx_av
Definition: mod_atmos_vars.F90:92
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdz
z-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:43
mod_atmos_vars::rhoq_tp
real(rp), dimension(:,:,:,:), allocatable, public rhoq_tp
Definition: mod_atmos_vars.F90:121
scale_index
module Index
Definition: scale_index.F90:11
scale_tracer::tracer_mass
real(rp), dimension(qa_max), public tracer_mass
Definition: scale_tracer.F90:47
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdx
reciprocal of face-dx
Definition: scale_atmos_grid_cartesC.F90:68
mod_atmos_dyn_driver::atmos_dyn_tinteg_large_type
character(len=h_short), public atmos_dyn_tinteg_large_type
Definition: mod_atmos_dyn_driver.F90:40
mod_atmos_admin::atmos_use_average
logical, public atmos_use_average
Definition: mod_atmos_admin.F90:49
mod_atmos_vars::qtrc_av
real(rp), dimension(:,:,:,:), pointer, public qtrc_av
Definition: mod_atmos_vars.F90:95
scale_atmos_refstate::atmos_refstate_dens
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_dens
refernce density [kg/m3]
Definition: scale_atmos_refstate.F90:39
mod_atmos_dyn_driver::atmos_dyn_fvm_flux_tracer_type
character(len=h_short), public atmos_dyn_fvm_flux_tracer_type
Definition: mod_atmos_dyn_driver.F90:50
mod_atmos_dyn_vars::prog
real(rp), dimension(:,:,:,:), allocatable, public prog
Definition: mod_atmos_dyn_vars.F90:59
scale_precision
module PRECISION
Definition: scale_precision.F90:14
mod_atmos_vars::rhov_tp
real(rp), dimension(:,:,:), allocatable, public rhov_tp
Definition: mod_atmos_vars.F90:118
scale_atmos_grid_cartesc_index::ka
integer, public ka
Definition: scale_atmos_grid_cartesC_index.F90:47
scale_atmos_dyn
module Atmosphere / Dynamics FENT + FCT
Definition: scale_atmos_dyn.F90:13
mod_atmos_admin
module ATMOS admin
Definition: mod_atmos_admin.F90:11
scale_coriolis
module Coriolis
Definition: scale_coriolis.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdx
reciprocal of center-dx
Definition: scale_atmos_grid_cartesC.F90:66
mod_atmos_vars::cptot
real(rp), dimension(:,:,:), allocatable, target, public cptot
Definition: mod_atmos_vars.F90:143
mod_atmos_admin::atmos_sw_dyn
logical, public atmos_sw_dyn
Definition: mod_atmos_admin.F90:51
scale_atmos_refstate::atmos_refstate_pres
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_pres
refernce pressure [Pa]
Definition: scale_atmos_refstate.F90:37
scale_atmos_grid_cartesc_metric
module Atmosphere Grid CartesianC metirc
Definition: scale_atmos_grid_cartesC_metric.F90:12
mod_atmos_vars::rhot_av
real(rp), dimension(:,:,:), pointer, public rhot_av
Definition: mod_atmos_vars.F90:94
scale_atmos_hydrometeor
module atmosphere / hydrometeor
Definition: scale_atmos_hydrometeor.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdy
y-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:62
scale_atmos_refstate
module atmosphere / reference state
Definition: scale_atmos_refstate.F90:12
scale_atmos_refstate::atmos_refstate_pott
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_pott
refernce potential temperature [K]
Definition: scale_atmos_refstate.F90:40
mod_atmos_bnd_driver::atmos_boundary_dens
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_dens
Definition: mod_atmos_bnd_driver.F90:45
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdy
y-length of grid(j+1) to grid(j) [m]
Definition: scale_atmos_grid_cartesC.F90:64
mod_atmos_bnd_driver::atmos_boundary_alpha_dens
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_dens
Definition: mod_atmos_bnd_driver.F90:52
mod_atmos_bnd_driver::atmos_boundary_qtrc
real(rp), dimension(:,:,:,:), allocatable, public atmos_boundary_qtrc
Definition: mod_atmos_bnd_driver.F90:50
mod_atmos_bnd_driver::atmos_boundary_velx
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_velx
Definition: mod_atmos_bnd_driver.F90:47
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdz
reciprocal of center-dz
Definition: scale_atmos_grid_cartesC.F90:45
mod_atmos_vars::rhot
real(rp), dimension(:,:,:), allocatable, target, public rhot
Definition: mod_atmos_vars.F90:80
mod_atmos_vars::atmos_vars_check
subroutine, public atmos_vars_check(force)
Check variables for atmosphere.
Definition: mod_atmos_vars.F90:1475
scale_atmos_grid_cartesc_real
module Atmosphere GRID CartesC Real(real space)
Definition: scale_atmos_grid_cartesC_real.F90:11
mod_atmos_vars::qtrc
real(rp), dimension(:,:,:,:), allocatable, target, public qtrc
Definition: mod_atmos_vars.F90:81
scale_atmos_refstate::atmos_refstate_qv
real(rp), dimension(:,:,:), allocatable, public atmos_refstate_qv
refernce vapor [kg/kg]
Definition: scale_atmos_refstate.F90:41
scale_atmos_grid_cartesc::atmos_grid_cartesc_domain_center_y
real(rp), public atmos_grid_cartesc_domain_center_y
center position of global domain [m]: y
Definition: scale_atmos_grid_cartesC.F90:92
scale_prc
module PROCESS
Definition: scale_prc.F90:11
mod_atmos_bnd_driver::atmos_boundary_pott
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_pott
Definition: mod_atmos_bnd_driver.F90:49
scale_atmos_dyn::atmos_dyn_setup
subroutine, public atmos_dyn_setup(DYN_Tinteg_Short_TYPE, DYN_Tinteg_Tracer_TYPE, DYN_Tinteg_Large_TYPE, DYN_Tstep_Tracer_TYPE, DYN_Tstep_Large_TYPE, DYN_Tstep_Short_TYPE, DYN_FVM_FLUX_TYPE, DYN_FVM_FLUX_TYPE_TRACER, DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, CDZ, CDX, CDY, FDZ, FDX, FDY, wdamp_tau, wdamp_height, FZ, none)
Setup.
Definition: scale_atmos_dyn.F90:95
mod_atmos_vars::rhou_tp
real(rp), dimension(:,:,:), allocatable, public rhou_tp
Definition: mod_atmos_vars.F90:117
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_mapf
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_mapf
map factor
Definition: scale_atmos_grid_cartesC_metric.F90:35
scale_precision::rp
integer, parameter, public rp
Definition: scale_precision.F90:41
scale_atmos_grid_cartesc_index::ie
integer, public ie
end point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:54
scale_io
module STDIO
Definition: scale_io.F90:10
scale_atmos_grid_cartesc::atmos_grid_cartesc_rcdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rcdy
reciprocal of center-dy
Definition: scale_atmos_grid_cartesC.F90:67
mod_atmos_vars::dens
real(rp), dimension(:,:,:), allocatable, target, public dens
Definition: mod_atmos_vars.F90:76
mod_atmos_dyn_driver
module Atmosphere / Dynamics
Definition: mod_atmos_dyn_driver.F90:11
mod_atmos_bnd_driver::atmos_boundary_alpha_vely
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_vely
Definition: mod_atmos_bnd_driver.F90:55
scale_tracer::k
real(rp), public k
Definition: scale_tracer.F90:45
scale_atmos_grid_cartesc_index
module atmosphere / grid / cartesC index
Definition: scale_atmos_grid_cartesC_index.F90:12
mod_atmos_vars::momz
real(rp), dimension(:,:,:), allocatable, target, public momz
Definition: mod_atmos_vars.F90:77
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_gsqrt
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_gsqrt
transformation metrics from Z to Xi, {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:38
scale_atmos_grid_cartesc_index::kmax
integer, public kmax
Definition: scale_atmos_grid_cartesC_index.F90:36
mod_atmos_bnd_driver::bnd_iq
integer, dimension(:), allocatable, public bnd_iq
Definition: mod_atmos_bnd_driver.F90:43
mod_atmos_bnd_driver::bnd_qa
integer, public bnd_qa
Definition: mod_atmos_bnd_driver.F90:42
mod_atmos_vars::momy_av
real(rp), dimension(:,:,:), pointer, public momy_av
Definition: mod_atmos_vars.F90:93
scale_tracer::tracer_cv
real(rp), dimension(qa_max), public tracer_cv
Definition: scale_tracer.F90:42
scale_prc_cartesc
module process / cartesC
Definition: scale_prc_cartesC.F90:11
mod_atmos_bnd_driver
module ATMOSPHERE / Boundary treatment
Definition: mod_atmos_bnd_driver.F90:13
mod_atmos_bnd_driver::atmos_boundary_smoother_fact
real(rp), public atmos_boundary_smoother_fact
Definition: mod_atmos_bnd_driver.F90:62
scale_prof
module profiler
Definition: scale_prof.F90:11
mod_atmos_vars::momz_tp
real(rp), dimension(:,:,:), allocatable, public momz_tp
Definition: mod_atmos_vars.F90:116
mod_atmos_dyn_driver::atmos_dyn_tinteg_tracer_type
character(len=h_short), public atmos_dyn_tinteg_tracer_type
Definition: mod_atmos_dyn_driver.F90:46
mod_atmos_vars::momx
real(rp), dimension(:,:,:), allocatable, target, public momx
Definition: mod_atmos_vars.F90:78
mod_atmos_bnd_driver::atmos_boundary_alpha_velz
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_velz
Definition: mod_atmos_bnd_driver.F90:53
mod_atmos_vars::exner
real(rp), dimension(:,:,:), allocatable, target, public exner
Definition: mod_atmos_vars.F90:136
scale_atmos_grid_cartesc_index::is
integer, public is
start point of inner domain: x, local
Definition: scale_atmos_grid_cartesC_index.F90:53
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j23g
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j23g
(2,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:40
scale_time::time_dtsec_atmos_dyn
real(dp), public time_dtsec_atmos_dyn
time interval of dynamics [sec]
Definition: scale_time.F90:35
mod_atmos_vars::dens_tp
real(rp), dimension(:,:,:), allocatable, public dens_tp
Definition: mod_atmos_vars.F90:115
mod_atmos_vars::momy
real(rp), dimension(:,:,:), allocatable, target, public momy
Definition: mod_atmos_vars.F90:79
scale_time
module TIME
Definition: scale_time.F90:11
mod_atmos_dyn_driver::atmos_dyn_driver
subroutine, public atmos_dyn_driver(do_flag)
Dynamical Process (Wrapper)
Definition: mod_atmos_dyn_driver.F90:242
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdx
x-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:63
scale_tracer
module TRACER
Definition: scale_tracer.F90:12
mod_atmos_bnd_driver::atmos_boundary_alpha_qtrc
real(rp), dimension(:,:,:,:), allocatable, public atmos_boundary_alpha_qtrc
Definition: mod_atmos_bnd_driver.F90:57
mod_atmos_bnd_driver::atmos_boundary_vely
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_vely
Definition: mod_atmos_bnd_driver.F90:48
scale_atmos_hydrometeor::i_qv
integer, public i_qv
Definition: scale_atmos_hydrometeor.F90:93
mod_atmos_dyn_driver::atmos_dyn_driver_setup
subroutine, public atmos_dyn_driver_setup
Setup.
Definition: mod_atmos_dyn_driver.F90:91
mod_atmos_bnd_driver::atmos_boundary_velz
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_velz
Definition: mod_atmos_bnd_driver.F90:46
mod_atmos_vars::dens_av
real(rp), dimension(:,:,:), pointer, public dens_av
Definition: mod_atmos_vars.F90:90
scale_atmos_grid_cartesc_index::ks
integer, public ks
start point of inner domain: z, local
Definition: scale_atmos_grid_cartesC_index.F90:51
mod_atmos_dyn_driver::atmos_dyn_fvm_flux_type
character(len=h_short), public atmos_dyn_fvm_flux_type
Definition: mod_atmos_dyn_driver.F90:49
scale_tracer::tracer_cp
real(rp), dimension(qa_max), public tracer_cp
Definition: scale_tracer.F90:43
mod_atmos_vars::rhoh_p
real(rp), dimension(:,:,:), allocatable, public rhoh_p
Definition: mod_atmos_vars.F90:120
scale_atmos_grid_cartesc::atmos_grid_cartesc_fz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fz
face coordinate [m]: z, local
Definition: scale_atmos_grid_cartesC.F90:42
scale_coriolis::coriolis_f
real(rp), dimension(:,:), allocatable, public coriolis_f
Definition: scale_coriolis.F90:34
mod_atmos_admin::atmos_dyn_type
character(len=h_short), public atmos_dyn_type
Definition: mod_atmos_admin.F90:35
scale_time::time_dtsec
real(dp), public time_dtsec
time interval of model [sec]
Definition: scale_time.F90:33
mod_atmos_dyn_driver::atmos_dyn_tinteg_short_type
character(len=h_short), public atmos_dyn_tinteg_short_type
Definition: mod_atmos_dyn_driver.F90:42
scale_comm_cartesc
module COMMUNICATION
Definition: scale_comm_cartesC.F90:11
mod_atmos_vars
module ATMOSPHERIC Variables
Definition: mod_atmos_vars.F90:12
scale_tracer::tracer_r
real(rp), dimension(qa_max), public tracer_r
Definition: scale_tracer.F90:44
mod_atmos_vars::rhot_tp
real(rp), dimension(:,:,:), allocatable, public rhot_tp
Definition: mod_atmos_vars.F90:119
scale_atmos_grid_cartesc_index::js
integer, public js
start point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:55
scale_atmos_grid_cartesc::atmos_grid_cartesc_cy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cy
center coordinate [m]: y, local
Definition: scale_atmos_grid_cartesC.F90:57
mod_atmos_vars::momy_tp
real(rp), dimension(:,:,:), allocatable, public momy_tp
Definition: mod_atmos_vars.F90:125
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdz
reciprocal of face-dz
Definition: scale_atmos_grid_cartesC.F90:46
mod_atmos_vars::momx_tp
real(rp), dimension(:,:,:), allocatable, public momx_tp
Definition: mod_atmos_vars.F90:124
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j33g
real(rp), public atmos_grid_cartesc_metric_j33g
(3,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:41
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_phi
real(rp), dimension(:,:,:), allocatable, public atmos_grid_cartesc_real_phi
geopotential [m2/s2] (cell center)
Definition: scale_atmos_grid_cartesC_real.F90:64
mod_atmos_bnd_driver::atmos_boundary_mflux_offset_y
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_mflux_offset_y
Definition: mod_atmos_bnd_driver.F90:60
mod_atmos_bnd_driver::atmos_boundary_mflux_offset_x
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_mflux_offset_x
Definition: mod_atmos_bnd_driver.F90:59
scale_atmos_grid_cartesc::atmos_grid_cartesc_cdx
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_cdx
x-length of control volume [m]
Definition: scale_atmos_grid_cartesC.F90:61
scale_atmos_grid_cartesc_real::atmos_grid_cartesc_real_lat
real(rp), dimension(:,:), allocatable, public atmos_grid_cartesc_real_lat
latitude [rad,-pi,pi]
Definition: scale_atmos_grid_cartesC_real.F90:53
scale_atmos_grid_cartesc_metric::atmos_grid_cartesc_metric_j13g
real(rp), dimension(:,:,:,:), allocatable, public atmos_grid_cartesc_metric_j13g
(1,3) element of Jacobian matrix * {G}^1/2
Definition: scale_atmos_grid_cartesC_metric.F90:39
scale_atmos_grid_cartesc
module atmosphere / grid / cartesC
Definition: scale_atmos_grid_cartesC.F90:12
scale_atmos_grid_cartesc::atmos_grid_cartesc_fdz
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_fdz
z-length of grid(i+1) to grid(i) [m]
Definition: scale_atmos_grid_cartesC.F90:44
mod_atmos_bnd_driver::atmos_boundary_alpha_velx
real(rp), dimension(:,:,:), allocatable, public atmos_boundary_alpha_velx
Definition: mod_atmos_bnd_driver.F90:54
scale_atmos_grid_cartesc::atmos_grid_cartesc_rfdy
real(rp), dimension(:), allocatable, public atmos_grid_cartesc_rfdy
reciprocal of face-dy
Definition: scale_atmos_grid_cartesC.F90:69
scale_io::io_fid_conf
integer, public io_fid_conf
Config file ID.
Definition: scale_io.F90:57
scale_atmos_dyn::atmos_dyn
subroutine, public atmos_dyn(DENS, MOMZ, MOMX, MOMY, RHOT, QTRC, PROG, DENS_av, MOMZ_av, MOMX_av, MOMY_av, RHOT_av, QTRC_av, DENS_tp, MOMZ_tp, MOMX_tp, MOMY_tp, RHOT_tp, RHOQ_tp, CORIOLIS, CDZ, CDX, CDY, FDZ, FDX, FDY, RCDZ, RCDX, RCDY, RFDZ, RFDX, RFDY, PHI, GSQRT, J13G, J23G, J33G, MAPF, AQ_R, AQ_CV, AQ_CP, AQ_MASS, REF_dens, REF_pott, REF_qv, REF_pres, ND_COEF, ND_COEF_Q, ND_LAPLACIAN_NUM, ND_SFC_FACT, ND_USE_RS, BND_QA, BND_IQ, BND_SMOOTHER_FACT, DAMP_DENS, DAMP_VELZ, DAMP_VELX, DAMP_VELY, DAMP_POTT, DAMP_QTRC, DAMP_alpha_DENS, DAMP_alpha_VELZ, DAMP_alpha_VELX, DAMP_alpha_VELY, DAMP_alpha_POTT, DAMP_alpha_QTRC, MFLUX_OFFSET_X, MFLUX_OFFSET_Y, divdmp_coef, FLAG_TRACER_SPLIT_TEND, FLAG_FCT_MOMENTUM, FLAG_FCT_T, FLAG_FCT_TRACER, FLAG_FCT_ALONG_STREAM, USE_AVERAGE, I_QV, DTSEC, DTSEC_DYN)
Dynamical Process.
Definition: scale_atmos_dyn.F90:299
scale_atmos_grid_cartesc_index::je
integer, public je
end point of inner domain: y, local
Definition: scale_atmos_grid_cartesC_index.F90:56
scale_prc_cartesc::prc_twod
logical, public prc_twod
2D experiment
Definition: scale_prc_cartesC.F90:56
mod_atmos_dyn_driver::atmos_dyn_tstep_tracer_type
character(len=h_short), public atmos_dyn_tstep_tracer_type
Definition: mod_atmos_dyn_driver.F90:38
scale_atmos_dyn::atmos_dyn_finalize
subroutine, public atmos_dyn_finalize
finalize
Definition: scale_atmos_dyn.F90:236