FE-Project
Loading...
Searching...
No Matches
scale_atm_dyn_dgm_nonhydro3d_common Module Reference

module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / Common More...

Functions/Subroutines

subroutine, public atm_dyn_dgm_nonhydro3d_common_init (mesh)
subroutine, public atm_dyn_dgm_nonhydro3d_common_final ()
subroutine, public atm_dyn_dgm_nonhydro3d_common_get_varinfo (prgvar_info, auxvar_info, phytend_info)
subroutine, public atm_dyn_dgm_nonhydro3d_common_setup_variables (prgvars, qtrcvars, auxvars, phytends, prgvar_manager, qtrcvar_manager, auxvar_manager, phytend_manager, phytend_num_tot, mesh3d, prgvar_varinfo)
subroutine, public atm_dyn_dgm_nonhydro3d_common_calc_pressure (pres, dpres, ddens, momx, momy, momz, therm, pres_hyd, dens_hyd, therm_hyd, rtot, cvtot, cptot, mesh3d, entot_conserve_scheme_flag)
subroutine, public atm_dyn_dgm_nonhydro3d_common_calc_therm_phyd (therm_hyd, pres_hyd, dens_hyd, mesh3d, entot_conserve_scheme_flag)
subroutine, public atm_dyn_dgm_nonhydro3d_common_drhot2pres (pres, dpres, drhot, pres_hyd, therm_hyd, rtot, cvtot, cptot, lcmesh, elem3d)
subroutine, public atm_dyn_dgm_nonhydro3d_common_drhot2entot (entot, ddens, momx, momy, momz, drhot, dens_hyd, pres_hyd, therm_hyd, rtot, cvtot, cptot, lcmesh, elem3d)
subroutine, public atm_dyn_dgm_nonhydro3d_common_entot2pres (pres, dpres, ddens, momx, momy, momz, entot, pres_hyd, dens_hyd, rtot, cvtot, lcmesh, elem3d)
subroutine, public atm_dyn_dgm_nonhydro3d_common_calc_rhot_hyd (rhot_hyd, pres_hyd, lcmesh, elem3d)
subroutine, public atm_dyn_dgm_nonhydro3d_common_calc_phyd_hgrad_lc (dphyddx, dphyddy, pres_hyd, pres_hyd_ref, element3d_operation, lmesh, elem)
 Calculate horizontal graidient of hydrostatic pressure In this calculation, we assume that PRES_hyd_ref is continuous at element boundaries.

Variables

integer, parameter, public prgvar_ddens_id = 1
integer, parameter, public prgvar_therm_id = 2
integer, parameter, public prgvar_drhot_id = 2
integer, parameter, public prgvar_etot_id = 2
integer, parameter, public prgvar_momz_id = 3
integer, parameter, public prgvar_momx_id = 4
integer, parameter, public prgvar_momy_id = 5
integer, parameter, public prgvar_scalar_num = 3
integer, parameter, public prgvar_hvec_num = 1
integer, parameter, public prgvar_num = 5
integer, parameter, public phytend_dens_id = 1
integer, parameter, public phytend_momx_id = 2
integer, parameter, public phytend_momy_id = 3
integer, parameter, public phytend_momz_id = 4
integer, parameter, public phytend_rhot_id = 5
integer, parameter, public phytend_rhoh_id = 6
integer, parameter, public phytend_num = 6
integer, parameter, public auxvar_preshydro_id = 1
integer, parameter, public auxvar_denshydro_id = 2
integer, parameter, public auxvar_thermhydro_id = 3
integer, parameter, public auxvar_pres_id = 4
integer, parameter, public auxvar_pt_id = 5
integer, parameter, public auxvar_rtot_id = 6
integer, parameter, public auxvar_cvtot_id = 7
integer, parameter, public auxvar_cptot_id = 8
integer, parameter, public auxvar_qdry_id = 9
integer, parameter, public auxvar_preshydro_ref_id = 10
integer, parameter, public auxvar_num = 10
real(rp), dimension(:,:), allocatable, public intrpmat_vpordm1

Detailed Description

module FElib / Fluid dyn solver / Atmosphere / Nonhydrostatic model / Common

Description
A common model for atmospheric nonhydrostatic dynamical core
Author
Yuta Kawai, Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
PRES_hydhydrostatic part of pressurePaPRES_hyd
DENS_hydhydrostatic part of densitykg/m3DENS_hyd
THERM_hydhydrostatic part of THERMkg/m3THERM_hyd
PRESpressurePaPRES
PTpotential temperatureKPT
RTOTTotal gas constantJ/kg/KRTOT
CVTOTTotal heat capacityJ/kg/KCVTOT
CPTOTTotal heat capacityJ/kg/KCPTOT
QDRYdry airkg/kgQDRY
PRES_hyd_REFhydrostatic reference pressurePaPRES_hyd_REF
DENS_tpDENS_tpkg/m3/sDENS_tp
MOMX_tpMOMX_tpkg/m2/sMOMX_tp
MOMY_tpMOMY_tpkg/m2/sMOMY_tp
MOMZ_tpMOMZ_tpkg/m2/sMOMZ_tp
RHOT_tpRHOT_tpkg/m3.K/sRHOT_tp
RHOH_pRHOH_pkg/m3.J/sRHOH_p
DDENSdeviation of densitykg/m3DDENS
THERMTHERM-THERM
MOMZmomentum zkg/m2/sMOMZ
MOMXmomentum xkg/m2/sMOMX
MOMYmomentum ykg/m2/sMOMY
{TRACER_NAME}{TRACER_NAME};
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNITTRACER_NAME
{TRACER_NAME}_tptendency of physical process for {TRACER_NAME};
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNITTRACER_NAME

Function/Subroutine Documentation

◆ atm_dyn_dgm_nonhydro3d_common_init()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_init ( class(meshbase3d), intent(in) mesh)

Definition at line 111 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

112
113 implicit none
114 class(MeshBase3D), intent(in) :: mesh
115
116 integer :: p1, p2, p_
117 type(ElementBase3D), pointer :: elem
118 real(RP) :: invV_VPOrdM1(mesh%refElem3D%Np,mesh%refElem3D%Np)
119
120 !--------------------------------------------
121
122 elem => mesh%refElem3D
123 allocate( intrpmat_vpordm1(elem%Np,elem%Np) )
124
125 invv_vpordm1(:,:) = elem%invV
126 do p2=1, elem%Nnode_h1D
127 do p1=1, elem%Nnode_h1D
128 p_ = p1 + (p2-1)*elem%Nnode_h1D + (elem%Nnode_v-1)*elem%Nnode_h1D**2
129 invv_vpordm1(p_,:) = 0.0_rp
130 end do
131 end do
132 intrpmat_vpordm1(:,:) = matmul(elem%V, invv_vpordm1)
133
134 return

References intrpmat_vpordm1.

Referenced by scale_atm_dyn_dgm_globalnonhydro3d_etot_heve::atm_dyn_dgm_globalnonhydro3d_etot_heve_init(), scale_atm_dyn_dgm_globalnonhydro3d_etot_hevi::atm_dyn_dgm_globalnonhydro3d_etot_hevi_init(), scale_atm_dyn_dgm_globalnonhydro3d_rhot_heve::atm_dyn_dgm_globalnonhydro3d_rhot_heve_init(), scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_rhot_hevi_init(), scale_atm_dyn_dgm_nonhydro3d_etot_heve::atm_dyn_dgm_nonhydro3d_etot_heve_init(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi::atm_dyn_dgm_nonhydro3d_etot_hevi_init(), scale_atm_dyn_dgm_nonhydro3d_rhot_heve::atm_dyn_dgm_nonhydro3d_rhot_heve_init(), scale_atm_dyn_dgm_nonhydro3d_rhot_heve_splitform::atm_dyn_dgm_nonhydro3d_rhot_heve_splitform_init(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_init(), and scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform::atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform_init().

◆ atm_dyn_dgm_nonhydro3d_common_final()

◆ atm_dyn_dgm_nonhydro3d_common_get_varinfo()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_get_varinfo ( type(variableinfo), dimension(prgvar_num), intent(out) prgvar_info,
type(variableinfo), dimension(auxvar_num), intent(out) auxvar_info,
type(variableinfo), dimension(phytend_num), intent(out), optional phytend_info )

Definition at line 147 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

149
150 implicit none
151 type(VariableInfo), intent(out) :: prgvar_info(PRGVAR_NUM)
152 type(VariableInfo), intent(out) :: auxvar_info(AUXVAR_NUM)
153 type(VariableInfo), intent(out), optional :: phytend_info(PHYTEND_NUM)
154
155 type(VariableInfo) :: PRGVAR_VARINFO(PRGVAR_NUM)
156 DATA prgvar_varinfo / &
157 variableinfo( prgvar_ddens_id, 'DDENS', 'deviation of density', &
158 'kg/m3', 3, 'XYZ', 'air_density' ), &
159 variableinfo( prgvar_therm_id, 'THERM', 'THERM', &
160 '-', 3, 'XYZ', '' ), &
161 variableinfo( prgvar_momz_id , 'MOMZ', 'momentum z', &
162 'kg/m2/s', 3, 'XYZ', 'northward_mass_flux_of_air' ), &
163 variableinfo( prgvar_momx_id , 'MOMX', 'momentum x', &
164 'kg/m2/s', 3, 'XYZ', 'upward_mass_flux_of_air' ), &
165 variableinfo( prgvar_momy_id , 'MOMY', 'momentum y', &
166 'kg/m2/s', 3, 'XYZ', 'eastward_mass_flux_of_air' ) /
167
168 type(VariableInfo) :: AUXVAR_VARINFO(AUXVAR_NUM)
169 DATA auxvar_varinfo / &
170 variableinfo( auxvar_preshydro_id, 'PRES_hyd', 'hydrostatic part of pressure', &
171 'Pa', 3, 'XYZ', '' ), &
172 variableinfo( auxvar_denshydro_id, 'DENS_hyd', 'hydrostatic part of density', &
173 'kg/m3', 3, 'XYZ', '' ), &
174 variableinfo( auxvar_thermhydro_id, 'THERM_hyd', 'hydrostatic part of THERM', &
175 'kg/m3', 3, 'XYZ', '' ), &
176 variableinfo( auxvar_pres_id , 'PRES', 'pressure', &
177 'Pa', 3, 'XYZ', 'air_pressure' ), &
178 variableinfo( auxvar_pt_id , 'PT', 'potential temperature', &
179 'K', 3, 'XYZ', 'potential_temperature' ), &
180 variableinfo( auxvar_rtot_id , 'RTOT', 'Total gas constant', &
181 'J/kg/K', 3, 'XYZ', '' ), &
182 variableinfo( auxvar_cvtot_id , 'CVTOT', 'Total heat capacity', &
183 'J/kg/K', 3, 'XYZ', '' ), &
184 variableinfo( auxvar_cptot_id , 'CPTOT', 'Total heat capacity', &
185 'J/kg/K', 3, 'XYZ', '' ), &
186 variableinfo( auxvar_qdry_id , 'QDRY', 'dry air', &
187 'kg/kg', 3, 'XYZ', '' ), &
188 variableinfo( auxvar_preshydro_ref_id, 'PRES_hyd_REF', 'hydrostatic reference pressure', &
189 'Pa', 3, 'XYZ', '' )/
190
191 type(VariableInfo) :: PHYTEND_VARINFO(PHYTEND_NUM)
192 DATA phytend_varinfo / &
193 variableinfo( phytend_dens_id, 'DENS_tp', 'DENS_tp', &
194 'kg/m3/s', 3, 'XYZ', 'tendency of physical process for DENS' ), &
195 variableinfo( phytend_momx_id, 'MOMX_tp', 'MOMX_tp', &
196 'kg/m2/s', 3, 'XYZ', 'tendency of physical process for MOMX' ), &
197 variableinfo( phytend_momy_id, 'MOMY_tp', 'MOMY_tp', &
198 'kg/m2/s', 3, 'XYZ', 'tendency of physical process for MOMY' ), &
199 variableinfo( phytend_momz_id, 'MOMZ_tp', 'MOMZ_tp', &
200 'kg/m2/s', 3, 'XYZ', 'tendency of physical process for MOMZ' ), &
201 variableinfo( phytend_rhot_id, 'RHOT_tp', 'RHOT_tp', &
202 'kg/m3.K/s', 3, 'XYZ', 'tendency of physical process for RHOT' ), &
203 variableinfo( phytend_rhoh_id, 'RHOH_p', 'RHOH_p', &
204 'kg/m3.J/s', 3, 'XYZ', 'heating of physical process for THERM' ) /
205
206 !----------------------------------------------------------
207
208 prgvar_info(:) = prgvar_varinfo
209 auxvar_info(:) = auxvar_varinfo
210 if ( present(phytend_info) ) phytend_info(:) = phytend_varinfo
211
212 return

References auxvar_cptot_id, auxvar_cvtot_id, auxvar_denshydro_id, auxvar_pres_id, auxvar_preshydro_id, auxvar_preshydro_ref_id, auxvar_pt_id, auxvar_qdry_id, auxvar_rtot_id, auxvar_thermhydro_id, phytend_dens_id, phytend_momx_id, phytend_momy_id, phytend_momz_id, phytend_rhoh_id, phytend_rhot_id, prgvar_ddens_id, prgvar_momx_id, prgvar_momy_id, prgvar_momz_id, and prgvar_therm_id.

Referenced by atm_dyn_dgm_nonhydro3d_common_setup_variables().

◆ atm_dyn_dgm_nonhydro3d_common_setup_variables()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_setup_variables ( type(meshfield3d), dimension(prgvar_num), intent(inout) prgvars,
type(meshfield3d), dimension(0:qa), intent(inout) qtrcvars,
type(meshfield3d), dimension(auxvar_num), intent(inout) auxvars,
type(meshfield3d), dimension(phytend_num_tot), intent(inout) phytends,
type(modelvarmanager), intent(inout) prgvar_manager,
type(modelvarmanager), intent(inout) qtrcvar_manager,
type(modelvarmanager), intent(inout) auxvar_manager,
type(modelvarmanager), intent(inout) phytend_manager,
integer, intent(in) phytend_num_tot,
class(meshbase3d), intent(in) mesh3d,
type(variableinfo), dimension(prgvar_num), intent(out) prgvar_varinfo )

Definition at line 215 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

220
221 use scale_atmos_hydrometeor, only: &
222 atmos_hydrometeor_dry
223 use scale_tracer, only: &
224 qa, tracer_name, tracer_desc, tracer_unit
227 implicit none
228 integer, intent(in) :: PHYTEND_NUM_TOT
229 type(MeshField3D), intent(inout) :: prgvars(PRGVAR_NUM)
230 type(MeshField3D), intent(inout) :: qtrcvars(0:QA)
231 type(MeshField3D), intent(inout) :: auxvars(AUXVAR_NUM)
232 type(MeshField3D), intent(inout) :: phytends(PHYTEND_NUM_TOT)
233 type(ModelVarManager), intent(inout) :: prgvar_manager
234 type(ModelVarManager), intent(inout) :: qtrcvar_manager
235 type(ModelVarManager), intent(inout) :: auxvar_manager
236 type(ModelVarManager), intent(inout) :: phytend_manager
237 class(MeshBase3D), intent(in) :: mesh3D
238
239 type(VariableInfo), intent(out) :: PRGVAR_VARINFO(PRGVAR_NUM)
240
241 type(VariableInfo) :: AUXVAR_VARINFO(AUXVAR_NUM)
242 type(VariableInfo) :: PHYTEND_VARINFO(PHYTEND_NUM)
243
244 integer :: iv
245 integer :: iq
246 logical :: reg_file_hist
247
248 type(VariableInfo) :: qtrc_dry_vinfo_tmp
249 type(VariableInfo) :: qtrc_dry_tp_vinfo_tmp
250 type(VariableInfo) :: qtrc_vinfo_tmp
251 type(VariableInfo) :: qtrc_tp_vinfo_tmp
252 !----------------------------------------------------------
253
254 call atm_dyn_dgm_nonhydro3d_common_get_varinfo( prgvar_varinfo, auxvar_varinfo, phytend_varinfo ) ! (out)
255
256 !- Initialize prognostic variables
257
258 reg_file_hist = .true.
259 do iv = 1, prgvar_num
260 call prgvar_manager%Regist( &
261 prgvar_varinfo(iv), mesh3d, & ! (in)
262 prgvars(iv), & ! (inout)
263 reg_file_hist, monitor_flag=.true., fill_zero=.true. ) ! (out)
264 end do
265
266 !- Initialize tracer variables
267
268 reg_file_hist = .true.
269
270 if ( qa == 0 ) then
271 ! Dummy
272 qtrc_dry_vinfo_tmp%ndims = 3
273 qtrc_dry_vinfo_tmp%dim_type = 'XYZ'
274 qtrc_dry_vinfo_tmp%STDNAME = ''
275
276 qtrc_dry_vinfo_tmp%keyID = 0
277 qtrc_dry_vinfo_tmp%NAME = "QV"
278 qtrc_dry_vinfo_tmp%DESC = "Ratio of Water Vapor mass to total mass (Specific humidity)"
279 qtrc_dry_vinfo_tmp%UNIT = "kg/kg"
280 call qtrcvar_manager%Regist( &
281 qtrc_dry_vinfo_tmp, mesh3d, & ! (in)
282 qtrcvars(0), & ! (inout)
283 .false., monitor_flag=.false., fill_zero=.true. ) ! (in)
284 else
285 qtrc_vinfo_tmp%ndims = 3
286 qtrc_vinfo_tmp%dim_type = 'XYZ'
287 qtrc_vinfo_tmp%STDNAME = ''
288
289 do iq = 1, qa
290 qtrc_vinfo_tmp%keyID = iq
291 qtrc_vinfo_tmp%NAME = tracer_name(iq)
292 qtrc_vinfo_tmp%DESC = tracer_desc(iq)
293 qtrc_vinfo_tmp%UNIT = tracer_unit(iq)
294 call qtrcvar_manager%Regist( &
295 qtrc_vinfo_tmp, mesh3d, & ! (in)
296 qtrcvars(iq), & ! (inout)
297 reg_file_hist, monitor_flag=.true., fill_zero=.true. ) ! (in)
298 end do
299 end if
300
301 !- Initialize auxiliary variables
302
303 reg_file_hist = .true.
304 do iv = 1, auxvar_num
305 call auxvar_manager%Regist( &
306 auxvar_varinfo(iv), mesh3d, & ! (in)
307 auxvars(iv), & ! (inout)
308 reg_file_hist, fill_zero=.true. ) ! (in)
309 end do
310
311 !- Initialize the tendency of physical processes
312
313 reg_file_hist = .true.
314 do iv = 1, phytend_num
315 call phytend_manager%Regist( &
316 phytend_varinfo(iv), mesh3d, & ! (in)
317 phytends(iv), & ! (inout)
318 reg_file_hist, fill_zero=.true. ) ! (in)
319 end do
320
321 if ( qa == 0 ) then
322 ! Dummy
323 qtrc_dry_tp_vinfo_tmp%ndims = 3
324 qtrc_dry_tp_vinfo_tmp%dim_type = 'XYZ'
325 qtrc_dry_tp_vinfo_tmp%STDNAME = ''
326
327 iv = phytend_num + 1
328 qtrc_dry_tp_vinfo_tmp%keyID = iv
329 qtrc_dry_tp_vinfo_tmp%NAME = "QV_tp"
330 qtrc_dry_tp_vinfo_tmp%DESC = "tendency of physical process for QV"
331 qtrc_dry_tp_vinfo_tmp%UNIT = "kg/m3/s"
332 call phytend_manager%Regist( &
333 qtrc_dry_tp_vinfo_tmp, mesh3d, & ! (in)
334 phytends(iv), & ! (inout)
335 .false., fill_zero=.true. ) ! (in)
336 else
337 qtrc_tp_vinfo_tmp%ndims = 3
338 qtrc_tp_vinfo_tmp%dim_type = 'XYZ'
339 qtrc_tp_vinfo_tmp%STDNAME = ''
340
341 do iq = 1, qa
342 iv = phytend_num + iq
343 qtrc_tp_vinfo_tmp%keyID = iv
344 qtrc_tp_vinfo_tmp%NAME = trim(tracer_name(iq))//'_tp'
345 qtrc_tp_vinfo_tmp%DESC = 'tendency of physical process for '//trim(tracer_desc(iq))
346 qtrc_tp_vinfo_tmp%UNIT = trim(tracer_unit(iq))//'/s'
347
348 call phytend_manager%Regist( &
349 qtrc_tp_vinfo_tmp, mesh3d, & ! (in)
350 phytends(iv), & ! (inout)
351 reg_file_hist, fill_zero=.true. ) ! (in)
352 end do
353 end if
354
355 return
module FElib / Mesh / Base 3D
module FElib / Data / base
Derived type representing a field with 3D mesh.

References atm_dyn_dgm_nonhydro3d_common_get_varinfo(), auxvar_num, phytend_num, and prgvar_num.

◆ atm_dyn_dgm_nonhydro3d_common_calc_pressure()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_calc_pressure ( class(meshfield3d), intent(inout) pres,
class(meshfield3d), intent(inout) dpres,
class(meshfield3d), intent(in) ddens,
class(meshfield3d), intent(in) momx,
class(meshfield3d), intent(in) momy,
class(meshfield3d), intent(in) momz,
class(meshfield3d), intent(in) therm,
class(meshfield3d), intent(in) pres_hyd,
class(meshfield3d), intent(in) dens_hyd,
class(meshfield3d), intent(in) therm_hyd,
class(meshfield3d), intent(in) rtot,
class(meshfield3d), intent(in) cvtot,
class(meshfield3d), intent(in) cptot,
class(meshbase3d), intent(in), target mesh3d,
logical, intent(in) entot_conserve_scheme_flag )

Definition at line 359 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

364
365 implicit none
366 class(MeshField3D), intent(inout) :: PRES
367 class(MeshField3D), intent(inout) :: DPRES
368 class(MeshField3D), intent(in) :: DDENS
369 class(MeshField3D), intent(in) :: MOMX
370 class(MeshField3D), intent(in) :: MOMY
371 class(MeshField3D), intent(in) :: MOMZ
372 class(MeshField3D), intent(in) :: THERM
373 class(MeshField3D), intent(in) :: PRES_hyd
374 class(MeshField3D), intent(in) :: DENS_hyd
375 class(MeshField3D), intent(in) :: THERM_hyd
376 class(MeshField3D), intent(in) :: Rtot
377 class(MeshField3D), intent(in) :: CVtot
378 class(MeshField3D), intent(in) :: CPtot
379 class(MeshBase3D), intent(in), target :: mesh3D
380 logical, intent(in) :: ENTOT_CONSERVE_SCHEME_FLAG
381
382 integer :: n
383 class(LocalMesh3D), pointer :: lcmesh3D
384 !---------------------------
385
386 do n=1, mesh3d%LOCAL_MESH_NUM
387 lcmesh3d => mesh3d%lcmesh_list(n)
388
389 if ( entot_conserve_scheme_flag ) then
390 call atm_dyn_dgm_nonhydro3d_common_entot2pres( pres%local(n)%val, dpres%local(n)%val, &
391 ddens%local(n)%val, momx%local(n)%val, momy%local(n)%val, momz%local(n)%val, therm%local(n)%val, &
392 pres_hyd%local(n)%val, dens_hyd%local(n)%val, rtot%local(n)%val, cvtot%local(n)%val, &
393 lcmesh3d, lcmesh3d%refElem3D )
394 else
395 call atm_dyn_dgm_nonhydro3d_common_drhot2pres( pres%local(n)%val, dpres%local(n)%val, &
396 therm%local(n)%val, pres_hyd%local(n)%val, therm_hyd%local(n)%val, rtot%local(n)%val, cvtot%local(n)%val, cptot%local(n)%val, &
397 lcmesh3d, lcmesh3d%refElem3D )
398 end if
399 end do
400
401 return

References atm_dyn_dgm_nonhydro3d_common_drhot2pres(), and atm_dyn_dgm_nonhydro3d_common_entot2pres().

◆ atm_dyn_dgm_nonhydro3d_common_calc_therm_phyd()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_calc_therm_phyd ( class(meshfield3d), intent(inout) therm_hyd,
class(meshfield3d), intent(in) pres_hyd,
class(meshfield3d), intent(in) dens_hyd,
class(meshbase3d), intent(in), target mesh3d,
logical, intent(in) entot_conserve_scheme_flag )

Definition at line 405 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

409
410 implicit none
411 class(MeshField3D), intent(inout) :: THERM_hyd
412 class(MeshField3D), intent(in) :: PRES_hyd
413 class(MeshField3D), intent(in) :: DENS_hyd
414 class(MeshBase3D), intent(in), target :: mesh3D
415 logical, intent(in) :: ENTOT_CONSERVE_SCHEME_FLAG
416
417 integer :: n
418 class(LocalMesh3D), pointer :: lcmesh3D
419 !---------------------------
420
421 do n=1, mesh3d%LOCAL_MESH_NUM
422 lcmesh3d => mesh3d%lcmesh_list(n)
423
424 if ( entot_conserve_scheme_flag ) then
425 ! ToDO
426 else
427 call atm_dyn_dgm_nonhydro3d_common_calc_rhot_hyd( therm_hyd%local(n)%val, &
428 pres_hyd%local(n)%val, lcmesh3d, lcmesh3d%refElem3D )
429 end if
430 end do
431
432 return

References atm_dyn_dgm_nonhydro3d_common_calc_rhot_hyd().

◆ atm_dyn_dgm_nonhydro3d_common_drhot2pres()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_drhot2pres ( real(rp), dimension(elem3d%np,lcmesh%nea), intent(out) pres,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(out) dpres,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) drhot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) therm_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) rtot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) cvtot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) cptot,
class(localmesh3d), intent(in) lcmesh,
class(elementbase3d), intent(in) elem3d )

Definition at line 436 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

439
440 use scale_const, only: &
441 rdry => const_rdry, &
442 cpdry => const_cpdry, &
443 cvdry => const_cvdry, &
444 pres00 => const_pre00
445
446 implicit none
447 class(LocalMesh3D), intent(in) :: lcmesh
448 class(ElementBase3D), intent(in) :: elem3D
449 real(RP), intent(out) :: PRES(elem3D%Np,lcmesh%NeA)
450 real(RP), intent(out) :: DPRES(elem3D%Np,lcmesh%NeA)
451 real(RP), intent(in) :: DRHOT(elem3D%Np,lcmesh%NeA)
452 real(RP), intent(in) :: PRES_hyd(elem3D%Np,lcmesh%NeA)
453 real(RP), intent(in) :: THERM_hyd(elem3D%Np,lcmesh%NeA)
454 real(RP), intent(in) :: Rtot(elem3D%Np,lcmesh%NeA)
455 real(RP), intent(in) :: CVtot(elem3D%Np,lcmesh%NeA)
456 real(RP), intent(in) :: CPtot(elem3D%Np,lcmesh%NeA)
457
458 integer :: ke
459 real(RP) :: RHOT(elem3D%Np)
460 real(RP) :: rP0
461 !---------------------------------------------------------------
462
463 rp0 = 1.0_rp / pres00
464 !$omp parallel do private( RHOT )
465 do ke=lcmesh%NeS, lcmesh%NeE
466! RHOT(:) = PRES00 / Rdry * ( PRES_hyd(:,ke) / PRES00 )**(CvDry/CpDry) + DRHOT(:,ke)
467 rhot(:) = therm_hyd(:,ke) + drhot(:,ke)
468
469 pres(:,ke) = pres00 * ( rtot(:,ke) * rp0 * rhot(:) )**( cptot(:,ke) / cvtot(:,ke) )
470 dpres(:,ke) = pres(:,ke) - pres_hyd(:,ke)
471 end do
472
473 return

Referenced by atm_dyn_dgm_nonhydro3d_common_calc_pressure(), and atm_dyn_dgm_nonhydro3d_common_drhot2entot().

◆ atm_dyn_dgm_nonhydro3d_common_drhot2entot()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_drhot2entot ( real(rp), dimension(elem3d%np,lcmesh%nea), intent(out) entot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) ddens,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) momx,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) momy,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) momz,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) drhot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) dens_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) therm_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) rtot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) cvtot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) cptot,
class(localmesh3d), intent(in) lcmesh,
class(elementbase3d), intent(in) elem3d )

Definition at line 477 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

481
482 use scale_const, only: &
483 grav => const_grav
484
485 implicit none
486 class(LocalMesh3D), intent(in) :: lcmesh
487 class(ElementBase3D), intent(in) :: elem3D
488 real(RP), intent(out) :: EnTot(elem3D%Np,lcmesh%NeA)
489 real(RP), intent(in) :: DDENS(elem3D%Np,lcmesh%NeA)
490 real(RP), intent(in) :: MOMX(elem3D%Np,lcmesh%NeA)
491 real(RP), intent(in) :: MOMY(elem3D%Np,lcmesh%NeA)
492 real(RP), intent(in) :: MOMZ(elem3D%Np,lcmesh%NeA)
493 real(RP), intent(in) :: DRHOT(elem3D%Np,lcmesh%NeA)
494 real(RP), intent(in) :: DENS_hyd(elem3D%Np,lcmesh%NeA)
495 real(RP), intent(in) :: PRES_hyd(elem3D%Np,lcmesh%NeA)
496 real(RP), intent(in) :: THERM_hyd(elem3D%Np,lcmesh%NeA)
497 real(RP), intent(in) :: Rtot(elem3D%Np,lcmesh%NeA)
498 real(RP), intent(in) :: CVtot(elem3D%Np,lcmesh%NeA)
499 real(RP), intent(in) :: CPtot(elem3D%Np,lcmesh%NeA)
500
501 integer :: ke, ke2D
502
503 real(RP) :: DENS(elem3D%Np)
504 real(RP) :: mom_u1(elem3D%Np), mom_u2(elem3D%Np)
505
506 real(RP) :: PRES(elem3D%Np,lcmesh%NeA)
507 real(RP) :: DPRES(elem3D%Np,lcmesh%NeA)
508 !---------------------------------------------------------------
509
510 call atm_dyn_dgm_nonhydro3d_common_drhot2pres( pres, dpres, &
511 drhot, pres_hyd, therm_hyd, rtot, cvtot, cptot, &
512 lcmesh, elem3d )
513
514 !$omp parallel do private( ke2D, DENS, mom_u1, mom_u2 )
515 do ke=lcmesh%NeS, lcmesh%NeE
516 ke2d = lcmesh%EMap3Dto2D(ke)
517
518 dens(:) = dens_hyd(:,ke) + ddens(:,ke)
519 mom_u1(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,1,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momy(:,ke)
520 mom_u2(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,2) * momy(:,ke)
521
522 entot(:,ke) = pres(:,ke) * cvtot(:,ke) / rtot(:,ke) &
523 + 0.5_rp * ( momx(:,ke) * mom_u1(:) + momy(:,ke) * mom_u2(:) + momz(:,ke)**2 ) / dens(:) &
524 + grav * dens(:) * lcmesh%zlev(:,ke)
525 end do
526
527 return

References atm_dyn_dgm_nonhydro3d_common_drhot2pres().

Referenced by mod_dg_prep::dg_prep().

◆ atm_dyn_dgm_nonhydro3d_common_entot2pres()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_entot2pres ( real(rp), dimension(elem3d%np,lcmesh%nea), intent(out) pres,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(out) dpres,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) ddens,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) momx,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) momy,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) momz,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) entot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) dens_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) rtot,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) cvtot,
class(localmesh3d), intent(in) lcmesh,
class(elementbase3d), intent(in) elem3d )

Definition at line 531 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

535
536 use scale_const, only: &
537 grav => const_grav
538
539 implicit none
540 class(LocalMesh3D), intent(in) :: lcmesh
541 class(ElementBase3D), intent(in) :: elem3D
542 real(RP), intent(out) :: PRES(elem3D%Np,lcmesh%NeA)
543 real(RP), intent(out) :: DPRES(elem3D%Np,lcmesh%NeA)
544 real(RP), intent(in) :: DDENS(elem3D%Np,lcmesh%NeA)
545 real(RP), intent(in) :: MOMX(elem3D%Np,lcmesh%NeA)
546 real(RP), intent(in) :: MOMY(elem3D%Np,lcmesh%NeA)
547 real(RP), intent(in) :: MOMZ(elem3D%Np,lcmesh%NeA)
548 real(RP), intent(in) :: EnTot(elem3D%Np,lcmesh%NeA)
549 real(RP), intent(in) :: PRES_hyd(elem3D%Np,lcmesh%NeA)
550 real(RP), intent(in) :: DENS_hyd(elem3D%Np,lcmesh%NeA)
551 real(RP), intent(in) :: Rtot(elem3D%Np,lcmesh%NeA)
552 real(RP), intent(in) :: CVtot(elem3D%Np,lcmesh%NeA)
553
554 integer :: ke, ke2D
555
556 real(RP) :: DENS(elem3D%Np)
557 real(RP) :: mom_u1(elem3D%Np), mom_u2(elem3D%Np)
558 !---------------------------------------------------------------
559
560 !$omp parallel do private( ke2D, DENS, mom_u1, mom_u2 )
561 do ke=lcmesh%NeS, lcmesh%NeE
562 ke2d = lcmesh%EMap3Dto2D(ke)
563
564 dens(:) = dens_hyd(:,ke) + ddens(:,ke)
565 mom_u1(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,1,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momy(:,ke)
566 mom_u2(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,2) * momy(:,ke)
567
568 pres(:,ke) = ( entot(:,ke) - grav * dens(:) * lcmesh%zlev(:,ke) &
569 - 0.5_rp * ( momx(:,ke) * mom_u1(:) + momy(:,ke) * mom_u2(:) + momz(:,ke)**2 ) / dens(:) &
570 ) * rtot(:,ke) / cvtot(:,ke)
571
572 dpres(:,ke) = pres(:,ke) - pres_hyd(:,ke)
573 end do
574
575 return

Referenced by atm_dyn_dgm_nonhydro3d_common_calc_pressure().

◆ atm_dyn_dgm_nonhydro3d_common_calc_rhot_hyd()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_calc_rhot_hyd ( real(rp), dimension(elem3d%np,lcmesh%nea), intent(out) rhot_hyd,
real(rp), dimension(elem3d%np,lcmesh%nea), intent(in) pres_hyd,
class(localmesh3d), intent(in) lcmesh,
class(elementbase3d), intent(in) elem3d )

Definition at line 579 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

582
583 use scale_const, only: &
584 rdry => const_rdry, &
585 cpdry => const_cpdry, &
586 cvdry => const_cvdry, &
587 pres00 => const_pre00
588
589 implicit none
590 class(LocalMesh3D), intent(in) :: lcmesh
591 class(ElementBase3D), intent(in) :: elem3D
592 real(RP), intent(out) :: RHOT_hyd(elem3D%Np,lcmesh%NeA)
593 real(RP), intent(in) :: PRES_hyd(elem3D%Np,lcmesh%NeA)
594
595 integer :: ke
596 real(RP) :: rP0
597 !---------------------------------------------------------------
598
599 rp0 = 1.0_rp / pres00
600
601 !$omp parallel do
602 do ke=lcmesh%NeS, lcmesh%NeE
603 rhot_hyd(:,ke) = pres00 / rdry * ( pres_hyd(:,ke) / pres00 )**(cvdry/cpdry)
604 end do
605
606 return

Referenced by atm_dyn_dgm_nonhydro3d_common_calc_therm_phyd(), and mod_dg_prep::dg_prep().

◆ atm_dyn_dgm_nonhydro3d_common_calc_phyd_hgrad_lc()

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::atm_dyn_dgm_nonhydro3d_common_calc_phyd_hgrad_lc ( real(rp), dimension(elem%np,lmesh%nea), intent(out) dphyddx,
real(rp), dimension(elem%np,lmesh%nea), intent(out) dphyddy,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd,
real(rp), dimension(elem%np,lmesh%nea), intent(in) pres_hyd_ref,
class(elementoperationbase3d), intent(in) element3d_operation,
class(localmesh3d), intent(in) lmesh,
class(elementbase3d), intent(in) elem )

Calculate horizontal graidient of hydrostatic pressure In this calculation, we assume that PRES_hyd_ref is continuous at element boundaries.

Definition at line 612 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

615 implicit none
616 class(LocalMesh3D), intent(in) :: lmesh
617 class(ElementBase3D), intent(in) :: elem
618 real(RP), intent(out) :: DPhydDx(elem%Np,lmesh%NeA)
619 real(RP), intent(out) :: DPhydDy(elem%Np,lmesh%NeA)
620 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
621 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
622 class(ElementOperationBase3D), intent(in) :: element3D_operation
623
624 integer :: ke, ke2D, p
625
626 real(RP) :: Flux(elem%Np,3), Fz(elem%Np), DFlux(elem%Np,4,2)
627 real(RP) :: del_flux_hyd(elem%NfpTot,2,lmesh%Ne)
628 real(RP) :: GsqrtV, RGsqrtV(elem%Np)
629
630 real(RP) :: E33
631 real(RP) :: GradPhyd_x, GradPhyd_y
632 !-----------------------------------------
633
634 call get_phyd_hgrad_numflux_generalhvc( del_flux_hyd, & ! (out)
635 pres_hyd, pres_hyd_ref, & ! (in)
636 lmesh%Gsqrt, lmesh%GsqrtH, lmesh%gam, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), & ! (in)
637 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
638 lmesh%vmapM, lmesh%vmapP, elem%IndexH2Dto3D_bnd, & ! (in)
639 lmesh, elem, lmesh%lcmesh2D, lmesh%lcmesh2D%refElem2D ) ! (in)
640
641 !$omp parallel do private( ke, ke2D, p, &
642 !$omp Flux, Fz, DFlux, GsqrtV, RGsqrtV, E33, &
643 !$omp GradPhyd_x, GradPhyd_y )
644 do ke = lmesh%NeS, lmesh%NeE
645 ke2d = lmesh%EMap3Dto2D(ke)
646
647 do p=1, elem%Np
648 gsqrtv = lmesh%Gsqrt(p,ke) / ( lmesh%gam(p,ke)**2 * lmesh%GsqrtH(elem%IndexH2Dto3D(p),ke2d) )
649
650 rgsqrtv(p) = 1.0_rp / gsqrtv
651 flux(p,1) = gsqrtv * ( pres_hyd(p,ke) - pres_hyd_ref(p,ke) )
652 end do
653
654 do p=1, elem%Np
655 flux(p,2) = flux(p,1)
656 flux(p,3) = lmesh%GI3(p,ke,1) * flux(p,1)
657 fz(p) = lmesh%GI3(p,ke,2) * flux(p,1)
658 end do
659
660 call element3d_operation%Div( flux, del_flux_hyd(:,1,ke), &
661 dflux(:,:,1) )
662 call element3d_operation%Dz( fz, dflux(:,3,2) )
663 call element3d_operation%Lift( del_flux_hyd(:,2,ke), dflux(:,4,2) )
664
665 do p=1, elem%Np
666 e33 = lmesh%Escale(p,ke,3,3)
667
668 gradphyd_x = lmesh%Escale(p,ke,1,1) * dflux(p,1,1) &
669 + e33 * dflux(p,3,1) &
670 + dflux(p,4,1)
671
672 gradphyd_y = lmesh%Escale(p,ke,2,2) * dflux(p,2,1) &
673 + e33 * dflux(p,3,2) &
674 + dflux(p,4,2)
675
676 dphyddx(p,ke) = gradphyd_x * rgsqrtv(p)
677 dphyddy(p,ke) = gradphyd_y * rgsqrtv(p)
678 end do
679 end do
680
681 return

Variable Documentation

◆ prgvar_ddens_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_ddens_id = 1

Definition at line 68 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

68 integer, public, parameter :: PRGVAR_DDENS_ID = 1

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), mod_atmos_dyn::atmosdyn_setup(), and mod_atmos_phy_tb::atmosphytb_setup().

◆ prgvar_therm_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_therm_id = 2

Definition at line 69 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

69 integer, public, parameter :: PRGVAR_THERM_ID = 2 ! Variable associated with energy equation (DRHOT or ETOT)

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ prgvar_drhot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_drhot_id = 2

Definition at line 70 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

70 integer, public, parameter :: PRGVAR_DRHOT_ID = 2

◆ prgvar_etot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_etot_id = 2

Definition at line 71 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

71 integer, public, parameter :: PRGVAR_ETOT_ID = 2

◆ prgvar_momz_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_momz_id = 3

Definition at line 72 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

72 integer, public, parameter :: PRGVAR_MOMZ_ID = 3

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ prgvar_momx_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_momx_id = 4

Definition at line 73 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

73 integer, public, parameter :: PRGVAR_MOMX_ID = 4

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ prgvar_momy_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_momy_id = 5

Definition at line 74 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

74 integer, public, parameter :: PRGVAR_MOMY_ID = 5

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ prgvar_scalar_num

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_scalar_num = 3

Definition at line 75 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

75 integer, public, parameter :: PRGVAR_SCALAR_NUM = 3

◆ prgvar_hvec_num

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_hvec_num = 1

Definition at line 76 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

76 integer, public, parameter :: PRGVAR_HVEC_NUM = 1

◆ prgvar_num

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_num = 5

Definition at line 77 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

77 integer, public, parameter :: PRGVAR_NUM = 5

Referenced by atm_dyn_dgm_nonhydro3d_common_setup_variables(), and mod_atmos_dyn::atmosdyn_setup().

◆ phytend_dens_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::phytend_dens_id = 1

Definition at line 79 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

79 integer, public, parameter :: PHYTEND_DENS_ID = 1

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_atmos_component::atmos_setup().

◆ phytend_momx_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::phytend_momx_id = 2

Definition at line 80 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

80 integer, public, parameter :: PHYTEND_MOMX_ID = 2

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_atmos_component::atmos_setup().

◆ phytend_momy_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::phytend_momy_id = 3

Definition at line 81 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

81 integer, public, parameter :: PHYTEND_MOMY_ID = 3

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_atmos_component::atmos_setup().

◆ phytend_momz_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::phytend_momz_id = 4

Definition at line 82 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

82 integer, public, parameter :: PHYTEND_MOMZ_ID = 4

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_atmos_component::atmos_setup().

◆ phytend_rhot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::phytend_rhot_id = 5

Definition at line 83 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

83 integer, public, parameter :: PHYTEND_RHOT_ID = 5

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_atmos_component::atmos_setup().

◆ phytend_rhoh_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::phytend_rhoh_id = 6

Definition at line 84 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

84 integer, public, parameter :: PHYTEND_RHOH_ID = 6

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_atmos_component::atmos_setup().

◆ phytend_num

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::phytend_num = 6

Definition at line 85 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

85 integer, public, parameter :: PHYTEND_NUM = 6

Referenced by atm_dyn_dgm_nonhydro3d_common_setup_variables(), and mod_atmos_component::atmos_setup().

◆ auxvar_preshydro_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_preshydro_id = 1

Definition at line 88 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

88 integer, public, parameter :: AUXVAR_PRESHYDRO_ID = 1

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_experiment::exp_setinitcond_lc::exp_setinitcond_lc().

◆ auxvar_denshydro_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_denshydro_id = 2

Definition at line 89 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

89 integer, public, parameter :: AUXVAR_DENSHYDRO_ID = 2

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo(), and mod_experiment::exp_setinitcond_lc::exp_setinitcond_lc().

◆ auxvar_thermhydro_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_thermhydro_id = 3

Definition at line 90 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

90 integer, public, parameter :: AUXVAR_THERMHYDRO_ID = 3

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_pres_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_pres_id = 4

Definition at line 91 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

91 integer, public, parameter :: AUXVAR_PRES_ID = 4

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_pt_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_pt_id = 5

Definition at line 92 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

92 integer, public, parameter :: AUXVAR_PT_ID = 5

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_rtot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_rtot_id = 6

Definition at line 93 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

93 integer, public, parameter :: AUXVAR_Rtot_ID = 6

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_cvtot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_cvtot_id = 7

Definition at line 94 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

94 integer, public, parameter :: AUXVAR_CVtot_ID = 7

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_cptot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_cptot_id = 8

Definition at line 95 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

95 integer, public, parameter :: AUXVAR_CPtot_ID = 8

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_qdry_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_qdry_id = 9

Definition at line 96 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

96 integer, public, parameter :: AUXVAR_Qdry_ID = 9

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_preshydro_ref_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_preshydro_ref_id = 10

Definition at line 97 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

97 integer, public, parameter :: AUXVAR_PRESHYDRO_REF_ID = 10

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_num

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_num = 10

Definition at line 98 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

98 integer, public, parameter :: AUXVAR_NUM = 10

Referenced by atm_dyn_dgm_nonhydro3d_common_setup_variables().

◆ intrpmat_vpordm1

real(rp), dimension(:,:), allocatable, public scale_atm_dyn_dgm_nonhydro3d_common::intrpmat_vpordm1

Definition at line 102 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

102 real(RP), public, allocatable :: IntrpMat_VPOrdM1(:,:)

Referenced by scale_atm_dyn_dgm_globalnonhydro3d_rhot_heve::atm_dyn_dgm_globalnh3d_rhot_heve_cal_tend_shallow_atm_asis(), scale_atm_dyn_dgm_globalnonhydro3d_etot_heve::atm_dyn_dgm_globalnonhydro3d_etot_heve_cal_tend(), scale_atm_dyn_dgm_globalnonhydro3d_etot_hevi::atm_dyn_dgm_globalnonhydro3d_etot_hevi_cal_vi(), scale_atm_dyn_dgm_globalnonhydro3d_rhot_heve::atm_dyn_dgm_globalnonhydro3d_rhot_heve_cal_tend_deep_atm(), scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_rhot_hevi_cal_vi(), scale_atm_dyn_dgm_globalnonhydro3d_rhot_hevi::atm_dyn_dgm_globalnonhydro3d_rhot_hevi_cal_vi_asis(), atm_dyn_dgm_nonhydro3d_common_final(), atm_dyn_dgm_nonhydro3d_common_init(), scale_atm_dyn_dgm_nonhydro3d_etot_heve::atm_dyn_dgm_nonhydro3d_etot_heve_cal_tend(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi::atm_dyn_dgm_nonhydro3d_etot_hevi_cal_vi(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_construct_matbnd(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_construct_matbnd_uv(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_eval_ax(), scale_atm_dyn_dgm_nonhydro3d_etot_hevi_common::atm_dyn_dgm_nonhydro3d_etot_hevi_common_eval_ax_uv(), scale_atm_dyn_dgm_nonhydro3d_rhot_heve::atm_dyn_dgm_nonhydro3d_rhot_heve_cal_tend_asis(), scale_atm_dyn_dgm_nonhydro3d_rhot_heve_splitform::atm_dyn_dgm_nonhydro3d_rhot_heve_splitform_cal_tend(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi::atm_dyn_dgm_nonhydro3d_rhot_hevi_cal_vi_asis(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_construct_matbnd_uv(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_eval_ax_uv(), scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_common_2::atm_dyn_dgm_nonhydro3d_rhot_hevi_common_solve(), and scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform::atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform_cal_vi().