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

module Atmosphere / Physics / Cloud Microphysics More...

Data Types

type  atmosphympvars
 Derived type to manage variables with cloud microphysics component. More...

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_RAINprecipitation flux (liquid) in MP processkg/m2/sMP_SFLX_RAIN
MP_SFLX_SNOWprecipitation flux (solid) in MP processkg/m2/sMP_SFLX_SNOW
MP_SFLX_ENGIinternal energy flux flux in MP processJ/m2/sMP_SFLX_ENGI
MP_DENS_ttendency of x-momentum in MP processkg/m3/sMP_DENS_t
MP_MOMX_ttendency of x-momentum in MP processkg/m2/s2MP_MOMX_t
MP_MOMY_ttendency of y-momentum in MP processkg/m2/s2MP_MOMY_t
MP_MOMZ_ttendency of z-momentum in MP processkg/m2/s2MP_MOMZ_t
MP_RHOT_ttendency of rho*PT in MP processkg/m3.K/sMP_RHOT_t
mp_RHOHdiabatic heating rate in MP processJ/kg/smp_RHOH
mp_EVAPORATEnumber concentration of evaporated cloudm-3mp_EVAPORATE
MP_{TRACER_NAME}_ttendency of rho*{TRACER_NAME} in MP process;
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNITTRACER_NAME
Vterm_{TRACER_NAME}terminal velocity of {TRACER_NAME};
{TRACER_NAME} depends on the physics schemes, e.g., QV, QC, QR.
TRACER_UNITTRACER_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 281 of file mod_atmos_phy_mp_vars.F90.

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

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 353 of file mod_atmos_phy_mp_vars.F90.

355
356 use scale_mesh_base, only: meshbase
358 implicit none
359
360 integer, intent(in) :: domID
361 class(MeshBase), intent(in) :: mesh
362 class(ModelVarManager), intent(inout) :: sfcflx_list
363 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_rain
364 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_snow
365 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_engi
366
367 class(MeshFieldBase), pointer :: field
368 !-------------------------------------------------------
369
370 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_rain_id, field)
371 call field%GetLocalMeshField(domid, sflx_rain)
372
373 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_snow_id, field)
374 call field%GetLocalMeshField(domid, sflx_snow)
375
376 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_engi_id, field)
377 call field%GetLocalMeshField(domid, sflx_engi)
378
379 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 83 of file mod_atmos_phy_mp_vars.F90.

83 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 84 of file mod_atmos_phy_mp_vars.F90.

84 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 85 of file mod_atmos_phy_mp_vars.F90.

85 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 86 of file mod_atmos_phy_mp_vars.F90.

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

87 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 88 of file mod_atmos_phy_mp_vars.F90.

88 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 89 of file mod_atmos_phy_mp_vars.F90.

89 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 90 of file mod_atmos_phy_mp_vars.F90.

90 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 92 of file mod_atmos_phy_mp_vars.F90.

92 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 110 of file mod_atmos_phy_mp_vars.F90.

110 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 111 of file mod_atmos_phy_mp_vars.F90.

111 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 112 of file mod_atmos_phy_mp_vars.F90.

112 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 113 of file mod_atmos_phy_mp_vars.F90.

113 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 115 of file mod_atmos_phy_mp_vars.F90.

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