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

module Atmosphere / Physics / sub-grid scale turbulence More...

Data Types

type  atmosphytbvars
 Derived type to manage variables with atmospheric SGS turbulent component. More...

Functions/Subroutines

subroutine atmosphytbvars_init (this, model_mesh)
 Setup an object to manage variables with atmospheric SGS turbulent component.
subroutine, public atmosphytbvars_getlocalmeshfields_tend (domid, mesh, tb_tends_list, tb_momx_t, tb_momy_t, tb_momz_t, tb_rhot_t, tb_rhoq_t, lcmesh3d)

Detailed Description

module Atmosphere / Physics / sub-grid scale turbulence

Description
Container for variables with sub-grid scale turbulence model
Author
Yuta Kawai, Team SCALE

Function/Subroutine Documentation

◆ atmosphytbvars_init()

subroutine mod_atmos_phy_tb_vars::atmosphytbvars_init ( class(atmosphytbvars), intent(inout), target this,
class(modelmeshbase), intent(in), target model_mesh )

Setup an object to manage variables with atmospheric SGS turbulent component.

Parameters
model_meshObject to manage computational mesh of atmospheric model

Definition at line 96 of file mod_atmos_phy_tb_vars.F90.

97
98 use scale_tracer, only: &
99 tracer_name, tracer_desc, tracer_unit
100 use scale_atm_phy_tb_dgm_common, only: &
102
103 implicit none
104 class(AtmosPhyTbVars), target, intent(inout) :: this
105 class(ModelMeshBase), target, intent(in) :: model_mesh
106
107 integer :: iv
108 integer :: iq
109 logical :: reg_file_hist
110
111 class(AtmosMesh), pointer :: atm_mesh
112 class(MeshBase2D), pointer :: mesh2D
113 class(MeshBase3D), pointer :: mesh3D
114 !--------------------------------------------------
115
116 log_info('AtmosPhyTbVars_Init',*)
117
118 this%TENDS_NUM_TOT = atmos_phy_tb_tends_num1 + qa
119
120 !- Initialize auxiliary and diagnostic variables
121
122 nullify( atm_mesh )
123 select type(model_mesh)
124 class is (atmosmesh)
125 atm_mesh => model_mesh
126 end select
127 mesh3d => atm_mesh%ptr_mesh
128
129 call mesh3d%GetMesh2D( mesh2d )
130
131 !- Get variable information
132 ! call atm_phy_tb_dgm_common_get_varinfo( auxvar_info, diagvar_info, tbtend_info ) ! (out)
133
134 !- Initialize objects to variables with turbulent model
135
136 call this%tends_manager%Init()
137 call this%auxvars_manager%Init()
138 call this%diagvars_manager%Init()
139
140 allocate( this%diagvars(atmos_phy_tb_diag_num) )
141 allocate( this%tends(this%TENDS_NUM_TOT) )
142 allocate( this%auxvars(atmos_phy_tb_aux_num) )
143
145 this%tends, this%auxvars, this%diagvars, & ! (inout)
146 this%tends_manager, this%auxvars_manager, this%diagvars_manager, & ! (inout)
147 this%TENDS_NUM_TOT, mesh3d ) ! (in)
148
149 !- Setup communication
150
151 call atm_mesh%Create_communicator( &
154 this%auxvars_manager, & ! (inout)
155 this%auxvars(:), & ! (in)
156 this%auxvars_commid ) ! (out)
157
158 return
module FElib / Fluid dyn solver / Atmosphere / Physics turbulence / Common
integer, parameter, public atmos_phy_tb_aux_scalar_num
integer, parameter, public atmos_phy_tb_aux_hvec_num
integer, parameter, public atmos_phy_tb_diag_num
subroutine, public atm_phy_tb_dgm_common_setup_variables(tends, auxvars, diagvars, tends_manager, auxvar_manager, diagvar_manager, tends_num_tot, mesh3d)
integer, parameter, public atmos_phy_tb_tends_num1
integer, parameter, public atmos_phy_tb_aux_htensor_num
integer, parameter, public atmos_phy_tb_aux_num

References scale_atm_phy_tb_dgm_common::atm_phy_tb_dgm_common_setup_variables(), scale_atm_phy_tb_dgm_common::atmos_phy_tb_aux_htensor_num, scale_atm_phy_tb_dgm_common::atmos_phy_tb_aux_hvec_num, scale_atm_phy_tb_dgm_common::atmos_phy_tb_aux_num, scale_atm_phy_tb_dgm_common::atmos_phy_tb_aux_scalar_num, scale_atm_phy_tb_dgm_common::atmos_phy_tb_diag_num, and scale_atm_phy_tb_dgm_common::atmos_phy_tb_tends_num1.

◆ atmosphytbvars_getlocalmeshfields_tend()

subroutine, public mod_atmos_phy_tb_vars::atmosphytbvars_getlocalmeshfields_tend ( integer, intent(in) domid,
class(meshbase), intent(in) mesh,
class(modelvarmanager), intent(inout) tb_tends_list,
class(localmeshfieldbase), intent(out), pointer tb_momx_t,
class(localmeshfieldbase), intent(out), pointer tb_momy_t,
class(localmeshfieldbase), intent(out), pointer tb_momz_t,
class(localmeshfieldbase), intent(out), pointer tb_rhot_t,
type(localmeshfieldbaselist), dimension(:), intent(out), optional tb_rhoq_t,
class(localmesh3d), intent(out), optional, pointer lcmesh3d )

Definition at line 183 of file mod_atmos_phy_tb_vars.F90.

187
188 use scale_mesh_base, only: meshbase
190 implicit none
191
192 integer, intent(in) :: domID
193 class(MeshBase), intent(in) :: mesh
194 class(ModelVarManager), intent(inout) :: tb_tends_list
195 class(LocalMeshFieldBase), pointer, intent(out) :: tb_MOMX_t
196 class(LocalMeshFieldBase), pointer, intent(out) :: tb_MOMY_t
197 class(LocalMeshFieldBase), pointer, intent(out) :: tb_MOMZ_t
198 class(LocalMeshFieldBase), pointer, intent(out) :: tb_RHOT_t
199 type(LocalMeshFieldBaseList), intent(out), optional :: tb_RHOQ_t(:)
200 class(LocalMesh3D), pointer, intent(out), optional :: lcmesh3D
201
202 class(MeshFieldBase), pointer :: field
203 class(LocalMeshBase), pointer :: lcmesh
204
205 integer :: iq
206 !-------------------------------------------------------
207
208 !--
209 call tb_tends_list%Get(tb_momx_t_vid, field)
210 call field%GetLocalMeshField(domid, tb_momx_t)
211
212 call tb_tends_list%Get(tb_momy_t_vid, field)
213 call field%GetLocalMeshField(domid, tb_momy_t)
214
215 call tb_tends_list%Get(tb_momz_t_vid, field)
216 call field%GetLocalMeshField(domid, tb_momz_t)
217
218 call tb_tends_list%Get(tb_rhot_t_vid, field)
219 call field%GetLocalMeshField(domid, tb_rhot_t)
220
221 !---
222 if ( present(tb_rhoq_t) ) then
223 do iq = 1, size(tb_rhoq_t)
224 call tb_tends_list%Get(atmos_phy_tb_tends_num1 + iq, field)
225 call field%GetLocalMeshField(domid, tb_rhoq_t(iq)%ptr)
226 end do
227 end if
228
229 if (present(lcmesh3d)) then
230 call mesh%GetLocalMesh( domid, lcmesh )
231 nullify( lcmesh3d )
232
233 select type(lcmesh)
234 type is (localmesh3d)
235 if (present(lcmesh3d)) lcmesh3d => lcmesh
236 end select
237 end if
238
239 return
module FElib / Mesh / Base
module FElib / Data / base
Derived type representing a field (base type)

Referenced by mod_atmos_phy_tb::atmosphytb_setup(), and mod_atmos_phy_tb_vars::atmosphytbvars::history().