FE-Project
Loading...
Searching...
No Matches
Data Types | Functions/Subroutines | Variables
mod_atmos_phy_mp_vars Module Reference

module ATMOSPHERE physics / Cloud Microphysics More...

Data Types

type  atmosphympvars
 

Functions/Subroutines

subroutine, public atmosphympvars_getlocalmeshfields_tend (domid, mesh, mp_tends_list, mp_dens_t, mp_momx_t, mp_momy_t, mp_momz_t, mp_rhot_t, mp_rhoh, mp_evap, mp_rhoq_t, lcmesh3d)
 
subroutine, public atmosphympvars_getlocalmeshfields_sfcflx (domid, mesh, sfcflx_list, sflx_rain, sflx_snow, sflx_engi)
 

Variables

integer, parameter, public atmos_phy_mp_dens_t_id = 1
 
integer, parameter, public atmos_phy_mp_momx_t_id = 2
 
integer, parameter, public atmos_phy_mp_momy_t_id = 3
 
integer, parameter, public atmos_phy_mp_momz_t_id = 4
 
integer, parameter, public atmos_phy_mp_rhot_t_id = 5
 
integer, parameter, public atmos_phy_mp_rhoh_id = 6
 
integer, parameter, public atmos_phy_mp_evaporate_id = 7
 
integer, parameter, public atmos_phy_mp_tends_num1 = 7
 
type(variableinfo), dimension(atmos_phy_mp_tends_num1), public atmos_phy_mp_tend_vinfo
 
integer, parameter, public atmos_phy_mp_aux2d_sflx_rain_id = 1
 
integer, parameter, public atmos_phy_mp_aux2d_sflx_snow_id = 2
 
integer, parameter, public atmos_phy_mp_aux2d_sflx_engi_id = 3
 
integer, parameter, public atmos_phy_mp_aux2d_num = 3
 
type(variableinfo), dimension(atmos_phy_mp_aux2d_num), public atmos_phy_mp_aux2d_vinfo
 

Detailed Description

module ATMOSPHERE physics / Cloud Microphysics

Description
Container for variables with cloud microphysics component
Author
Yuta Kawai, Team SCALE
NAMELIST
  • No namelist group
History Output
namedescriptionunitvariable
MP_SFLX_RAIN precipitation flux (liquid) in MP process kg/m2/s MP_SFLX_RAIN
MP_SFLX_SNOW precipitation flux (solid) in MP process kg/m2/s MP_SFLX_SNOW
MP_SFLX_ENGI internal energy flux flux in MP process J/m2/s MP_SFLX_ENGI
MP_DENS_t tendency of x-momentum in MP process kg/m3/s MP_DENS_t
MP_MOMX_t tendency of x-momentum in MP process kg/m2/s2 MP_MOMX_t
MP_MOMY_t tendency of y-momentum in MP process kg/m2/s2 MP_MOMY_t
MP_MOMZ_t tendency of z-momentum in MP process kg/m2/s2 MP_MOMZ_t
MP_RHOT_t tendency of rho*PT in MP process kg/m3.K/s MP_RHOT_t
mp_RHOH diabatic heating rate in MP process J/kg/s mp_RHOH
mp_EVAPORATE number concentration of evaporated cloud m-3 mp_EVAPORATE
MP_{TRACER_NAME}_t tendency of rho*{TRACER_NAME} in MP process;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNIT TRACER_NAME
Vterm_{TRACER_NAME} terminal velocity of {TRACER_NAME};
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNIT TRACER_NAME

Function/Subroutine Documentation

◆ atmosphympvars_getlocalmeshfields_tend()

subroutine, public mod_atmos_phy_mp_vars::atmosphympvars_getlocalmeshfields_tend ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) mp_tends_list,
class(localmeshfieldbase), intent(out), pointer mp_dens_t,
class(localmeshfieldbase), intent(out), pointer mp_momx_t,
class(localmeshfieldbase), intent(out), pointer mp_momy_t,
class(localmeshfieldbase), intent(out), pointer mp_momz_t,
class(localmeshfieldbase), intent(out), pointer mp_rhot_t,
class(localmeshfieldbase), intent(out), pointer mp_rhoh,
class(localmeshfieldbase), intent(out), pointer mp_evap,
type(localmeshfieldbaselist), dimension(:), intent(out) mp_rhoq_t,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 277 of file mod_atmos_phy_mp_vars.F90.

282
283 use scale_mesh_base, only: meshbase
285 implicit none
286
287 integer, intent(in) :: domID
288 class(MeshBase), intent(in) :: mesh
289 class(ModelVarManager), intent(inout) :: mp_tends_list
290 class(LocalMeshFieldBase), pointer, intent(out) :: mp_DENS_t
291 class(LocalMeshFieldBase), pointer, intent(out) :: mp_MOMX_t
292 class(LocalMeshFieldBase), pointer, intent(out) :: mp_MOMY_t
293 class(LocalMeshFieldBase), pointer, intent(out) :: mp_MOMZ_t
294 class(LocalMeshFieldBase), pointer, intent(out) :: mp_RHOT_t
295 class(LocalMeshFieldBase), pointer, intent(out) :: mp_RHOH
296 class(LocalMeshFieldBase), pointer, intent(out) :: mp_EVAP
297 type(LocalMeshFieldBaseList), intent(out) :: mp_RHOQ_t(:)
298 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
299
300 class(MeshFieldBase), pointer :: field
301 class(LocalMeshBase), pointer :: lcmesh
302
303 integer :: iq
304 !-------------------------------------------------------
305
306 !--
307 call mp_tends_list%Get(atmos_phy_mp_dens_t_id, field)
308 call field%GetLocalMeshField(domid, mp_dens_t)
309
310 call mp_tends_list%Get(atmos_phy_mp_momx_t_id, field)
311 call field%GetLocalMeshField(domid, mp_momx_t)
312
313 call mp_tends_list%Get(atmos_phy_mp_momy_t_id, field)
314 call field%GetLocalMeshField(domid, mp_momy_t)
315
316 call mp_tends_list%Get(atmos_phy_mp_momz_t_id, field)
317 call field%GetLocalMeshField(domid, mp_momz_t)
318
319 call mp_tends_list%Get(atmos_phy_mp_rhot_t_id, field)
320 call field%GetLocalMeshField(domid, mp_rhot_t)
321
322 call mp_tends_list%Get(atmos_phy_mp_rhoh_id, field)
323 call field%GetLocalMeshField(domid, mp_rhoh)
324
325 call mp_tends_list%Get(atmos_phy_mp_evaporate_id, field)
326 call field%GetLocalMeshField(domid, mp_evap)
327
328 !---
329 do iq = 1, size(mp_rhoq_t)
330 call mp_tends_list%Get(atmos_phy_mp_tends_num1 + iq, field)
331 call field%GetLocalMeshField(domid, mp_rhoq_t(iq)%ptr)
332 end do
333
334
335 if (present(lcmesh3d)) then
336 call mesh%GetLocalMesh( domid, lcmesh )
337 nullify( lcmesh3d )
338
339 select type(lcmesh)
340 type is (localmesh3d)
341 if (present(lcmesh3d)) lcmesh3d => lcmesh
342 end select
343 end if
344
345 return
module FElib / Mesh / Base
module FElib / Data / base

References atmos_phy_mp_dens_t_id, atmos_phy_mp_evaporate_id, atmos_phy_mp_momx_t_id, atmos_phy_mp_momy_t_id, atmos_phy_mp_momz_t_id, atmos_phy_mp_rhoh_id, atmos_phy_mp_rhot_t_id, and atmos_phy_mp_tends_num1.

Referenced by mod_atmos_phy_mp_vars::atmosphympvars::history().

◆ atmosphympvars_getlocalmeshfields_sfcflx()

subroutine, public mod_atmos_phy_mp_vars::atmosphympvars_getlocalmeshfields_sfcflx ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) sfcflx_list,
class(localmeshfieldbase), intent(out), pointer sflx_rain,
class(localmeshfieldbase), intent(out), pointer sflx_snow,
class(localmeshfieldbase), intent(out), pointer sflx_engi )

Definition at line 349 of file mod_atmos_phy_mp_vars.F90.

351
352 use scale_mesh_base, only: meshbase
354 implicit none
355
356 integer, intent(in) :: domID
357 class(MeshBase), intent(in) :: mesh
358 class(ModelVarManager), intent(inout) :: sfcflx_list
359 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_rain
360 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_snow
361 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_engi
362
363 class(MeshFieldBase), pointer :: field
364 !-------------------------------------------------------
365
366 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_rain_id, field)
367 call field%GetLocalMeshField(domid, sflx_rain)
368
369 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_snow_id, field)
370 call field%GetLocalMeshField(domid, sflx_snow)
371
372 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_engi_id, field)
373 call field%GetLocalMeshField(domid, sflx_engi)
374
375 return

References atmos_phy_mp_aux2d_num, atmos_phy_mp_aux2d_sflx_engi_id, atmos_phy_mp_aux2d_sflx_rain_id, and atmos_phy_mp_aux2d_sflx_snow_id.

Referenced by mod_atmos_component::atmos_setup(), mod_atmos_vars::atmosvars_getlocalmeshqtrcphytend(), and mod_atmos_phy_mp_vars::atmosphympvars::history().

Variable Documentation

◆ atmos_phy_mp_dens_t_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_dens_t_id = 1

Definition at line 82 of file mod_atmos_phy_mp_vars.F90.

82 integer, public, parameter :: ATMOS_PHY_MP_DENS_t_ID = 1

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_momx_t_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_momx_t_id = 2

Definition at line 83 of file mod_atmos_phy_mp_vars.F90.

83 integer, public, parameter :: ATMOS_PHY_MP_MOMX_t_ID = 2

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_momy_t_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_momy_t_id = 3

Definition at line 84 of file mod_atmos_phy_mp_vars.F90.

84 integer, public, parameter :: ATMOS_PHY_MP_MOMY_t_ID = 3

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_momz_t_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_momz_t_id = 4

Definition at line 85 of file mod_atmos_phy_mp_vars.F90.

85 integer, public, parameter :: ATMOS_PHY_MP_MOMZ_t_ID = 4

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_rhot_t_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_rhot_t_id = 5

Definition at line 86 of file mod_atmos_phy_mp_vars.F90.

86 integer, public, parameter :: ATMOS_PHY_MP_RHOT_t_ID = 5

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_rhoh_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_rhoh_id = 6

Definition at line 87 of file mod_atmos_phy_mp_vars.F90.

87 integer, public, parameter :: ATMOS_PHY_MP_RHOH_ID = 6

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_evaporate_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_evaporate_id = 7

Definition at line 88 of file mod_atmos_phy_mp_vars.F90.

88 integer, public, parameter :: ATMOS_PHY_MP_EVAPORATE_ID = 7

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_tends_num1

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_tends_num1 = 7

Definition at line 89 of file mod_atmos_phy_mp_vars.F90.

89 integer, public, parameter :: ATMOS_PHY_MP_TENDS_NUM1 = 7

Referenced by atmosphympvars_getlocalmeshfields_tend().

◆ atmos_phy_mp_tend_vinfo

type(variableinfo), dimension(atmos_phy_mp_tends_num1), public mod_atmos_phy_mp_vars::atmos_phy_mp_tend_vinfo

Definition at line 91 of file mod_atmos_phy_mp_vars.F90.

91 type(VariableInfo), public :: ATMOS_PHY_MP_TEND_VINFO(ATMOS_PHY_MP_TENDS_NUM1)

◆ atmos_phy_mp_aux2d_sflx_rain_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_aux2d_sflx_rain_id = 1

Definition at line 109 of file mod_atmos_phy_mp_vars.F90.

109 integer, public, parameter :: ATMOS_PHY_MP_AUX2D_SFLX_RAIN_ID = 1

Referenced by atmosphympvars_getlocalmeshfields_sfcflx().

◆ atmos_phy_mp_aux2d_sflx_snow_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_aux2d_sflx_snow_id = 2

Definition at line 110 of file mod_atmos_phy_mp_vars.F90.

110 integer, public, parameter :: ATMOS_PHY_MP_AUX2D_SFLX_SNOW_ID = 2

Referenced by atmosphympvars_getlocalmeshfields_sfcflx().

◆ atmos_phy_mp_aux2d_sflx_engi_id

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_aux2d_sflx_engi_id = 3

Definition at line 111 of file mod_atmos_phy_mp_vars.F90.

111 integer, public, parameter :: ATMOS_PHY_MP_AUX2D_SFLX_ENGI_ID = 3

Referenced by atmosphympvars_getlocalmeshfields_sfcflx().

◆ atmos_phy_mp_aux2d_num

integer, parameter, public mod_atmos_phy_mp_vars::atmos_phy_mp_aux2d_num = 3

Definition at line 112 of file mod_atmos_phy_mp_vars.F90.

112 integer, public, parameter :: ATMOS_PHY_MP_AUX2D_NUM = 3

Referenced by atmosphympvars_getlocalmeshfields_sfcflx().

◆ atmos_phy_mp_aux2d_vinfo

type(variableinfo), dimension(atmos_phy_mp_aux2d_num), public mod_atmos_phy_mp_vars::atmos_phy_mp_aux2d_vinfo

Definition at line 114 of file mod_atmos_phy_mp_vars.F90.

114 type(VariableInfo), public :: ATMOS_PHY_MP_AUX2D_VINFO(ATMOS_PHY_MP_AUX2D_NUM)