FE-Project
Loading...
Searching...
No Matches
mod_atmos_phy_tb_vars.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> module Atmosphere / Physics / sub-grid scale turbulence
3!!
4!! @par Description
5!! Container for variables with sub-grid scale turbulence model
6!!
7!! @author Yuta Kawai, Team SCALE
8!!
9!<
10!-------------------------------------------------------------------------------
11#include "scaleFElib.h"
13 !-----------------------------------------------------------------------------
14 !
15 !++ Used modules
16 !
17 use scale_precision
18 use scale_io
19 use scale_prc
20 use scale_tracer, only: qa
21
23 use scale_mesh_base, only: meshbase
25 use scale_mesh_base3d, only: &
26 meshbase3d, &
27 dimtype_xyz => meshbase3d_dimtypeid_xyz
30 use scale_localmeshfield_base, only: &
32 use scale_meshfield_base, only: &
34
37
39
40 use scale_model_var_manager, only: &
41 modelvarmanager, variableinfo
42 use scale_model_mesh_manager, only: modelmeshbase
43
49 tb_momx_t_vid => atmos_phy_tb_momx_t_id,tb_momy_t_vid => atmos_phy_tb_momy_t_id, &
50 tb_momz_t_vid => atmos_phy_tb_momz_t_id, tb_rhot_t_vid => atmos_phy_tb_rhot_t_id
51
52 use mod_atmos_mesh, only: atmosmesh
53
54 !-----------------------------------------------------------------------------
55 implicit none
56 private
57
58
59 !-----------------------------------------------------------------------------
60 !
61 !++ Public type & procedures
62 !
63
64 !> Derived type to manage variables with atmospheric SGS turbulent component
65 type, public :: atmosphytbvars
66 type(meshfield3d), allocatable :: tends(:)
67 type(modelvarmanager) :: tends_manager
68
69 type(meshfield3d), allocatable :: auxvars(:)
70 type(modelvarmanager) :: auxvars_manager
71 integer :: auxvars_commid
72
73 type(meshfield3d), allocatable :: diagvars(:)
74 type(modelvarmanager) :: diagvars_manager
75
76 integer :: tends_num_tot
77 contains
78 procedure :: init => atmosphytbvars_init
79 procedure :: final => atmosphytbvars_final
80 procedure :: history => atmosphytbvars_history
81 end type atmosphytbvars
82
84
85 !-----------------------------------------------------------------------------
86 !
87 !++ Private procedures
88 !
89 !-------------------
90
91contains
92!> Setup an object to manage variables with atmospheric SGS turbulent component
93!!
94!! @param model_mesh Object to manage computational mesh of atmospheric model
95!OCL SERIAL
96 subroutine atmosphytbvars_init( this, model_mesh )
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
159 end subroutine atmosphytbvars_init
160
161!> Finalize an object to manage variables with atmospheric SGS turbulent component
162!OCL SERIAL
163 subroutine atmosphytbvars_final( this )
164 implicit none
165 class(atmosphytbvars), intent(inout) :: this
166 !--------------------------------------------------
167
168 log_info('AtmosPhyTbVars_Final',*)
169
170 call this%tends_manager%Final()
171 deallocate( this%tends )
172
173 call this%auxvars_manager%Final()
174 deallocate( this%auxvars )
175
176 call this%diagvars_manager%Final()
177 deallocate( this%diagvars )
178
179 return
180 end subroutine atmosphytbvars_final
181
182!OCL SERIAL
183 subroutine atmosphytbvars_getlocalmeshfields_tend( domID, mesh, tb_tends_list, &
184 tb_MOMX_t, tb_MOMY_t, tb_MOMZ_t, tb_RHOT_t, tb_RHOQ_t, &
185 lcmesh3D &
186 )
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
241
242!OCL SERIAL
243 subroutine atmosphytbvars_history( this )
245 implicit none
246 class(atmosphytbvars), intent(inout) :: this
247
248 integer :: v
249 integer :: hst_id
250 type(meshfield3d) :: tmp_field
251 class(meshbase3d), pointer :: mesh3d
252 !-------------------------------------------------------------------------
253
254 mesh3d => this%auxvars(1)%mesh
255
256 do v = 1, this%TENDS_NUM_TOT
257 hst_id = this%tends(v)%hist_id
258 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%tends(v) )
259 end do
260
261 do v = 1, atmos_phy_tb_aux_num
262 hst_id = this%auxvars(v)%hist_id
263 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%auxvars(v) )
264 end do
265
266 do v = 1, atmos_phy_tb_diag_num
267 hst_id = this%diagvars(v)%hist_id
268 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%diagvars(v) )
269 end do
270
271 return
272 end subroutine atmosphytbvars_history
273
274end module mod_atmos_phy_tb_vars
module Atmosphere / Mesh
module Atmosphere / Physics / sub-grid scale turbulence
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)
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_momy_t_id
integer, parameter, public atmos_phy_tb_aux_htensor_num
integer, parameter, public atmos_phy_tb_aux_num
integer, parameter, public atmos_phy_tb_momz_t_id
integer, parameter, public atmos_phy_tb_rhot_t_id
integer, parameter, public atmos_phy_tb_momx_t_id
module FElib / Element / Base
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
integer, public meshbase3d_dimtypeid_xyz
module FElib / Mesh / Base
module FElib / Data / base
module FElib / Data / Communication base
FElib / model framework / mesh manager.
FElib / model framework / variable manager.
Derived type to manage a computational mesh (base class)
Derived type to manage variables with atmospheric SGS turbulent component.
Derived type representing a 3D reference element.
Derived type to manage a local 3D computational domain.
Derived type to manage a local computational domain (base type)
Derived type representing a field with local mesh (base type)
Derived type representing a field with 2D mesh.
Derived type representing a field with 3D mesh.
Derived type representing a field (base type)
Container to save a pointer of MeshField(1D, 2D, 3D) object.