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

module Atmosphere / Physics / surface process More...

Data Types

type  atmosphysfcvars
 Derived type to manage variables with a surface component. More...

Functions/Subroutines

subroutine, public atmosphysfcvars_getlocalmeshfields (domid, mesh, svars_list, sflx_list, sfc_temp, sflx_mu, sflx_mv, sflx_mw, sflx_sh, sflx_lh, sflx_qv, lcmesh3d)

Variables

integer, parameter, public atmos_phy_sf_svar_temp_id = 1
integer, parameter, public atmos_phy_sf_svar_num = 1
type(variableinfo), dimension(atmos_phy_sf_svar_num), public atmos_phy_sf_svar_vinfo
integer, parameter, public atmos_phy_sf_sflx_mu_id = 1
integer, parameter, public atmos_phy_sf_sflx_mv_id = 2
integer, parameter, public atmos_phy_sf_sflx_mw_id = 3
integer, parameter, public atmos_phy_sf_sflx_sh_id = 4
integer, parameter, public atmos_phy_sf_sflx_lh_id = 5
integer, parameter, public atmos_phy_sf_sflx_qv_id = 6
integer, parameter, public atmos_phy_sf_sflx_num1 = 6
type(variableinfo), dimension(atmos_phy_sf_sflx_num1), public atmos_phy_sf_sflx_vinfo

Detailed Description

module Atmosphere / Physics / surface process

Description
Container for variables with surface model
Author
Yuta Kawai, Team SCALE
NAMELIST
  • PARAM_ATMOS_PHY_SFC_VARS
    nametypedefault valuecomment
    ATMOS_PHY_SFC_DEFAULT_SFC_TEMPreal(RP)300.0_RP

History Output
namedescriptionunitvariable
SFLX_MUx-momentum fluxm/s*kg/m2/sSFLX_MU
SFLX_MVy-momentum fluxm/s*kg/m2/sSFLX_MV
SFLX_MWz-momentum fluxm/s*kg/m2/sSFLX_MW
SFLX_SHsensible heat fluxJ/m2/sSFLX_SH
SFLX_LHlatent heat fluxJ/m2/sSFLX_LH
SFLX_QVwater vapor fluxkg/m2/sSFLX_QV
SFC_TEMPsurface skin temperatureKSFC_TEMP

Function/Subroutine Documentation

◆ atmosphysfcvars_getlocalmeshfields()

subroutine, public mod_atmos_phy_sfc_vars::atmosphysfcvars_getlocalmeshfields ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) svars_list,
class(modelvarmanager), intent(inout) sflx_list,
class(localmeshfieldbase), intent(out), pointer sfc_temp,
class(localmeshfieldbase), intent(out), pointer sflx_mu,
class(localmeshfieldbase), intent(out), pointer sflx_mv,
class(localmeshfieldbase), intent(out), pointer sflx_mw,
class(localmeshfieldbase), intent(out), pointer sflx_sh,
class(localmeshfieldbase), intent(out), pointer sflx_lh,
class(localmeshfieldbase), intent(out), pointer sflx_qv,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 231 of file mod_atmos_phy_sfc_vars.F90.

235
236 use scale_mesh_base, only: meshbase
238 implicit none
239
240 integer, intent(in) :: domID
241 class(MeshBase), intent(in) :: mesh
242 class(ModelVarManager), intent(inout) :: svars_list
243 class(ModelVarManager), intent(inout) :: sflx_list
244 class(LocalMeshFieldBase), pointer, intent(out) :: SFC_TEMP
245 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_MU
246 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_MV
247 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_MW
248 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_SH
249 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_LH
250 class(LocalMeshFieldBase), pointer, intent(out) :: SFLX_QV
251 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
252
253 class(MeshFieldBase), pointer :: field
254 class(LocalMeshBase), pointer :: lcmesh
255 !-------------------------------------------------------
256
257 !--
258 call svars_list%Get(atmos_phy_sf_svar_temp_id, field)
259 call field%GetLocalMeshField(domid, sfc_temp)
260
261 call sflx_list%Get(atmos_phy_sf_sflx_mu_id, field)
262 call field%GetLocalMeshField(domid, sflx_mu)
263
264 call sflx_list%Get(atmos_phy_sf_sflx_mv_id, field)
265 call field%GetLocalMeshField(domid, sflx_mv)
266
267 call sflx_list%Get(atmos_phy_sf_sflx_mw_id, field)
268 call field%GetLocalMeshField(domid, sflx_mw)
269
270 call sflx_list%Get(atmos_phy_sf_sflx_sh_id, field)
271 call field%GetLocalMeshField(domid, sflx_sh)
272
273 call sflx_list%Get(atmos_phy_sf_sflx_lh_id, field)
274 call field%GetLocalMeshField(domid, sflx_lh)
275
276 call sflx_list%Get(atmos_phy_sf_sflx_qv_id, field)
277 call field%GetLocalMeshField(domid, sflx_qv)
278 !---
279
280 if (present(lcmesh3d)) then
281 call mesh%GetLocalMesh( domid, lcmesh )
282 nullify( lcmesh3d )
283
284 select type(lcmesh)
285 type is (localmesh3d)
286 if (present(lcmesh3d)) lcmesh3d => lcmesh
287 end select
288 end if
289
290 return
module FElib / Mesh / Base
module FElib / Data / base
Derived type representing a field (base type)

References atmos_phy_sf_sflx_lh_id, atmos_phy_sf_sflx_mu_id, atmos_phy_sf_sflx_mv_id, atmos_phy_sf_sflx_mw_id, atmos_phy_sf_sflx_qv_id, atmos_phy_sf_sflx_sh_id, and atmos_phy_sf_svar_temp_id.

Variable Documentation

◆ atmos_phy_sf_svar_temp_id

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_svar_temp_id = 1

Definition at line 68 of file mod_atmos_phy_sfc_vars.F90.

68 integer, public, parameter :: ATMOS_PHY_SF_SVAR_TEMP_ID = 1

Referenced by atmosphysfcvars_getlocalmeshfields().

◆ atmos_phy_sf_svar_num

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_svar_num = 1

Definition at line 69 of file mod_atmos_phy_sfc_vars.F90.

69 integer, public, parameter :: ATMOS_PHY_SF_SVAR_NUM = 1

◆ atmos_phy_sf_svar_vinfo

type(variableinfo), dimension(atmos_phy_sf_svar_num), public mod_atmos_phy_sfc_vars::atmos_phy_sf_svar_vinfo

Definition at line 70 of file mod_atmos_phy_sfc_vars.F90.

70 type(VariableInfo), public :: ATMOS_PHY_SF_SVAR_VINFO(ATMOS_PHY_SF_SVAR_NUM)

◆ atmos_phy_sf_sflx_mu_id

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_mu_id = 1

Definition at line 75 of file mod_atmos_phy_sfc_vars.F90.

75 integer, public, parameter :: ATMOS_PHY_SF_SFLX_MU_ID = 1

Referenced by atmosphysfcvars_getlocalmeshfields().

◆ atmos_phy_sf_sflx_mv_id

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_mv_id = 2

Definition at line 76 of file mod_atmos_phy_sfc_vars.F90.

76 integer, public, parameter :: ATMOS_PHY_SF_SFLX_MV_ID = 2

Referenced by atmosphysfcvars_getlocalmeshfields().

◆ atmos_phy_sf_sflx_mw_id

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_mw_id = 3

Definition at line 77 of file mod_atmos_phy_sfc_vars.F90.

77 integer, public, parameter :: ATMOS_PHY_SF_SFLX_MW_ID = 3

Referenced by atmosphysfcvars_getlocalmeshfields().

◆ atmos_phy_sf_sflx_sh_id

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_sh_id = 4

Definition at line 78 of file mod_atmos_phy_sfc_vars.F90.

78 integer, public, parameter :: ATMOS_PHY_SF_SFLX_SH_ID = 4

Referenced by atmosphysfcvars_getlocalmeshfields().

◆ atmos_phy_sf_sflx_lh_id

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_lh_id = 5

Definition at line 79 of file mod_atmos_phy_sfc_vars.F90.

79 integer, public, parameter :: ATMOS_PHY_SF_SFLX_LH_ID = 5

Referenced by atmosphysfcvars_getlocalmeshfields().

◆ atmos_phy_sf_sflx_qv_id

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_qv_id = 6

Definition at line 80 of file mod_atmos_phy_sfc_vars.F90.

80 integer, public, parameter :: ATMOS_PHY_SF_SFLX_QV_ID = 6

Referenced by atmosphysfcvars_getlocalmeshfields().

◆ atmos_phy_sf_sflx_num1

integer, parameter, public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_num1 = 6

Definition at line 81 of file mod_atmos_phy_sfc_vars.F90.

81 integer, public, parameter :: ATMOS_PHY_SF_SFLX_NUM1 = 6

◆ atmos_phy_sf_sflx_vinfo

type(variableinfo), dimension(atmos_phy_sf_sflx_num1), public mod_atmos_phy_sfc_vars::atmos_phy_sf_sflx_vinfo

Definition at line 83 of file mod_atmos_phy_sfc_vars.F90.

83 type(VariableInfo), public :: ATMOS_PHY_SF_SFLX_VINFO(ATMOS_PHY_SF_SFLX_NUM1)