FE-Project
Loading...
Searching...
No Matches
mod_atmos_phy_mp_vars.F90
Go to the documentation of this file.
1!-------------------------------------------------------------------------------
2!> module Atmosphere / Physics / Cloud Microphysics
3!!
4!! @par Description
5!! Container for variables with cloud microphysics component
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
22 use scale_mesh_base, only: meshbase
24 use scale_mesh_base3d, only: &
25 meshbase3d, &
26 dimtype_xyz => meshbase3d_dimtypeid_xyz
29 use scale_localmeshfield_base, only: &
31 use scale_meshfield_base, only: &
33
36
38
39 use scale_model_var_manager, only: &
40 modelvarmanager, variableinfo
41 use scale_model_mesh_manager, only: modelmeshbase
42
43 use mod_atmos_mesh, only: atmosmesh
44
45 !-----------------------------------------------------------------------------
46 implicit none
47 private
48
49 !-----------------------------------------------------------------------------
50 !
51 !++ Public type & procedures
52 !
53
54 !> Derived type to manage variables with cloud microphysics component
55 type, public :: atmosphympvars
56 type(meshfield3d), allocatable :: tends(:) !< Array of tendency variables
57 type(modelvarmanager) :: tends_manager !< Object to manage tendencies
58
59 type(meshfield2d), allocatable :: auxvars2d(:) !< Array of 2D auxiliary variables
60 type(modelvarmanager) :: auxvars2d_manager !< Object to manage 2D auxiliary variables
61
62 integer :: qs !< Start index of tracer variables with cloud microphysics
63 integer :: qe !< End index of tracer variables with cloud microphysics
64 integer :: qa !< Number of tracer variables with cloud microphysics
65
66 integer :: tends_num_tot !< Number of tendency variables with cloud microphysics
67 integer, allocatable :: vterm_hist_id(:)
68 type(meshfield3d), allocatable :: vterm_hist(:)
69 contains
70 procedure :: init => atmosphympvars_init
71 procedure :: final => atmosphympvars_final
72 procedure :: history => atmosphympvars_history
73 end type atmosphympvars
74
77
78 !-----------------------------------------------------------------------------
79 !
80 !++ Public variables
81 !
82
83 integer, public, parameter :: atmos_phy_mp_dens_t_id = 1
84 integer, public, parameter :: atmos_phy_mp_momx_t_id = 2
85 integer, public, parameter :: atmos_phy_mp_momy_t_id = 3
86 integer, public, parameter :: atmos_phy_mp_momz_t_id = 4
87 integer, public, parameter :: atmos_phy_mp_rhot_t_id = 5
88 integer, public, parameter :: atmos_phy_mp_rhoh_id = 6
89 integer, public, parameter :: atmos_phy_mp_evaporate_id = 7
90 integer, public, parameter :: atmos_phy_mp_tends_num1 = 7
91
92 type(variableinfo), public :: atmos_phy_mp_tend_vinfo(atmos_phy_mp_tends_num1)
94 variableinfo( atmos_phy_mp_dens_t_id, 'MP_DENS_t', 'tendency of x-momentum in MP process', &
95 'kg/m3/s', 3, 'XYZ', '' ), &
96 variableinfo( atmos_phy_mp_momx_t_id, 'MP_MOMX_t', 'tendency of x-momentum in MP process', &
97 'kg/m2/s2', 3, 'XYZ', '' ), &
98 variableinfo( atmos_phy_mp_momy_t_id, 'MP_MOMY_t', 'tendency of y-momentum in MP process', &
99 'kg/m2/s2', 3, 'XYZ', '' ), &
100 variableinfo( atmos_phy_mp_momz_t_id, 'MP_MOMZ_t', 'tendency of z-momentum in MP process', &
101 'kg/m2/s2', 3, 'XYZ', '' ), &
102 variableinfo( atmos_phy_mp_rhot_t_id, 'MP_RHOT_t', 'tendency of rho*PT in MP process', &
103 'kg/m3.K/s', 3, 'XYZ', '' ), &
104 variableinfo( atmos_phy_mp_rhoh_id , 'mp_RHOH', 'diabatic heating rate in MP process', &
105 'J/kg/s', 3, 'XYZ', '' ), &
106 variableinfo( atmos_phy_mp_evaporate_id, 'mp_EVAPORATE', 'number concentration of evaporated cloud', &
107 'm-3' , 3, 'XYZ', '' ) /
108
109
110 integer, public, parameter :: atmos_phy_mp_aux2d_sflx_rain_id = 1
111 integer, public, parameter :: atmos_phy_mp_aux2d_sflx_snow_id = 2
112 integer, public, parameter :: atmos_phy_mp_aux2d_sflx_engi_id = 3
113 integer, public, parameter :: atmos_phy_mp_aux2d_num = 3
114
117 variableinfo( atmos_phy_mp_aux2d_sflx_rain_id, 'MP_SFLX_RAIN', 'precipitation flux (liquid) in MP process', &
118 'kg/m2/s', 2, 'XY', '' ), &
119 variableinfo( atmos_phy_mp_aux2d_sflx_snow_id, 'MP_SFLX_SNOW', 'precipitation flux (solid) in MP process', &
120 'kg/m2/s', 2, 'XY', '' ), &
121 variableinfo( atmos_phy_mp_aux2d_sflx_engi_id, 'MP_SFLX_ENGI', 'internal energy flux flux in MP process', &
122 'J/m2/s', 2, 'XY', '' ) /
123
124
125 !-----------------------------------------------------------------------------
126 !
127 !++ Private procedures
128 !
129 !-------------------
130
131contains
132
133!> Setup an object to manage variables with a cloud microphysics component
134!OCL SERIAL
135 subroutine atmosphympvars_init( this, model_mesh, &
136 QS_MP, QE_MP, QA_MP )
137
138 use scale_tracer, only: &
139 tracer_name, tracer_desc, tracer_unit
140 use scale_file_history, only: &
141 file_history_reg
142
143 implicit none
144 class(atmosphympvars), target, intent(inout) :: this
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
149
150 integer :: iv
151 integer :: iq
152 integer :: n
153 logical :: reg_file_hist
154
155 class(atmosmesh), pointer :: atm_mesh
156 class(meshbase2d), pointer :: mesh2d
157 class(meshbase3d), pointer :: mesh3d
158
159 type(variableinfo) :: qtrc_tp_vinfo_tmp
160 type(variableinfo) :: qtrc_vterm_vinfo_tmp
161 !--------------------------------------------------
162
163 log_info('AtmosPhyMpVars_Init',*)
164
165 this%QS = qs_mp
166 this%QE = qe_mp
167 this%QA = qa_mp
168 this%TENDS_NUM_TOT = atmos_phy_mp_tends_num1 + qe_mp - qs_mp + 1
169
170 !- Initialize auxiliary and diagnostic variables
171
172 nullify( atm_mesh )
173 select type(model_mesh)
174 class is (atmosmesh)
175 atm_mesh => model_mesh
176 end select
177 mesh3d => atm_mesh%ptr_mesh
178
179 call mesh3d%GetMesh2D( mesh2d )
180
181 !----
182
183 call this%tends_manager%Init()
184 allocate( this%tends(this%TENDS_NUM_TOT) )
185
186 reg_file_hist = .true.
187 do iv = 1, atmos_phy_mp_tends_num1
188 call this%tends_manager%Regist( &
189 atmos_phy_mp_tend_vinfo(iv), mesh3d, &
190 this%tends(iv), reg_file_hist )
191
192 do n = 1, mesh3d%LOCAL_MESH_NUM
193 this%tends(iv)%local(n)%val(:,:) = 0.0_rp
194 end do
195 end do
196
197 qtrc_tp_vinfo_tmp%ndims = 3
198 qtrc_tp_vinfo_tmp%dim_type = 'XYZ'
199 qtrc_tp_vinfo_tmp%STDNAME = ''
200
201 do iq = 1, this%QA
202 iv = atmos_phy_mp_tends_num1 + iq
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'
207
208 reg_file_hist = .true.
209 call this%tends_manager%Regist( &
210 qtrc_tp_vinfo_tmp, mesh3d, &
211 this%tends(iv), reg_file_hist )
212
213 do n = 1, mesh3d%LOCAL_MESH_NUM
214 this%tends(iv)%local(n)%val(:,:) = 0.0_rp
215 end do
216 end do
217
218 !--
219
220 allocate( this%vterm_hist_id(this%QS+1:this%QE) )
221 allocate( this%vterm_hist (this%QS+1:this%QE) )
222
223 qtrc_vterm_vinfo_tmp%ndims = 3
224 qtrc_vterm_vinfo_tmp%dim_type = 'XYZ'
225 qtrc_vterm_vinfo_tmp%STDNAME = ''
226
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 )
234 end do
235
236 !--
237
238 call this%auxvars2D_manager%Init()
239 allocate( this%auxvars2D(atmos_phy_mp_aux2d_num) )
240
241 reg_file_hist = .true.
242 do iv = 1, atmos_phy_mp_aux2d_num
243 call this%auxvars2D_manager%Regist( &
244 atmos_phy_mp_aux2d_vinfo(iv), mesh2d, & ! (in)
245 this%auxvars2D(iv), reg_file_hist ) ! (out)
246
247 do n = 1, mesh3d%LOCAL_MESH_NUM
248 this%auxvars2D(iv)%local(n)%val(:,:) = 0.0_rp
249 end do
250 end do
251
252 return
253 end subroutine atmosphympvars_init
254
255!> Finalize an object to manage variables with cloud microphysics component
256!OCL SERIAL
257 subroutine atmosphympvars_final( this )
258 implicit none
259 class(atmosphympvars), intent(inout) :: this
260
261 integer :: iq
262 !--------------------------------------------------
263
264 log_info('AtmosPhyMpVars_Final',*)
265
266 call this%tends_manager%Final()
267 deallocate( this%tends )
268
269 call this%auxvars2D_manager%Final()
270 deallocate( this%auxvars2D )
271
272 do iq = this%QS+1, this%QE
273 if ( this%vterm_hist_id(iq) > 0 ) call this%vterm_hist(iq)%Final()
274 end do
275 deallocate( this%vterm_hist_id )
276
277 return
278 end subroutine atmosphympvars_final
279
280!OCL SERIAL
281 subroutine atmosphympvars_getlocalmeshfields_tend( domID, mesh, mp_tends_list, &
282 mp_DENS_t, mp_MOMX_t, mp_MOMY_t, mp_MOMZ_t, mp_RHOT_t, mp_RHOH, mp_EVAP, &
283 mp_RHOQ_t, &
284 lcmesh3D &
285 )
286
287 use scale_mesh_base, only: meshbase
289 implicit none
290
291 integer, intent(in) :: domid
292 class(meshbase), intent(in) :: mesh
293 class(modelvarmanager), intent(inout) :: mp_tends_list
294 class(localmeshfieldbase), pointer, intent(out) :: mp_dens_t
295 class(localmeshfieldbase), pointer, intent(out) :: mp_momx_t
296 class(localmeshfieldbase), pointer, intent(out) :: mp_momy_t
297 class(localmeshfieldbase), pointer, intent(out) :: mp_momz_t
298 class(localmeshfieldbase), pointer, intent(out) :: mp_rhot_t
299 class(localmeshfieldbase), pointer, intent(out) :: mp_rhoh
300 class(localmeshfieldbase), pointer, intent(out) :: mp_evap
301 type(localmeshfieldbaselist), intent(out) :: mp_rhoq_t(:)
302 class(localmesh3d), pointer, intent(out), optional :: lcmesh3d
303
304 class(meshfieldbase), pointer :: field
305 class(localmeshbase), pointer :: lcmesh
306
307 integer :: iq
308 !-------------------------------------------------------
309
310 !--
311 call mp_tends_list%Get(atmos_phy_mp_dens_t_id, field)
312 call field%GetLocalMeshField(domid, mp_dens_t)
313
314 call mp_tends_list%Get(atmos_phy_mp_momx_t_id, field)
315 call field%GetLocalMeshField(domid, mp_momx_t)
316
317 call mp_tends_list%Get(atmos_phy_mp_momy_t_id, field)
318 call field%GetLocalMeshField(domid, mp_momy_t)
319
320 call mp_tends_list%Get(atmos_phy_mp_momz_t_id, field)
321 call field%GetLocalMeshField(domid, mp_momz_t)
322
323 call mp_tends_list%Get(atmos_phy_mp_rhot_t_id, field)
324 call field%GetLocalMeshField(domid, mp_rhot_t)
325
326 call mp_tends_list%Get(atmos_phy_mp_rhoh_id, field)
327 call field%GetLocalMeshField(domid, mp_rhoh)
328
329 call mp_tends_list%Get(atmos_phy_mp_evaporate_id, field)
330 call field%GetLocalMeshField(domid, mp_evap)
331
332 !---
333 do iq = 1, size(mp_rhoq_t)
334 call mp_tends_list%Get(atmos_phy_mp_tends_num1 + iq, field)
335 call field%GetLocalMeshField(domid, mp_rhoq_t(iq)%ptr)
336 end do
337
338
339 if (present(lcmesh3d)) then
340 call mesh%GetLocalMesh( domid, lcmesh )
341 nullify( lcmesh3d )
342
343 select type(lcmesh)
344 type is (localmesh3d)
345 if (present(lcmesh3d)) lcmesh3d => lcmesh
346 end select
347 end if
348
349 return
351
352!OCL SERIAL
353 subroutine atmosphympvars_getlocalmeshfields_sfcflx( domID, mesh, sfcflx_list, &
354 SFLX_rain, SFLX_snow, SFLX_engi )
355
356 use scale_mesh_base, only: meshbase
358 implicit none
359
360 integer, intent(in) :: domid
361 class(meshbase), intent(in) :: mesh
362 class(modelvarmanager), intent(inout) :: sfcflx_list
363 class(localmeshfieldbase), pointer, intent(out) :: sflx_rain
364 class(localmeshfieldbase), pointer, intent(out) :: sflx_snow
365 class(localmeshfieldbase), pointer, intent(out) :: sflx_engi
366
367 class(meshfieldbase), pointer :: field
368 !-------------------------------------------------------
369
370 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_rain_id, field)
371 call field%GetLocalMeshField(domid, sflx_rain)
372
373 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_snow_id, field)
374 call field%GetLocalMeshField(domid, sflx_snow)
375
376 call sfcflx_list%Get(atmos_phy_mp_aux2d_sflx_engi_id, field)
377 call field%GetLocalMeshField(domid, sflx_engi)
378
379 return
381
382!OCL SERIAL
383 subroutine atmosphympvars_history( this )
385 implicit none
386 class(atmosphympvars), intent(inout) :: this
387
388 integer :: v
389 integer :: iq
390 integer :: hst_id
391 type(meshfield3d) :: tmp_field
392 class(meshbase3d), pointer :: mesh3d
393 !-------------------------------------------------------------------------
394
395 mesh3d => this%tends(1)%mesh
396
397 do v = 1, this%TENDS_NUM_TOT
398 hst_id = this%tends(v)%hist_id
399 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%tends(v) )
400 end do
401
402 do v = 1, atmos_phy_mp_aux2d_num
403 hst_id = this%auxvars2D(v)%hist_id
404 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%auxvars2D(v) )
405 end do
406
407 do iq = this%QS+1, this%QE
408 hst_id = this%vterm_hist_id(iq)
409 if ( hst_id > 0 ) call file_history_meshfield_put( hst_id, this%vterm_hist(iq) )
410 end do
411
412 return
413 end subroutine atmosphympvars_history
414
415end module mod_atmos_phy_mp_vars
module Atmosphere / Mesh
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 / 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 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.