11#include "scaleFElib.h"
66 integer :: tends_num_tot
67 integer,
allocatable :: vterm_hist_id(:)
70 procedure :: init => atmosphympvars_init
71 procedure :: final => atmosphympvars_final
72 procedure :: history => atmosphympvars_history
95 'kg/m3/s', 3,
'XYZ',
'' ), &
97 'kg/m2/s2', 3,
'XYZ',
'' ), &
99 'kg/m2/s2', 3,
'XYZ',
'' ), &
101 'kg/m2/s2', 3,
'XYZ',
'' ), &
103 'kg/m3.K/s', 3,
'XYZ',
'' ), &
105 'J/kg/s', 3,
'XYZ',
'' ), &
107 'm-3' , 3,
'XYZ',
'' ) /
118 'kg/m2/s', 2,
'XY',
'' ), &
120 'kg/m2/s', 2,
'XY',
'' ), &
122 'J/m2/s', 2,
'XY',
'' ) /
135 subroutine atmosphympvars_init( this, model_mesh, &
136 QS_MP, QE_MP, QA_MP )
138 use scale_tracer,
only: &
139 tracer_name, tracer_desc, tracer_unit
140 use scale_file_history,
only: &
145 class(modelmeshbase),
target,
intent(in) :: model_mesh
146 integer,
intent(in) :: qs_mp
147 integer,
intent(in) :: qe_mp
148 integer,
intent(in) :: qa_mp
153 logical :: reg_file_hist
159 type(variableinfo) :: qtrc_tp_vinfo_tmp
160 type(variableinfo) :: qtrc_vterm_vinfo_tmp
163 log_info(
'AtmosPhyMpVars_Init',*)
173 select type(model_mesh)
175 atm_mesh => model_mesh
177 mesh3d => atm_mesh%ptr_mesh
179 call mesh3d%GetMesh2D( mesh2d )
183 call this%tends_manager%Init()
184 allocate( this%tends(this%TENDS_NUM_TOT) )
186 reg_file_hist = .true.
188 call this%tends_manager%Regist( &
190 this%tends(iv), reg_file_hist )
192 do n = 1, mesh3d%LOCAL_MESH_NUM
193 this%tends(iv)%local(n)%val(:,:) = 0.0_rp
197 qtrc_tp_vinfo_tmp%ndims = 3
198 qtrc_tp_vinfo_tmp%dim_type =
'XYZ'
199 qtrc_tp_vinfo_tmp%STDNAME =
''
203 qtrc_tp_vinfo_tmp%keyID = iv
204 qtrc_tp_vinfo_tmp%NAME =
'MP_'//trim(tracer_name(this%QS+iq-1))//
'_t'
205 qtrc_tp_vinfo_tmp%DESC =
'tendency of rho*'//trim(tracer_name(this%QS+iq-1))//
' in MP process'
206 qtrc_tp_vinfo_tmp%UNIT =
'kg/m3/s'
208 reg_file_hist = .true.
209 call this%tends_manager%Regist( &
210 qtrc_tp_vinfo_tmp, mesh3d, &
211 this%tends(iv), reg_file_hist )
213 do n = 1, mesh3d%LOCAL_MESH_NUM
214 this%tends(iv)%local(n)%val(:,:) = 0.0_rp
220 allocate( this%vterm_hist_id(this%QS+1:this%QE) )
221 allocate( this%vterm_hist (this%QS+1:this%QE) )
223 qtrc_vterm_vinfo_tmp%ndims = 3
224 qtrc_vterm_vinfo_tmp%dim_type =
'XYZ'
225 qtrc_vterm_vinfo_tmp%STDNAME =
''
227 do iq = this%QS+1, this%QE
228 qtrc_vterm_vinfo_tmp%NAME =
'Vterm_'//trim(tracer_name(this%QS+iq-1))
229 qtrc_vterm_vinfo_tmp%DESC =
'terminal velocity of '//trim(tracer_name(this%QS+iq-1))
230 qtrc_vterm_vinfo_tmp%UNIT =
'm/s'
231 call file_history_reg( qtrc_vterm_vinfo_tmp%NAME, qtrc_vterm_vinfo_tmp%DESC, qtrc_vterm_vinfo_tmp%UNIT, &
232 this%vterm_hist_id(iq), dim_type=
'XYZ' )
233 if ( this%vterm_hist_id(iq) > 0 )
call this%vterm_hist(iq)%Init( qtrc_vterm_vinfo_tmp%NAME, qtrc_vterm_vinfo_tmp%UNIT, mesh3d )
238 call this%auxvars2D_manager%Init()
241 reg_file_hist = .true.
243 call this%auxvars2D_manager%Regist( &
245 this%auxvars2D(iv), reg_file_hist )
247 do n = 1, mesh3d%LOCAL_MESH_NUM
248 this%auxvars2D(iv)%local(n)%val(:,:) = 0.0_rp
253 end subroutine atmosphympvars_init
257 subroutine atmosphympvars_final( this )
264 log_info(
'AtmosPhyMpVars_Final',*)
266 call this%tends_manager%Final()
267 deallocate( this%tends )
269 call this%auxvars2D_manager%Final()
270 deallocate( this%auxvars2D )
272 do iq = this%QS+1, this%QE
273 if ( this%vterm_hist_id(iq) > 0 )
call this%vterm_hist(iq)%Final()
275 deallocate( this%vterm_hist_id )
278 end subroutine atmosphympvars_final
282 mp_DENS_t, mp_MOMX_t, mp_MOMY_t, mp_MOMZ_t, mp_RHOT_t, mp_RHOH, mp_EVAP, &
291 integer,
intent(in) :: domid
302 class(
localmesh3d),
pointer,
intent(out),
optional :: lcmesh3d
312 call field%GetLocalMeshField(domid, mp_dens_t)
315 call field%GetLocalMeshField(domid, mp_momx_t)
318 call field%GetLocalMeshField(domid, mp_momy_t)
321 call field%GetLocalMeshField(domid, mp_momz_t)
324 call field%GetLocalMeshField(domid, mp_rhot_t)
327 call field%GetLocalMeshField(domid, mp_rhoh)
330 call field%GetLocalMeshField(domid, mp_evap)
333 do iq = 1,
size(mp_rhoq_t)
335 call field%GetLocalMeshField(domid, mp_rhoq_t(iq)%ptr)
339 if (
present(lcmesh3d))
then
340 call mesh%GetLocalMesh( domid, lcmesh )
345 if (
present(lcmesh3d)) lcmesh3d => lcmesh
354 SFLX_rain, SFLX_snow, SFLX_engi )
360 integer,
intent(in) :: domid
371 call field%GetLocalMeshField(domid, sflx_rain)
374 call field%GetLocalMeshField(domid, sflx_snow)
377 call field%GetLocalMeshField(domid, sflx_engi)
383 subroutine atmosphympvars_history( this )
391 type(meshfield3d) :: tmp_field
395 mesh3d => this%tends(1)%mesh
397 do v = 1, this%TENDS_NUM_TOT
398 hst_id = this%tends(v)%hist_id
403 hst_id = this%auxvars2D(v)%hist_id
407 do iq = this%QS+1, this%QE
408 hst_id = this%vterm_hist_id(iq)
413 end subroutine atmosphympvars_history
module Atmosphere / Physics / Cloud Microphysics
integer, parameter, public atmos_phy_mp_tends_num1
integer, parameter, public atmos_phy_mp_momx_t_id
integer, parameter, public atmos_phy_mp_momz_t_id
integer, parameter, public atmos_phy_mp_evaporate_id
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)
integer, parameter, public atmos_phy_mp_aux2d_sflx_snow_id
integer, parameter, public atmos_phy_mp_aux2d_sflx_engi_id
subroutine, public atmosphympvars_getlocalmeshfields_sfcflx(domid, mesh, sfcflx_list, sflx_rain, sflx_snow, sflx_engi)
type(variableinfo), dimension(atmos_phy_mp_aux2d_num), public atmos_phy_mp_aux2d_vinfo
integer, parameter, public atmos_phy_mp_aux2d_num
integer, parameter, public atmos_phy_mp_momy_t_id
integer, parameter, public atmos_phy_mp_aux2d_sflx_rain_id
type(variableinfo), dimension(atmos_phy_mp_tends_num1), public atmos_phy_mp_tend_vinfo
integer, parameter, public atmos_phy_mp_dens_t_id
integer, parameter, public atmos_phy_mp_rhot_t_id
integer, parameter, public atmos_phy_mp_rhoh_id
module FElib / Element / Base
module FElib / File / History
module FElib / File / Restart
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Data / 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 cloud microphysics 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.