FE-Project
Loading...
Searching...
No Matches
Functions/Subroutines | Variables
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, rtot, cvtot, cptot, mesh3d, entot_conserve_scheme_flag)
 
subroutine, public atm_dyn_dgm_nonhydro3d_common_drhot2pres (pres, dpres, drhot, pres_hyd, rtot, cvtot, cptot, lcmesh, elem3d)
 
subroutine, public atm_dyn_dgm_nonhydro3d_common_drhot2entot (entot, ddens, momx, momy, momz, drhot, dens_hyd, pres_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_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_pres_id = 3
 
integer, parameter, public auxvar_pt_id = 4
 
integer, parameter, public auxvar_rtot_id = 5
 
integer, parameter, public auxvar_cvtot_id = 6
 
integer, parameter, public auxvar_cptot_id = 7
 
integer, parameter, public auxvar_qdry_id = 8
 
integer, parameter, public auxvar_preshydro_ref_id = 9
 
integer, parameter, public auxvar_num = 9
 
real(rp), dimension(:,:), allocatable, public intrpmat_vpordm1
 

Detailed Description

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

Description
A coomon model for atmospheric nonhydrostatic dynamical core
Author
Yuta Kawai, Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
PRES_hyd hydrostatic part of pressure Pa PRES_hyd
DENS_hyd hydrostatic part of density kg/m3 DENS_hyd
PRES pressure Pa PRES
PT potential temperature K PT
RTOT Total gas constant J/kg/K RTOT
CVTOT Total heat capacity J/kg/K CVTOT
CPTOT Total heat capacity J/kg/K CPTOT
QDRY dry air kg/kg QDRY
PRES_hyd_REF hydrostatic reference pressure Pa PRES_hyd_REF
DENS_tp DENS_tp kg/m3/s DENS_tp
MOMX_tp MOMX_tp kg/m2/s MOMX_tp
MOMY_tp MOMY_tp kg/m2/s MOMY_tp
MOMZ_tp MOMZ_tp kg/m2/s MOMZ_tp
RHOT_tp RHOT_tp kg/m3.K/s RHOT_tp
RHOH_p RHOH_p kg/m3.J/s RHOH_p
DDENS deviation of density kg/m3 DDENS
THERM THERM - THERM
MOMZ momentum z kg/m2/s MOMZ
MOMX momentum x kg/m2/s MOMX
MOMY momentum y kg/m2/s MOMY
{TRACER_NAME} {TRACER_NAME};
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNIT TRACER_NAME
{TRACER_NAME}_tp tendency of physical process for {TRACER_NAME};
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNIT TRACER_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 108 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

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

subroutine, public scale_atm_dyn_dgm_nonhydro3d_common::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 144 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

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

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

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) 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 354 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

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

References atm_dyn_dgm_nonhydro3d_common_drhot2pres(), and atm_dyn_dgm_nonhydro3d_common_entot2pres().

◆ 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) 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 399 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

402
403 use scale_const, only: &
404 rdry => const_rdry, &
405 cpdry => const_cpdry, &
406 cvdry => const_cvdry, &
407 pres00 => const_pre00
408
409 implicit none
410 class(LocalMesh3D), intent(in) :: lcmesh
411 class(elementbase3D), intent(in) :: elem3D
412 real(RP), intent(out) :: PRES(elem3D%Np,lcmesh%NeA)
413 real(RP), intent(out) :: DPRES(elem3D%Np,lcmesh%NeA)
414 real(RP), intent(in) :: DRHOT(elem3D%Np,lcmesh%NeA)
415 real(RP), intent(in) :: PRES_hyd(elem3D%Np,lcmesh%NeA)
416 real(RP), intent(in) :: Rtot(elem3D%Np,lcmesh%NeA)
417 real(RP), intent(in) :: CVtot(elem3D%Np,lcmesh%NeA)
418 real(RP), intent(in) :: CPtot(elem3D%Np,lcmesh%NeA)
419
420 integer :: ke
421 real(RP) :: RHOT(elem3D%Np)
422 real(RP) :: rP0
423 !---------------------------------------------------------------
424
425 rp0 = 1.0_rp / pres00
426 !$omp parallel do private( RHOT )
427 do ke=lcmesh%NeS, lcmesh%NeE
428 rhot(:) = pres00 / rdry * ( pres_hyd(:,ke) / pres00 )**(cvdry/cpdry) + drhot(:,ke)
429
430 pres(:,ke) = pres00 * ( rtot(:,ke) * rp0 * rhot(:) )**( cptot(:,ke) / cvtot(:,ke) )
431 dpres(:,ke) = pres(:,ke) - pres_hyd(:,ke)
432 end do
433
434 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) 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 438 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

442
443 use scale_const, only: &
444 grav => const_grav
445
446 implicit none
447 class(LocalMesh3D), intent(in) :: lcmesh
448 class(elementbase3D), intent(in) :: elem3D
449 real(RP), intent(out) :: EnTot(elem3D%Np,lcmesh%NeA)
450 real(RP), intent(in) :: DDENS(elem3D%Np,lcmesh%NeA)
451 real(RP), intent(in) :: MOMX(elem3D%Np,lcmesh%NeA)
452 real(RP), intent(in) :: MOMY(elem3D%Np,lcmesh%NeA)
453 real(RP), intent(in) :: MOMZ(elem3D%Np,lcmesh%NeA)
454 real(RP), intent(in) :: DRHOT(elem3D%Np,lcmesh%NeA)
455 real(RP), intent(in) :: DENS_hyd(elem3D%Np,lcmesh%NeA)
456 real(RP), intent(in) :: PRES_hyd(elem3D%Np,lcmesh%NeA)
457 real(RP), intent(in) :: Rtot(elem3D%Np,lcmesh%NeA)
458 real(RP), intent(in) :: CVtot(elem3D%Np,lcmesh%NeA)
459 real(RP), intent(in) :: CPtot(elem3D%Np,lcmesh%NeA)
460
461 integer :: ke, ke2D
462
463 real(RP) :: DENS(elem3D%Np)
464 real(RP) :: mom_u1(elem3D%Np), mom_u2(elem3D%Np)
465
466 real(RP) :: PRES(elem3D%Np,lcmesh%NeA)
467 real(RP) :: DPRES(elem3D%Np,lcmesh%NeA)
468 !---------------------------------------------------------------
469
470 call atm_dyn_dgm_nonhydro3d_common_drhot2pres( pres, dpres, &
471 drhot, pres_hyd, rtot, cvtot, cptot, &
472 lcmesh, elem3d )
473
474 !$omp parallel do private( ke2D, DENS, mom_u1, mom_u2 )
475 do ke=lcmesh%NeS, lcmesh%NeE
476 ke2d = lcmesh%EMap3Dto2D(ke)
477
478 dens(:) = dens_hyd(:,ke) + ddens(:,ke)
479 mom_u1(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,1,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momy(:,ke)
480 mom_u2(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,2) * momy(:,ke)
481
482 entot(:,ke) = pres(:,ke) * cvtot(:,ke) / rtot(:,ke) &
483 + 0.5_rp * ( momx(:,ke) * mom_u1(:) + momy(:,ke) * mom_u2(:) + momz(:,ke)**2 ) / dens(:) &
484 + grav * dens(:) * lcmesh%zlev(:,ke)
485 end do
486
487 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 491 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

495
496 use scale_const, only: &
497 grav => const_grav
498
499 implicit none
500 class(LocalMesh3D), intent(in) :: lcmesh
501 class(elementbase3D), intent(in) :: elem3D
502 real(RP), intent(out) :: PRES(elem3D%Np,lcmesh%NeA)
503 real(RP), intent(out) :: DPRES(elem3D%Np,lcmesh%NeA)
504 real(RP), intent(in) :: DDENS(elem3D%Np,lcmesh%NeA)
505 real(RP), intent(in) :: MOMX(elem3D%Np,lcmesh%NeA)
506 real(RP), intent(in) :: MOMY(elem3D%Np,lcmesh%NeA)
507 real(RP), intent(in) :: MOMZ(elem3D%Np,lcmesh%NeA)
508 real(RP), intent(in) :: EnTot(elem3D%Np,lcmesh%NeA)
509 real(RP), intent(in) :: PRES_hyd(elem3D%Np,lcmesh%NeA)
510 real(RP), intent(in) :: DENS_hyd(elem3D%Np,lcmesh%NeA)
511 real(RP), intent(in) :: Rtot(elem3D%Np,lcmesh%NeA)
512 real(RP), intent(in) :: CVtot(elem3D%Np,lcmesh%NeA)
513
514 integer :: ke, ke2D
515
516 real(RP) :: DENS(elem3D%Np)
517 real(RP) :: mom_u1(elem3D%Np), mom_u2(elem3D%Np)
518 !---------------------------------------------------------------
519
520 !$omp parallel do private( ke2D, DENS, mom_u1, mom_u2 )
521 do ke=lcmesh%NeS, lcmesh%NeE
522 ke2d = lcmesh%EMap3Dto2D(ke)
523
524 dens(:) = dens_hyd(:,ke) + ddens(:,ke)
525 mom_u1(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,1,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momy(:,ke)
526 mom_u2(:) = lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,1) * momx(:,ke) + lcmesh%G_ij(lcmesh%refElem3D%IndexH2Dto3D,ke2d,2,2) * momy(:,ke)
527
528 pres(:,ke) = ( entot(:,ke) - grav * dens(:) * lcmesh%zlev(:,ke) &
529 - 0.5_rp * ( momx(:,ke) * mom_u1(:) + momy(:,ke) * mom_u2(:) + momz(:,ke)**2 ) / dens(:) &
530 ) * rtot(:,ke) / cvtot(:,ke)
531
532 dpres(:,ke) = pres(:,ke) - pres_hyd(:,ke)
533 end do
534
535 return

Referenced by atm_dyn_dgm_nonhydro3d_common_calc_pressure().

◆ 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 541 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

544 implicit none
545 class(LocalMesh3D), intent(in) :: lmesh
546 class(ElementBase3D), intent(in) :: elem
547 real(RP), intent(out) :: DPhydDx(elem%Np,lmesh%NeA)
548 real(RP), intent(out) :: DPhydDy(elem%Np,lmesh%NeA)
549 real(RP), intent(in) :: PRES_hyd(elem%Np,lmesh%NeA)
550 real(RP), intent(in) :: PRES_hyd_ref(elem%Np,lmesh%NeA)
551 class(ElementOperationBase3D), intent(in) :: element3D_operation
552
553 integer :: ke, ke2D, p
554
555 real(RP) :: Fx(elem%Np), Fy(elem%Np), Fz(elem%Np,2), DFlux(elem%Np,4,2)
556 real(RP) :: del_flux_hyd(elem%NfpTot,2,lmesh%Ne)
557 real(RP) :: GsqrtV(elem%Np), RGsqrtV(elem%Np)
558
559 real(RP) :: E33
560 real(RP) :: GradPhyd_x, GradPhyd_y
561 !-----------------------------------------
562
563 call get_phyd_hgrad_numflux_generalhvc( del_flux_hyd, & ! (out)
564 pres_hyd, pres_hyd_ref, & ! (in)
565 lmesh%Gsqrt, lmesh%GsqrtH, lmesh%gam, lmesh%GI3(:,:,1), lmesh%GI3(:,:,2), & ! (in)
566 lmesh%normal_fn(:,:,1), lmesh%normal_fn(:,:,2), lmesh%normal_fn(:,:,3), & ! (in)
567 lmesh%vmapM, lmesh%vmapP, elem%IndexH2Dto3D_bnd, & ! (in)
568 lmesh, elem, lmesh%lcmesh2D, lmesh%lcmesh2D%refElem2D ) ! (in)
569
570 !$omp parallel do private( ke, ke2D, p, &
571 !$omp Fx, Fy, Fz, DFlux, GsqrtV, RGsqrtV, E33, &
572 !$omp GradPhyd_x, GradPhyd_y )
573 do ke = lmesh%NeS, lmesh%NeE
574 ke2d = lmesh%EMap3Dto2D(ke)
575
576 do p=1, elem%Np
577 gsqrtv(p) = lmesh%Gsqrt(p,ke) / ( lmesh%gam(p,ke)**2 * lmesh%GsqrtH(elem%IndexH2Dto3D(p),ke2d) )
578 rgsqrtv(p) = 1.0_rp / gsqrtv(p)
579 end do
580
581 do p=1, elem%Np
582 fx(p) = gsqrtv(p) * ( pres_hyd(p,ke) - pres_hyd_ref(p,ke) )
583 fy(p) = fx(p)
584 fz(p,1) = lmesh%GI3(p,ke,1) * fx(p)
585 fz(p,2) = lmesh%GI3(p,ke,2) * fx(p)
586 end do
587
588 call element3d_operation%Div( fx, fy, fz(:,1), del_flux_hyd(:,1,ke), &
589 dflux(:,1,1), dflux(:,2,1), dflux(:,3,1), dflux(:,4,1) )
590 call element3d_operation%Dz( fz(:,2), dflux(:,3,2) )
591 call element3d_operation%Lift( del_flux_hyd(:,2,ke), dflux(:,4,2) )
592
593 do p=1, elem%Np
594 e33 = lmesh%Escale(p,ke,3,3)
595
596 gradphyd_x = lmesh%Escale(p,ke,1,1) * dflux(p,1,1) &
597 + e33 * dflux(p,3,1) &
598 + dflux(p,4,1)
599
600 gradphyd_y = lmesh%Escale(p,ke,2,2) * dflux(p,2,1) &
601 + e33 * dflux(p,3,2) &
602 + dflux(p,4,2)
603
604 dphyddx(p,ke) = gradphyd_x * rgsqrtv(p)
605 dphyddy(p,ke) = gradphyd_y * rgsqrtv(p)
606 end do
607 end do
608
609 return

Variable Documentation

◆ prgvar_ddens_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_ddens_id = 1

Definition at line 66 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

66 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 67 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

67 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 68 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

68 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 69 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

69 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 70 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

70 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 71 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

71 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 72 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

72 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 73 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

73 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 74 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

74 integer, public, parameter :: PRGVAR_HVEC_NUM = 1

◆ prgvar_num

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::prgvar_num = 5

Definition at line 75 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

75 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 77 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

77 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 78 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

78 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 79 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

79 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 80 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

80 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 81 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

81 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 82 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

82 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 83 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

83 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 86 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

86 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 87 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

87 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_pres_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_pres_id = 3

Definition at line 88 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

88 integer, public, parameter :: AUXVAR_PRES_ID = 3

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_pt_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_pt_id = 4

Definition at line 89 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

89 integer, public, parameter :: AUXVAR_PT_ID = 4

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_rtot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_rtot_id = 5

Definition at line 90 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

90 integer, public, parameter :: AUXVAR_Rtot_ID = 5

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_cvtot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_cvtot_id = 6

Definition at line 91 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

91 integer, public, parameter :: AUXVAR_CVtot_ID = 6

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_cptot_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_cptot_id = 7

Definition at line 92 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

92 integer, public, parameter :: AUXVAR_CPtot_ID = 7

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_qdry_id

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_qdry_id = 8

Definition at line 93 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

93 integer, public, parameter :: AUXVAR_Qdry_ID = 8

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 = 9

Definition at line 94 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

94 integer, public, parameter :: AUXVAR_PRESHYDRO_REF_ID = 9

Referenced by atm_dyn_dgm_nonhydro3d_common_get_varinfo().

◆ auxvar_num

integer, parameter, public scale_atm_dyn_dgm_nonhydro3d_common::auxvar_num = 9

Definition at line 95 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

95 integer, public, parameter :: AUXVAR_NUM = 9

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 99 of file scale_atm_dyn_dgm_nonhydro3d_common.F90.

99 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(), 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_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(), and scale_atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform::atm_dyn_dgm_nonhydro3d_rhot_hevi_splitform_cal_vi().