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

module ATMOSPHERE physics / sub-grid scale turbulence More...

Data Types

type  atmosphytbvars
 

Functions/Subroutines

subroutine atmosphytbvars_init (this, model_mesh)
 
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 )

Definition at line 91 of file mod_atmos_phy_tb_vars.F90.

92
93 use scale_tracer, only: &
94 tracer_name, tracer_desc, tracer_unit
97
98 implicit none
99 class(AtmosPhyTbVars), target, intent(inout) :: this
100 class(ModelMeshBase), target, intent(in) :: model_mesh
101
102 integer :: iv
103 integer :: iq
104 logical :: reg_file_hist
105
106 class(AtmosMesh), pointer :: atm_mesh
107 class(MeshBase2D), pointer :: mesh2D
108 class(MeshBase3D), pointer :: mesh3D
109 !--------------------------------------------------
110
111 log_info('AtmosPhyTbVars_Init',*)
112
113 this%TENDS_NUM_TOT = atmos_phy_tb_tends_num1 + qa
114
115 !- Initialize auxiliary and diagnostic variables
116
117 nullify( atm_mesh )
118 select type(model_mesh)
119 class is (atmosmesh)
120 atm_mesh => model_mesh
121 end select
122 mesh3d => atm_mesh%ptr_mesh
123
124 call mesh3d%GetMesh2D( mesh2d )
125
126 !- Get variable information
127 ! call atm_phy_tb_dgm_common_get_varinfo( auxvar_info, diagvar_info, tbtend_info ) ! (out)
128
129 !- Initialize objects to variables with turbulent model
130
131 call this%tends_manager%Init()
132 call this%auxvars_manager%Init()
133 call this%diagvars_manager%Init()
134
135 allocate( this%diagvars(atmos_phy_tb_diag_num) )
136 allocate( this%tends(this%TENDS_NUM_TOT) )
137 allocate( this%auxvars(atmos_phy_tb_aux_num) )
138
140 this%tends, this%auxvars, this%diagvars, & ! (inout)
141 this%tends_manager, this%auxvars_manager, this%diagvars_manager, & ! (inout)
142 this%TENDS_NUM_TOT, mesh3d ) ! (in)
143
144 !- Setup communication
145
146 call atm_mesh%Create_communicator( &
149 this%auxvars_manager, & ! (inout)
150 this%auxvars(:), & ! (in)
151 this%auxvars_commid ) ! (out)
152
153 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 177 of file mod_atmos_phy_tb_vars.F90.

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

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