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