FE-Project
Loading...
Searching...
No Matches
mod_atmos_phy_mp.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_prc
19 use scale_io
20 use scale_prof
21 use scale_const, only: &
22 undef8 => const_undef8
23 use scale_tracer, only: qa
24
25 use scale_sparsemat, only: sparsemat
27
28 use scale_mesh_base, only: meshbase
31
37
38 use scale_meshfield_base, only: &
40 use scale_localmeshfield_base, only: &
42
43 use scale_model_mesh_manager, only: modelmeshbase
46
48
49
50 !-----------------------------------------------------------------------------
51 implicit none
52 private
53 !-----------------------------------------------------------------------------
54 !
55 !++ Public type & procedure
56 !
57
60 type, extends(modelcomponentproc), public :: atmosphymp
61 integer :: mp_typeid
62 type(atmosphympvars) :: vars
63
64 logical, private :: do_precipitation
65 logical, private :: do_negative_fixer
66 real(rp), private :: limit_negative
67 integer, private :: ntmax_sedimentation
68 real(rp), private :: max_term_vel
69 real(rp), private :: cldfrac_thleshold
70
71 real(dp), private :: dtsec
72 integer, private :: nstep_sedmientation
73 real(rp), private :: rnstep_sedmientation
74 real(dp), private :: dtsec_sedmientation
75
76 type(sparsemat) :: dz, lift
77 type(hexahedralelement) :: elem
78 contains
79 procedure, public :: setup => atmosphymp_setup
80 procedure, public :: calc_tendency => atmosphymp_calc_tendency
81 procedure, public :: update => atmosphymp_update
82 procedure, public :: finalize => atmosphymp_finalize
83 procedure, private :: calc_tendency_core => atmosphymp_calc_tendency_core
84 end type atmosphymp
85
86 !-----------------------------------------------------------------------------
87 !++ Public parameters & variables
88 !-----------------------------------------------------------------------------
89
90 integer, parameter :: mp_typeid_kessler = 1
91
92 !-----------------------------------------------------------------------------
93 !
94 !++ Private procedure
95 !
96 !-----------------------------------------------------------------------------
97 !
98 !++ Private parameters & variables
99 !
100 !-----------------------------------------------------------------------------
101
102contains
103
109 subroutine atmosphymp_setup( this, model_mesh, tm_parent_comp )
110 use scale_const, only: &
111 eps => const_eps
112 use scale_atmos_hydrometeor, only: &
113 atmos_hydrometeor_regist
114 use scale_atmos_phy_mp_kessler, only: &
115 atmos_phy_mp_kessler_setup, &
116 atmos_phy_mp_kessler_ntracers, &
117 atmos_phy_mp_kessler_nwaters, &
118 atmos_phy_mp_kessler_nices, &
119 atmos_phy_mp_kessler_tracer_names, &
120 atmos_phy_mp_kessler_tracer_descriptions, &
121 atmos_phy_mp_kessler_tracer_units
122 !------------------------------
124 use mod_atmos_mesh, only: atmosmesh
125
126 implicit none
127 class(atmosphymp), intent(inout) :: this
128 class(modelmeshbase), target, intent(in) :: model_mesh
129 class(time_manager_component), intent(inout) :: tm_parent_comp
130
131 real(dp) :: time_dt = undef8
132 character(len=H_SHORT) :: time_dt_unit = 'SEC'
133
134 character(len=H_MID) :: mp_type = 'KESSLER'
135
136 logical :: do_precipitation
137 logical :: do_negative_fixer
138 real(rp) :: limit_negative
139 integer :: ntmax_sedimentation
140 real(rp) :: max_term_vel
141 real(rp) :: cldfrac_thleshold
142
143 namelist /param_atmos_phy_mp/ &
144 time_dt, &
145 time_dt_unit, &
146 mp_type, &
147 do_precipitation, &
148 do_negative_fixer, &
149 limit_negative, &
150 ntmax_sedimentation, &
151 max_term_vel, &
152 cldfrac_thleshold
153
154 class(atmosmesh), pointer :: atm_mesh
155 class(meshbase), pointer :: ptr_mesh
156
157 integer :: ierr
158
159 integer :: qs_mp, qe_mp, qa_mp
160
161 integer :: nstep_max
162 !--------------------------------------------------
163
164 if (.not. this%IsActivated()) return
165
166 log_newline
167 log_info("ATMOS_PHY_MP_setup",*) 'Setup'
168
169 cldfrac_thleshold = eps
170
171 do_precipitation = .true.
172 do_negative_fixer = .true.
173 limit_negative = 0.1_rp
174 ntmax_sedimentation = 1
175 max_term_vel = 10.0_rp
176
177 !--- read namelist
178 rewind(io_fid_conf)
179 read(io_fid_conf,nml=param_atmos_phy_mp,iostat=ierr)
180 if( ierr < 0 ) then !--- missing
181 log_info("ATMOS_PHY_MP_setup",*) 'Not found namelist. Default used.'
182 elseif( ierr > 0 ) then !--- fatal error
183 log_error("ATMOS_PHY_MP_setup",*) 'Not appropriate names in namelist PARAM_ATMOS_PHY_MP. Check!'
184 call prc_abort
185 endif
186 log_nml(param_atmos_phy_mp)
187
188 this%cldfrac_thleshold = cldfrac_thleshold
189 this%do_precipitation = do_precipitation
190 this%do_negative_fixer = do_negative_fixer
191 this%ntmax_sedimentation = ntmax_sedimentation
192 this%max_term_vel = max_term_vel
193
194 log_newline
195 log_info("ATMOS_PHY_MP_setup",*) 'Enable negative fixer? : ', this%do_negative_fixer
196 log_info("ATMOS_PHY_MP_setup",*) 'Value limit of negative fixer (abs) : ', abs(this%limit_negative)
197 log_info("ATMOS_PHY_MP_setup",*) 'Enable sedimentation (precipitation)? : ', this%do_precipitation
198
199 !- get mesh --------------------------------------------------
200
201 call model_mesh%GetModelMesh( ptr_mesh )
202 select type(model_mesh)
203 class is (atmosmesh)
204 atm_mesh => model_mesh
205 end select
206
207 !--- Regist this compoent in the time manager
208
209 call tm_parent_comp%Regist_process( 'ATMOS_PHY_MP', time_dt, time_dt_unit, & ! (in)
210 this%tm_process_id ) ! (out)
211
212 this%dtsec = tm_parent_comp%process_list(this%tm_process_id)%dtsec
213 nstep_max = 0 ! ceiling( this%dtsec * this%max_term_vel / maxval( CDZ) )
214 this%ntmax_sedimentation = max( ntmax_sedimentation, nstep_max )
215
216 this%nstep_sedmientation = ntmax_sedimentation
217 this%rnstep_sedmientation = 1.0_rp / real(ntmax_sedimentation,kind=rp)
218 this%dtsec_sedmientation = this%dtsec * this%rnstep_sedmientation
219
220 log_newline
221 log_info("ATMOS_PHY_MP_setup",*) 'Enable negative fixer? : ', this%do_negative_fixer
222 !LOG_INFO("ATMOS_PHY_MP_setup",*) 'Value limit of negative fixer (abs) : ', abs(MP_limit_negative)
223 log_info("ATMOS_PHY_MP_setup",*) 'Enable sedimentation (precipitation)? : ', this%do_precipitation
224 log_info("ATMOS_PHY_MP_setup",*) 'Timestep of sedimentation is divided into : ', this%ntmax_sedimentation, 'step'
225 log_info("ATMOS_PHY_MP_setup",*) 'DT of sedimentation : ', this%dtsec_sedmientation, '[s]'
226
227 !--- Set the type of microphysics scheme
228
229 select case( mp_type )
230 case ('KESSLER')
231 this%MP_TYPEID = mp_typeid_kessler
232 call atmos_hydrometeor_regist( &
233 atmos_phy_mp_kessler_nwaters, & ! (in)
234 atmos_phy_mp_kessler_nices, & ! (in)
235 atmos_phy_mp_kessler_tracer_names(:), & ! (in)
236 atmos_phy_mp_kessler_tracer_descriptions(:), & ! (in)
237 atmos_phy_mp_kessler_tracer_units(:), & ! (in)
238 qs_mp ) ! (out)
239 qa_mp = atmos_phy_mp_kessler_ntracers
240
241 case default
242 log_error("ATMOS_PHY_MP_setup",*) 'Not appropriate names of MP_TYPE in namelist PARAM_ATMOS_PHY_MP. Check!'
243 call prc_abort
244 end select
245
246 qe_mp = qs_mp + qa_mp - 1
247
248 !- initialize the variables
249 call this%vars%Init( model_mesh, qs_mp, qe_mp, qa_mp )
250
251 !- Setup a module for microcloud physics
252
253 select case( this%MP_TYPEID )
254 case( mp_typeid_kessler )
255 call atmos_phy_mp_kessler_setup()
256 end select
257
258 !
259 call this%elem%Init( atm_mesh%ptr_mesh%refElem3D%PolyOrder_h, atm_mesh%ptr_mesh%refElem3D%PolyOrder_v, .true. )
260 call this%Dz%Init( this%elem%Dx3, storage_format='ELL' )
261 call this%Lift%Init( this%elem%Lift, storage_format='ELL' )
262
263 return
264 end subroutine atmosphymp_setup
265
275!OCL SERIAL
276 subroutine atmosphymp_calc_tendency( &
277 this, model_mesh, prgvars_list, trcvars_list, &
278 auxvars_list, forcing_list, is_update )
279
280 use mod_atmos_vars, only: &
285 use mod_atmos_phy_mp_vars, only: &
288
289 implicit none
290
291 class(atmosphymp), intent(inout) :: this
292 class(modelmeshbase), intent(in) :: model_mesh
293 class(modelvarmanager), intent(inout) :: prgvars_list
294 class(modelvarmanager), intent(inout) :: trcvars_list
295 class(modelvarmanager), intent(inout) :: auxvars_list
296 class(modelvarmanager), intent(inout) :: forcing_list
297 logical, intent(in) :: is_update
298
299 class(meshbase), pointer :: mesh
300 class(localmesh3d), pointer :: lcmesh
301
302 class(localmeshfieldbase), pointer :: ddens, momx, momy, momz, drhot
303 class(localmeshfieldbase), pointer :: rtot, cvtot, cptot
304 type(localmeshfieldbaselist) :: qtrc(this%vars%qs:this%vars%qe)
305 class(localmeshfieldbase), pointer :: dens_hyd, pres_hyd
306 class(localmeshfieldbase), pointer :: pres, pt
307 class(localmeshfieldbase), pointer :: dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp, rhoh_p
308 type(localmeshfieldbaselist) :: rhoq_tp(qa)
309 class(localmeshfieldbase), pointer :: mp_dens_t, mp_momx_t, mp_momy_t, mp_momz_t, mp_rhot_t, mp_rhoh, mp_evap
310 type(localmeshfieldbaselist) :: mp_rhoq_t(this%vars%qs:this%vars%qe)
311 class(localmeshfieldbase), pointer :: sflx_rain, sflx_snow, sflx_engi
312
313 integer :: n
314 integer :: ke
315 integer :: iq
316 !--------------------------------------------------
317
318 if (.not. this%IsActivated()) return
319
320 log_progress(*) 'atmosphere / physics / cloud microphysics'
321
322 call model_mesh%GetModelMesh( mesh )
323
324 do n=1, mesh%LOCAL_MESH_NUM
325 call prof_rapstart('ATM_PHY_MP_get_localmesh_ptr', 2)
327 mesh, prgvars_list, auxvars_list, &
328 ddens, momx, momy, momz, drhot, &
329 dens_hyd, pres_hyd, rtot, cvtot, cptot, &
330 lcmesh )
331
333 mesh, trcvars_list, this%vars%QS, qtrc(:) )
334
336 mesh, auxvars_list, &
337 pres, pt )
338
340 mesh, forcing_list, &
341 dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp, &
342 rhoh_p, rhoq_tp )
343
345 mesh, this%vars%tends_manager, &
346 mp_dens_t, mp_momx_t, mp_momy_t, mp_momz_t, mp_rhot_t, mp_rhoh, mp_evap, &
347 mp_rhoq_t )
348
350 mesh, this%vars%auxvars2D_manager, &
351 sflx_rain, sflx_snow, sflx_engi )
352
353 call prof_rapend('ATM_PHY_MP_get_localmesh_ptr', 2)
354
355 call prof_rapstart('ATM_PHY_MP_cal_tend', 2)
356 if (is_update) then
357 call this%calc_tendency_core( &
358 mp_dens_t%val, mp_momx_t%val, mp_momy_t%val, mp_momz_t%val, mp_rhoq_t, & ! (out)
359 mp_rhoh%val, mp_evap%val, sflx_rain%val, sflx_snow%val, sflx_engi%val, & ! (out)
360 ddens%val, momx%val, momy%val, momz%val, pt%val, qtrc, & ! (in)
361 pres%val, pres_hyd%val, dens_hyd%val, rtot%val, cvtot%val, cptot%val, & ! (in)
362 model_mesh%DOptrMat(3), model_mesh%LiftOptrMat, & ! (in)
363 lcmesh, lcmesh%refElem3D, lcmesh%lcmesh2D, lcmesh%lcmesh2D%refElem2D ) ! (in)
364 end if
365
366 !$omp parallel private(ke, iq)
367 !$omp do
368 do ke = lcmesh%NeS, lcmesh%NeE
369 dens_tp%val(:,ke) = dens_tp%val(:,ke) + mp_dens_t%val(:,ke)
370 momx_tp%val(:,ke) = momx_tp%val(:,ke) + mp_momx_t%val(:,ke)
371 momy_tp%val(:,ke) = momy_tp%val(:,ke) + mp_momy_t%val(:,ke)
372 momz_tp%val(:,ke) = momz_tp%val(:,ke) + mp_momz_t%val(:,ke)
373 rhoh_p %val(:,ke) = rhoh_p %val(:,ke) + mp_rhoh %val(:,ke)
374 end do
375 !$omp end do
376 !$omp do collapse(2)
377 do iq = this%vars%QS, this%vars%QE
378 do ke = lcmesh%NeS, lcmesh%NeE
379 rhoq_tp(iq)%ptr%val(:,ke) = rhoq_tp(iq)%ptr%val(:,ke) &
380 + mp_rhoq_t(iq)%ptr%val(:,ke)
381 end do
382 end do
383 !$omp end do
384 !$omp end parallel
385
386 call prof_rapend('ATM_PHY_MP_cal_tend', 2)
387 end do
388
389 return
390 end subroutine atmosphymp_calc_tendency
391
400!OCL SERIAL
401 subroutine atmosphymp_update( this, model_mesh, &
402 prgvars_list, trcvars_list, &
403 auxvars_list, forcing_list, is_update )
404
405 implicit none
406 class(atmosphymp), intent(inout) :: this
407 class(modelmeshbase), intent(in) :: model_mesh
408 class(modelvarmanager), intent(inout) :: prgvars_list
409 class(modelvarmanager), intent(inout) :: trcvars_list
410 class(modelvarmanager), intent(inout) :: auxvars_list
411 class(modelvarmanager), intent(inout) :: forcing_list
412 logical, intent(in) :: is_update
413 !--------------------------------------------------
414
415 return
416 end subroutine atmosphymp_update
417
420!OCL SERIAL
421 subroutine atmosphymp_finalize( this )
422 implicit none
423 class(atmosphymp), intent(inout) :: this
424
425 !--------------------------------------------------
426 if (.not. this%IsActivated()) return
427
428 select case ( this%MP_TYPEID )
429 case( mp_typeid_kessler )
430 end select
431
432 call this%vars%Final()
433
434 return
435 end subroutine atmosphymp_finalize
436
437!- private ------------------------------------------------
438
439!OCL SERIAL
440 subroutine atmosphymp_calc_tendency_core( this, &
441 DENS_t_MP, RHOU_t_MP, RHOV_t_MP, MOMZ_t_MP, RHOQ_t_MP, & ! (out)
442 rhoh_mp, evaporate, sflx_rain, sflx_snow, sflx_engi, & ! (out)
443 ddens, rhou, rhov, momz, pt, qtrc, & ! (in)
444 pres, pres_hyd, dens_hyd, rtot, cvtot, cptot, & ! (in)
445 dz, lift, & ! (in)
446 lcmesh, elem3d, lcmesh2d, elem2d ) ! (in)
447
448 use scale_const, only: &
449 pre00 => const_pre00
450 use scale_const, only: &
451 cvdry => const_cvdry, &
452 cpdry => const_cpdry
453 use scale_atmos_hydrometeor, only: &
454 lhf, &
455 qha, qla, qia
456 use scale_atmos_phy_mp_kessler, only: &
457 atmos_phy_mp_kessler_adjustment, &
458 atmos_phy_mp_kessler_terminal_velocity
459
460 use scale_sparsemat, only: sparsemat
461 use scale_atm_phy_mp_dgm_common, only: &
466
467 implicit none
468
469 class(atmosphymp), intent(inout) :: this
470 class(localmesh3d), intent(in) :: lcmesh
471 class(elementbase3d), intent(in) :: elem3d
472 class(localmesh2d), intent(in) :: lcmesh2d
473 class(elementbase2d), intent(in) :: elem2d
474 real(rp), intent(out) :: dens_t_mp(elem3d%np,lcmesh%nea)
475 real(rp), intent(out) :: rhou_t_mp(elem3d%np,lcmesh%nea)
476 real(rp), intent(out) :: rhov_t_mp(elem3d%np,lcmesh%nea)
477 real(rp), intent(out) :: momz_t_mp(elem3d%np,lcmesh%nea)
478 type(localmeshfieldbaselist), intent(inout) :: rhoq_t_mp(this%vars%qs:this%vars%qe)
479 real(rp), intent(out) :: rhoh_mp(elem3d%np,lcmesh%nea)
480 real(rp), intent(out) :: evaporate(elem3d%np,lcmesh%nea)
481 real(rp), intent(out) :: sflx_rain(elem2d%np,lcmesh2d%nea)
482 real(rp), intent(out) :: sflx_snow(elem2d%np,lcmesh2d%nea)
483 real(rp), intent(out) :: sflx_engi(elem2d%np,lcmesh2d%nea)
484 real(rp), intent(in) :: ddens(elem3d%np,lcmesh%nea)
485 real(rp), intent(in) :: rhou(elem3d%np,lcmesh%nea)
486 real(rp), intent(in) :: rhov(elem3d%np,lcmesh%nea)
487 real(rp), intent(in) :: momz(elem3d%np,lcmesh%nea)
488 real(rp), intent(in) :: pt (elem3d%np,lcmesh%nea)
489 type(localmeshfieldbaselist), intent(in) :: qtrc(this%vars%qs:this%vars%qe)
490 real(rp), intent(in) :: pres(elem3d%np,lcmesh%nea)
491 real(rp), intent(in) :: pres_hyd(elem3d%np,lcmesh%nea)
492 real(rp), intent(in) :: dens_hyd(elem3d%np,lcmesh%nea)
493 real(rp), intent(in) :: rtot (elem3d%np,lcmesh%nea)
494 real(rp), intent(in) :: cvtot(elem3d%np,lcmesh%nea)
495 real(rp), intent(in) :: cptot(elem3d%np,lcmesh%nea)
496 class(sparsemat), intent(in) :: dz
497 class(sparsemat), intent(in) :: lift
498
499 real(rp) :: rhoe_t(elem3d%np,lcmesh%nea)
500
501 real(rp) :: dens (elem3d%np,lcmesh%nea)
502 real(rp) :: temp1(elem3d%np,lcmesh%nea)
503 real(rp) :: cptot1(elem3d%np,lcmesh%nea)
504 real(rp) :: cvtot1(elem3d%np,lcmesh%nea)
505 real(rp) :: qtrc1(elem3d%np,lcmesh%nea,this%vars%qs:this%vars%qe)
506
507 real(rp) :: cptot_t(elem3d%np,lcmesh%nea)
508 real(rp) :: cvtot_t(elem3d%np,lcmesh%nea)
509
510 real(rp) :: vterm(elem3d%np,lcmesh%nez,lcmesh%ne2d,this%vars%qs+1:this%vars%qe)
511 real(rp) :: dens0(elem3d%np,lcmesh%nez,lcmesh%ne2d)
512 real(rp) :: dens2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
513 real(rp) :: ref_dens(elem3d%np,lcmesh%nez,lcmesh%ne2d)
514 real(rp) :: rhou2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
515 real(rp) :: rhov2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
516 real(rp) :: momz2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
517 real(rp) :: temp2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
518 real(rp) :: pres2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
519 real(rp) :: cptot2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
520 real(rp) :: cvtot2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
521 real(rp) :: rhoe (elem3d%np,lcmesh%nez,lcmesh%ne2d)
522 real(rp) :: rhoe2(elem3d%np,lcmesh%nez,lcmesh%ne2d)
523 real(rp) :: rhoq (elem3d%np,lcmesh%nez,lcmesh%ne2d,this%vars%qs+1:this%vars%qe)
524 real(rp) :: rhoq2(elem3d%np,lcmesh%nez,lcmesh%ne2d,this%vars%qs+1:this%vars%qe)
525 real(rp) :: rhoq2_tmp(elem3d%nnode_v,lcmesh%nez,this%vars%qs+1:this%vars%qe)
526 real(rp) :: dens2_tmp(elem3d%nnode_v,lcmesh%nez)
527 real(rp) :: ref_dens_tmp(elem3d%nnode_v,lcmesh%nez)
528 real(rp) :: vterm_tmp(elem3d%nnode_v,lcmesh%nez,this%vars%qs+1:this%vars%qe)
529 real(rp) :: flx_hydro(elem3d%np,lcmesh%nez,lcmesh%ne2d)
530 real(rp) :: cp_t(elem3d%np), cv_t(elem3d%np)
531
532 integer :: iq
533 integer :: domid
534 integer :: ke
535 integer :: ke2d, ke_z
536 integer :: p2d, pv1d, p
537 integer :: colmask(elem3d%nnode_v)
538
539 integer :: step
540 real(rp) :: rdt_mp
541
542 integer :: vmapm(elem3d%nfptot,lcmesh%nez)
543 integer :: vmapp(elem3d%nfptot,lcmesh%nez)
544 real(rp) :: intweight(elem3d%nfaces,elem3d%nfptot)
545 real(rp) :: nz(elem3d%nfptot,lcmesh%nez,lcmesh%ne2d)
546 !----------------------------------------------
547
548 rdt_mp = 1.0_rp / this%dtsec
549 domid = lcmesh%lcdomID
550
551 call atm_phy_mp_dgm_common_gen_vmap( vmapm, vmapp, & ! (out)
552 lcmesh, elem3d ) ! (in)
553
554 call atm_phy_mp_dgm_common_gen_intweight( intweight, & ! (out)
555 lcmesh ) ! (in)
556
557 select case( this%MP_TYPEID )
558 case( mp_typeid_kessler )
559
560 !$omp parallel private(ke, iq)
561 !$omp do
562 do ke = lcmesh%NeS, lcmesh%NeE
563 dens(:,ke) = dens_hyd(:,ke) + ddens(:,ke)
564 temp1(:,ke) = pres(:,ke) / ( dens(:,ke) * rtot(:,ke) )
565 cptot1(:,ke) = cptot(:,ke)
566 cvtot1(:,ke) = cvtot(:,ke)
567 end do
568 !$omp end do
569 !$omp do collapse(2)
570 do iq = this%vars%QS, this%vars%QE
571 do ke = lcmesh%NeS, lcmesh%NeE
572 qtrc1(:,ke,iq) = qtrc(iq)%ptr%val(:,ke)
573 end do
574 end do
575 !$omp end do
576 !$omp end parallel
577
578 call atmos_phy_mp_kessler_adjustment( &
579 elem3d%Np, 1, elem3d%Np, lcmesh%NeA, lcmesh%NeS, lcmesh%NeE, 1, 1, 1, & ! (in)
580 dens, pres, this%dtsec, & ! (in)
581 temp1, qtrc1, cptot1, cvtot1, & ! (inout)
582 rhoe_t, evaporate ) ! (out)
583
584 !$omp parallel private(ke, iq)
585 !$omp do collapse(2)
586 do iq = this%vars%QS, this%vars%QE
587 do ke = lcmesh%NeS, lcmesh%NeE
588 rhoq_t_mp(iq)%ptr%val(:,ke) = ( qtrc1(:,ke,iq) - qtrc(iq)%ptr%val(:,ke) ) * dens(:,ke) * rdt_mp
589 end do
590 end do
591 !$omp end do
592 !$omp do
593 do ke = lcmesh%NeS, lcmesh%NeE
594 cptot_t(:,ke) = ( cptot1(:,ke) - cptot(:,ke) ) * rdt_mp
595 cvtot_t(:,ke) = ( cvtot1(:,ke) - cvtot(:,ke) ) * rdt_mp
596 end do
597 !$omp end do
598 !$omp end parallel
599
600 end select
601
602 !$omp parallel do
603 do ke = lcmesh%NeS, lcmesh%NeE
604 rhoh_mp(:,ke) = rhoe_t(:,ke) &
605 - ( cptot_t(:,ke) + log( pres(:,ke) / pre00 ) * ( cvtot(:,ke) / cptot(:,ke) * cptot_t(:,ke) - cvtot_t(:,ke) ) ) &
606 * pres(:,ke) / rtot(:,ke)
607 end do
608
609 if ( this%do_precipitation ) then
610
611 !$omp parallel private(ke,ke2D,ke_z,iq)
612 !$omp do collapse(2)
613 do ke2d = 1, lcmesh%Ne2D
614 do ke_z = 1, lcmesh%NeZ
615 ke = ke2d + (ke_z-1)*lcmesh%Ne2D
616
617 dens0(:,ke_z,ke2d) = dens(:,ke)
618 dens2(:,ke_z,ke2d) = dens(:,ke)
619 ref_dens(:,ke_z,ke2d) = dens_hyd(:,ke)
620 rhou2(:,ke_z,ke2d) = rhou(:,ke)
621 rhov2(:,ke_z,ke2d) = rhov(:,ke)
622 momz2(:,ke_z,ke2d) = momz(:,ke)
623 temp2(:,ke_z,ke2d) = pres(:,ke) / ( dens(:,ke) * rtot(:,ke) )
624 pres2(:,ke_z,ke2d) = pres(:,ke)
625 cptot2(:,ke_z,ke2d) = cptot(:,ke)
626 cvtot2(:,ke_z,ke2d) = cvtot(:,ke)
627 rhoe(:,ke_z,ke2d) = temp2(:,ke_z,ke2d) * cvtot(:,ke) * dens(:,ke)
628 rhoe2(:,ke_z,ke2d) = rhoe(:,ke_z,ke2d)
629
630 nz(:,ke_z,ke2d) = lcmesh%normal_fn(:,ke,3)
631
632 flx_hydro(:,ke_z,ke2d) = 0.0_rp
633 end do
634 end do
635 !$omp end do
636 !$omp do collapse(3)
637 do iq = this%vars%QS+1, this%vars%QE
638 do ke2d = 1, lcmesh%Ne2D
639 do ke_z = 1, lcmesh%NeZ
640 ke = ke2d + (ke_z-1)*lcmesh%Ne2D
641 rhoq(:,ke_z,ke2d,iq) = dens2(:,ke_z,ke2d) * qtrc(iq)%ptr%val(:,ke) &
642 + rhoq_t_mp(iq)%ptr%val(:,ke) * this%dtsec
643 rhoq2(:,ke_z,ke2d,iq) = rhoq(:,ke_z,ke2d,iq)
644 end do
645 end do
646 end do
647 !$omp end do
648 !$omp workshare
649 sflx_rain(:,:) = 0.0_rp
650 sflx_snow(:,:) = 0.0_rp
651 sflx_engi(:,:) = 0.0_rp
652 !$omp end workshare
653 !$omp end parallel
654
655 do iq = this%vars%QS+1, this%vars%QE
656 if ( this%vars%vterm_hist_id(iq) > 0 ) then
657 !$omp parallel do
658 do ke = lcmesh%NeS, lcmesh%NeE
659 this%vars%vterm_hist(iq)%local(domid)%val(:,ke) = 0.0_rp
660 end do
661 end if
662 end do
663
664 do step = 1, this%nstep_sedmientation
665
666 !- Calculate terminal velocity
667
668 select case( this%MP_TYPEID )
669 case( mp_typeid_kessler )
670 do ke2d = 1, lcmesh%Ne2D
671 !$omp parallel do private( &
672 !$omp p2D, ke_z, iq, ColMask, &
673 !$omp DENS2_tmp, REF_DENS_tmp, RHOQ2_tmp, vterm_tmp )
674 do p2d = 1, elem3d%Nnode_h1D**2
675 colmask(:) = elem3d%Colmask(:,p2d)
676 do ke_z = 1, lcmesh%NeZ
677 dens2_tmp(:,ke_z) = dens2(colmask(:),ke_z,ke2d)
678 ref_dens_tmp(:,ke_z) = ref_dens(colmask(:),ke_z,ke2d)
679 end do
680 do iq = this%vars%QS+1, this%vars%QE
681 do ke_z = 1, lcmesh%NeZ
682 rhoq2_tmp(:,ke_z,iq) = rhoq2(colmask(:),ke_z,ke2d,iq)
683 end do
684 end do
685
686 call atmos_phy_mp_kessler_terminal_velocity( &
687 elem3d%Nnode_v*lcmesh%NeZ, 1, elem3d%Nnode_v*lcmesh%NeZ, & ! (in)
688 dens2_tmp(:,:), rhoq2_tmp(:,:,:), ref_dens_tmp(:,:), & ! (in)
689 vterm_tmp(:,:,:) ) ! (out)
690
691 do iq = this%vars%QS+1, this%vars%QE
692 do ke_z = 1, lcmesh%NeZ
693 vterm(colmask(:),ke_z,ke2d,iq) = vterm_tmp(:,ke_z,iq)
694 end do
695 end do
696 end do
697 end do
698
699 end select
700
701 do iq = this%vars%QS+1, this%vars%QE
702 if ( this%vars%vterm_hist_id(iq) > 0 ) then
703 !$omp parallel do collapse(2) private(ke2D, ke_z, ke)
704 do ke2d = 1, lcmesh%Ne2D
705 do ke_z = 1, lcmesh%NeZ
706 ke = ke2d + (ke_z-1)*lcmesh%Ne2D
707 this%vars%vterm_hist(iq)%local(domid)%val(:,ke) = &
708 this%vars%vterm_hist(iq)%local(domid)%val(:,ke) + vterm(:,ke_z,ke2d,iq) * this%rnstep_sedmientation
709 end do
710 end do
711 end if
712 end do
713
714 !- precipiation of hydrometers
715
717 dens2, rhoq2, cptot2, cvtot2, rhoe2, & ! (inout)
718 flx_hydro, sflx_rain, sflx_snow, sflx_engi, & ! (inout)
719 temp2, vterm, & ! (in)
720 this%dtsec_sedmientation, this%rnstep_sedmientation, & ! (in)
721! Dz, Lift, nz, vmapM, vmapP, IntWeight, & ! (in)
722 this%Dz, this%Lift, nz, vmapm, vmapp, intweight, & ! (in)
723 this%vars%QE - this%vars%QS, qla, qia, & ! (in)
724 lcmesh, elem3d )
725
726 !$omp parallel do private(ke2D, ke_z, ke) collapse(2)
727 do ke2d = 1, lcmesh%Ne2D
728 do ke_z = 1, lcmesh%NeZ
729 ke = ke2d + (ke_z-1)*lcmesh%Ne2D
730 temp2(:,ke_z,ke2d) = rhoe2(:,ke_z,ke2d) / ( dens2(:,ke_z,ke2d) * cvtot2(:,ke_z,ke2d) )
731 end do
732 end do
733
734 end do ! end loop for step
735
736 !$omp parallel private(ke2D, ke_z, iq, ke, CP_t, CV_t)
737 !$omp workshare
738 sflx_engi(:,:) = sflx_engi(:,:) - sflx_snow(:,:) * lhf ! moist internal energy
739 !$omp end workshare
740 !$omp do collapse(2)
741 do ke2d = 1, lcmesh%Ne2D
742 do ke_z = 1, lcmesh%NeZ
743 ke = ke2d + (ke_z-1)*lcmesh%Ne2D
744 dens_t_mp(:,ke) = ( dens2(:,ke_z,ke2d) - dens(:,ke) ) * rdt_mp
745
746 cp_t(:) = ( cptot2(:,ke_z,ke2d) - cptot(:,ke) ) * rdt_mp
747 cv_t(:) = ( cvtot2(:,ke_z,ke2d) - cvtot(:,ke) ) * rdt_mp
748 rhoh_mp(:,ke) = rhoh_mp(:,ke) &
749 + ( rhoe2(:,ke_z,ke2d) - rhoe(:,ke_z,ke2d) ) * rdt_mp &
750 - ( cp_t(:) &
751 + log( pres(:,ke) / pre00 ) * ( cvtot(:,ke) / cptot(:,ke) * cp_t(:) - cv_t(:) ) &
752 ) * pres(:,ke) / rtot(:,ke)
753 end do
754 end do
755 !$omp end do
756 !$omp do collapse(2)
757 do iq = this%vars%QS+1, this%vars%QE
758 do ke2d = 1, lcmesh%Ne2D
759 do ke_z = 1, lcmesh%NeZ
760 ke = ke2d + (ke_z-1)*lcmesh%Ne2D
761 rhoq_t_mp(iq)%ptr%val(:,ke) = rhoq_t_mp(iq)%ptr%val(:,ke) &
762 + ( rhoq2(:,ke_z,ke2d,iq) - rhoq(:,ke_z,ke2d,iq) ) * rdt_mp
763 end do
764 end do
765 end do
766 !$omp end do
767 !$omp end parallel
768
769 !- precipiation of momentum
770
772 rhou_t_mp, rhov_t_mp, momz_t_mp, & ! (out)
773 dens0, rhou2, rhov2, momz2, flx_hydro, & ! (in)
774! Dz, Lift, nz, vmapM, vmapP, & ! (in)
775 this%Dz, this%Lift, nz, vmapm, vmapp, & ! (in)
776 lcmesh, elem3d ) ! (in)
777
778 end if ! end if do_precipitation
779
780 return
781 end subroutine atmosphymp_calc_tendency_core
782
783end module mod_atmos_phy_mp
module Atmosphere / Mesh
module ATMOSPHERE physics / Cloud Microphysics
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)
subroutine, public atmosphympvars_getlocalmeshfields_sfcflx(domid, mesh, sfcflx_list, sflx_rain, sflx_snow, sflx_engi)
module ATMOSPHERE physics / cloud microphysics
integer, parameter mp_typeid_kessler
module ATMOSPHERE / Variables
subroutine, public atmosvars_getlocalmeshprgvars(domid, mesh, prgvars_list, auxvars_list, ddens, momx, momy, momz, therm, dens_hyd, pres_hyd, rtot, cvtot, cptot, lcmesh3d)
subroutine, public atmosvars_getlocalmeshqtrcvarlist(domid, mesh, trcvars_list, varid_s, var_list, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphytends(domid, mesh, phytends_list, dens_tp, momx_tp, momy_tp, momz_tp, rhot_tp, rhoh_p, rhoq_tp, lcmesh3d)
subroutine, public atmosvars_getlocalmeshphyauxvars(domid, mesh, phyauxvars_list, pres, pt, lcmesh3d)
module FElib / Atmosphere / Physics cloud microphysics / common
subroutine, public atm_phy_mp_dgm_common_gen_intweight(intweight, lcmesh)
subroutine, public atm_phy_mp_dgm_common_precipitation(dens, rhoq, cptot, cvtot, rhoe, flx_hydro, sflx_rain, sflx_snow, esflx, temp, vterm, dt, rnstep, dz, lift, nz, vmapm, vmapp, intweight, qha, qla, qia, lcmesh, elem)
subroutine, public atm_phy_mp_dgm_common_precipitation_momentum(momu_t, momv_t, momz_t, dens, momu, momv, momz, mflx, dz, lift, nz, vmapm, vmapp, lcmesh, elem)
subroutine, public atm_phy_mp_dgm_common_gen_vmap(vmapm, vmapp, lmesh, elem)
module FElib / Element / Base
module FElib / Element / hexahedron
module FElib / Mesh / Local 2D
module FElib / Mesh / Local 3D
module FElib / Mesh / Local, Base
module FElib / Mesh / Base 2D
module FElib / Mesh / Base 3D
module FElib / Mesh / Base
module FElib / Data / base
FElib / model framework / physics process.
FElib / model framework / mesh manager.
FElib / model framework / variable manager.
module common / sparsemat
module common / time
Derived type to manage a component of cloud microphysics.