18 subroutine tempo_init(is_aerosol_aware_in, merra2_aerosol_aware_in, is_hail_aware_in, &
19 mpicomm, mpirank, mpiroot, threads, errmsg, errflg)
21 logical,
intent(in) :: is_aerosol_aware_in
22 logical,
intent(in) :: merra2_aerosol_aware_in
23 logical,
intent(in),
optional :: is_hail_aware_in
24 type(mpi_comm),
intent(in) :: mpicomm
25 integer,
intent(in) :: mpirank, mpiroot
26 integer,
intent(in) :: threads
27 character(len=*),
intent(inout) :: errmsg
28 integer,
intent(inout) :: errflg
30 integer :: i, j, k, l, m, n
32 real(wp) :: stime, etime
33 logical,
parameter :: precomputed_tables = .false.
36 call mp_tempo_params_init()
39 configs%aerosol_aware = is_aerosol_aware_in
40 merra2_aerosol_aware = merra2_aerosol_aware_in
41 if (
present(is_hail_aware_in))
then
42 configs%hail_aware = is_hail_aware_in
44 configs%hail_aware = .false.
46 if (configs%aerosol_aware .and. merra2_aerosol_aware)
then
47 errmsg =
'Logic error in tempo_init: only one of the two options can be true, ' // &
48 'not both: is_aerosol_aware or merra2_aerosol_aware'
52 if (mpirank==mpiroot)
then
53 if (configs%aerosol_aware)
then
54 write (*,
'(a)')
'Using aerosol-aware version of TEMPO microphysics'
55 else if(merra2_aerosol_aware)
then
56 write (*,
'(a)')
'Using merra2 aerosol-aware version of TEMPO microphysics'
58 write (*,
'(a)')
'Using non-aerosol-aware version of TEMPO microphysics'
64 if (configs%hail_aware)
then
67 av_g(idx_bg1) = av_g_old
68 bv_g(idx_bg1) = bv_g_old
72 if (mpirank==mpiroot)
then
73 write (*,*)
'Hail-aware option is: ', configs%hail_aware
74 write (*,*)
'Hail-aware option dimNRHG is: ', dimnrhg
78 if (.not.
allocated(tcg_racg))
then
79 allocate(tcg_racg(ntb_g1,ntb_g,dimnrhg,ntb_r1,ntb_r))
84 if (.not.
allocated(tmr_racg))
allocate(tmr_racg(ntb_g1,ntb_g,dimnrhg,ntb_r1,ntb_r))
85 if (.not.
allocated(tcr_gacr))
allocate(tcr_gacr(ntb_g1,ntb_g,dimnrhg,ntb_r1,ntb_r))
86 if (.not.
allocated(tnr_racg))
allocate(tnr_racg(ntb_g1,ntb_g,dimnrhg,ntb_r1,ntb_r))
87 if (.not.
allocated(tnr_gacr))
allocate(tnr_gacr(ntb_g1,ntb_g,dimnrhg,ntb_r1,ntb_r))
90 if (.not.
allocated(tcs_racs1))
allocate(tcs_racs1(ntb_s,ntb_t,ntb_r1,ntb_r))
91 if (.not.
allocated(tmr_racs1))
allocate(tmr_racs1(ntb_s,ntb_t,ntb_r1,ntb_r))
92 if (.not.
allocated(tcs_racs2))
allocate(tcs_racs2(ntb_s,ntb_t,ntb_r1,ntb_r))
93 if (.not.
allocated(tmr_racs2))
allocate(tmr_racs2(ntb_s,ntb_t,ntb_r1,ntb_r))
94 if (.not.
allocated(tcr_sacr1))
allocate(tcr_sacr1(ntb_s,ntb_t,ntb_r1,ntb_r))
95 if (.not.
allocated(tms_sacr1))
allocate(tms_sacr1(ntb_s,ntb_t,ntb_r1,ntb_r))
96 if (.not.
allocated(tcr_sacr2))
allocate(tcr_sacr2(ntb_s,ntb_t,ntb_r1,ntb_r))
97 if (.not.
allocated(tms_sacr2))
allocate(tms_sacr2(ntb_s,ntb_t,ntb_r1,ntb_r))
98 if (.not.
allocated(tnr_racs1))
allocate(tnr_racs1(ntb_s,ntb_t,ntb_r1,ntb_r))
99 if (.not.
allocated(tnr_racs2))
allocate(tnr_racs2(ntb_s,ntb_t,ntb_r1,ntb_r))
100 if (.not.
allocated(tnr_sacr1))
allocate(tnr_sacr1(ntb_s,ntb_t,ntb_r1,ntb_r))
101 if (.not.
allocated(tnr_sacr2))
allocate(tnr_sacr2(ntb_s,ntb_t,ntb_r1,ntb_r))
104 if (.not.
allocated(tpi_qcfz))
allocate(tpi_qcfz(ntb_c,nbc,ntb_t1,ntb_in))
105 if (.not.
allocated(tni_qcfz))
allocate(tni_qcfz(ntb_c,nbc,ntb_t1,ntb_in))
108 if (.not.
allocated(tpi_qrfz))
allocate(tpi_qrfz(ntb_r,ntb_r1,ntb_t1,ntb_in))
109 if (.not.
allocated(tpg_qrfz))
allocate(tpg_qrfz(ntb_r,ntb_r1,ntb_t1,ntb_in))
110 if (.not.
allocated(tni_qrfz))
allocate(tni_qrfz(ntb_r,ntb_r1,ntb_t1,ntb_in))
111 if (.not.
allocated(tnr_qrfz))
allocate(tnr_qrfz(ntb_r,ntb_r1,ntb_t1,ntb_in))
114 if (.not.
allocated(tps_iaus))
allocate(tps_iaus(ntb_i,ntb_i1))
115 if (.not.
allocated(tni_iaus))
allocate(tni_iaus(ntb_i,ntb_i1))
116 if (.not.
allocated(tpi_ide))
allocate(tpi_ide(ntb_i,ntb_i1))
119 if (.not.
allocated(t_efrw))
allocate(t_efrw(nbr,nbc))
120 if (.not.
allocated(t_efsw))
allocate(t_efsw(nbs,nbc))
123 if (.not.
allocated(tnr_rev))
allocate(tnr_rev(nbr,ntb_r1,ntb_r))
124 if (.not.
allocated(tpc_wev))
allocate(tpc_wev(nbc,ntb_c,nbc))
125 if (.not.
allocated(tnc_wev))
allocate(tnc_wev(nbc,ntb_c,nbc))
128 if (.not.
allocated(tnccn_act))
allocate(tnccn_act(ntb_arc,ntb_arw,ntb_art,ntb_arr,ntb_ark))
131 if_micro_init:
if (micro_init)
then
139 mu_c_l = min(15.0_wp, (1000.e6_wp/nt_c_l + 2.))
140 mu_c_o = min(15.0_wp, (1000.e6_wp/nt_c_o + 2.))
146 d0i = (xm0i/am_i)**(1./bm_i)
147 xm0s = am_s * d0s**bm_s
148 xm0g = am_g(nrhg) * d0g**bm_g
154 cce(2,n) = bm_r + n + 1.
155 cce(3,n) = bm_r + n + 4.
156 cce(4,n) = n + bv_c + 1.
157 cce(5,n) = bm_r + n + bv_c + 1.
158 ccg(1,n) = gamma(cce(1,n))
159 ccg(2,n) = gamma(cce(2,n))
160 ccg(3,n) = gamma(cce(3,n))
161 ccg(4,n) = gamma(cce(4,n))
162 ccg(5,n) = gamma(cce(5,n))
163 ocg1(n) = 1.0 / ccg(1,n)
164 ocg2(n) = 1.0 / ccg(2,n)
168 cie(2) = bm_i + mu_i + 1.
169 cie(3) = bm_i + mu_i + bv_i + 1.
170 cie(4) = mu_i + bv_i + 1.
172 cie(6) = bm_i*0.5 + mu_i + bv_i + 1.
173 cie(7) = bm_i*0.5 + mu_i + 1.
174 cig(1) = gamma(cie(1))
175 cig(2) = gamma(cie(2))
176 cig(3) = gamma(cie(3))
177 cig(4) = gamma(cie(4))
178 cig(5) = gamma(cie(5))
179 cig(6) = gamma(cie(6))
180 cig(7) = gamma(cie(7))
187 cre(3) = bm_r + mu_r + 1.
188 cre(4) = bm_r*2. + mu_r + 1.
189 cre(5) = mu_r + bv_r + 1.
190 cre(6) = bm_r + mu_r + bv_r + 1.
191 cre(7) = bm_r*0.5 + mu_r + bv_r + 1.
192 cre(8) = bm_r + mu_r + bv_r + 3.
193 cre(9) = mu_r + bv_r + 3.
195 cre(11) = 0.5*(bv_r + 5. + 2.*mu_r)
196 cre(12) = bm_r*0.5 + mu_r + 1.
197 cre(13) = bm_r*2. + mu_r + bv_r + 1.
200 crg(n) = gamma(cre(n))
212 cse(4) = bm_s + bv_s + 1.
213 cse(5) = bm_s*2. + bv_s + 1.
214 cse(6) = bm_s*2. + 1.
215 cse(7) = bm_s + mu_s + 1.
216 cse(8) = bm_s + mu_s + 2.
217 cse(9) = bm_s + mu_s + 3.
218 cse(10) = bm_s + mu_s + bv_s + 1.
219 cse(11) = bm_s*2. + mu_s + bv_s + 1.
220 cse(12) = bm_s*2. + mu_s + 1.
222 cse(14) = bm_s + bv_s
224 cse(16) = 1.0 + (1.0 + bv_s)/2.
226 if (original_thompson)
then
227 cse(17) = cse(16) + mu_s + 1.
228 cse(18) = bv_s + mu_s + 3.
230 csg(n) = gamma(cse(n))
233 cse(17) = bm_s + bv_s + 2.
235 csg(n) = gamma(cse(n))
245 cge(3,:) = bm_g + mu_g + 1.
246 cge(4,:) = bm_g*2. + mu_g + 1.
247 cge(10,:) = mu_g + 2.
248 cge(12,:) = bm_g*0.5 + mu_g + 1.
251 cge(5,m) = bm_g*2. + mu_g + bv_g(m) + 1.
252 cge(6,m) = bm_g + mu_g + bv_g(m) + 1.
253 cge(7,m) = bm_g*0.5 + mu_g + bv_g(m) + 1.
254 cge(8,m) = mu_g + bv_g(m) + 1.
255 cge(9,m) = mu_g + bv_g(m) + 3.
256 cge(11,m) = 0.5*(bv_g(m) + 5. + 2.*mu_g)
261 cgg(n,m) = gamma(cge(n,m))
269 oamg(m) = 1.0 / am_g(m)
270 ocmg(m) = oamg(m)**obmg
273 oge1 = 1.0 / cge(1,1)
274 ogg1 = 1.0 / cgg(1,1)
275 ogg2 = 1.0 / cgg(2,1)
276 ogg3 = 1.0 / cgg(3,1)
282 t1_qr_qc = pi * 0.25 * av_r * crg(9)
283 t1_qr_qi = pi * 0.25 * av_r * crg(9)
284 t2_qr_qi = pi * 0.25 * am_r*av_r * crg(8)
290 t1_qs_qc = pi * 0.25 * av_s
293 t1_qs_qi = pi * 0.25 * av_s
296 t1_qr_ev = 0.78 * crg(10)
297 t2_qr_ev = 0.308 * sc3 * sqrt(av_r) * crg(11)
301 t2_qs_sd = 0.28 * sc3 * sqrt(av_s)
304 t1_qs_me = pi * 4. *c_sqrd * olfus * 0.86
305 t2_qs_me = pi * 4. *c_sqrd * olfus * 0.28 * sc3 * sqrt(av_s)
308 t1_qg_sd = 0.86 * cgg(10,1)
312 t1_qg_me = pi * 4. * c_cube * olfus * 0.86 * cgg(10,1)
317 nic2 = nint(log10(r_c(1)))
318 nii2 = nint(log10(r_i(1)))
319 nii3 = nint(log10(nt_i(1)))
320 nir2 = nint(log10(r_r(1)))
321 nir3 = nint(log10(n0r_exp(1)))
322 nis2 = nint(log10(r_s(1)))
323 nig2 = nint(log10(r_g(1)))
324 nig3 = nint(log10(n0g_exp(1)))
325 niin2 = nint(log10(nt_in(1)))
331 dc(n) = dc(n-1) + 1.0e-6_dp
332 dtc(n) = (dc(n) - dc(n-1))
336 call create_bins(numbins=nbi, lowbin=d0i*1.0_dp, highbin=d0s*2.0_dp, &
337 bins=di, deltabins=dti)
340 call create_bins(numbins=nbr, lowbin=d0r*1.0_dp, highbin=0.005_dp, &
341 bins=dr, deltabins=dtr)
344 call create_bins(numbins=nbs, lowbin=d0s*1.0_dp, highbin=0.02_dp, &
345 bins=ds, deltabins=dts)
348 call create_bins(numbins=nbg, lowbin=d0g*1.0_dp, highbin=0.05_dp, &
349 bins=dg, deltabins=dtg)
352 call create_bins(numbins=nbc, lowbin=1.0_dp, highbin=3000.0_dp, &
354 t_nc = t_nc * 1.0e6_dp
355 nic1 = log(t_nc(nbc)/t_nc(1))
361 mpi_communicator = mpicomm
366 if (mpirank==mpiroot)
then
367 thompson_table_writer = .true.
369 thompson_table_writer = .false.
372 precomputed_tables_1:
if (.not.precomputed_tables)
then
381 tcg_racg(i,j,n,k,m) = 0.0_dp
382 tmr_racg(i,j,n,k,m) = 0.0_dp
383 tcr_gacr(i,j,n,k,m) = 0.0_dp
384 tnr_racg(i,j,n,k,m) = 0.0_dp
385 tnr_gacr(i,j,n,k,m) = 0.0_dp
396 tcs_racs1(i,j,k,m) = 0.0_dp
397 tmr_racs1(i,j,k,m) = 0.0_dp
398 tcs_racs2(i,j,k,m) = 0.0_dp
399 tmr_racs2(i,j,k,m) = 0.0_dp
400 tcr_sacr1(i,j,k,m) = 0.0_dp
401 tms_sacr1(i,j,k,m) = 0.0_dp
402 tcr_sacr2(i,j,k,m) = 0.0_dp
403 tms_sacr2(i,j,k,m) = 0.0_dp
404 tnr_racs1(i,j,k,m) = 0.0_dp
405 tnr_racs2(i,j,k,m) = 0.0_dp
406 tnr_sacr1(i,j,k,m) = 0.0_dp
407 tnr_sacr2(i,j,k,m) = 0.0_dp
417 tpi_qrfz(i,j,k,m) = 0.0_dp
418 tni_qrfz(i,j,k,m) = 0.0_dp
419 tpg_qrfz(i,j,k,m) = 0.0_dp
420 tnr_qrfz(i,j,k,m) = 0.0_dp
425 tpi_qcfz(i,j,k,m) = 0.0_dp
426 tni_qcfz(i,j,k,m) = 0.0_dp
434 tps_iaus(i,j) = 0.0_dp
435 tni_iaus(i,j) = 0.0_dp
436 tpi_ide(i,j) = 0.0_dp
452 tnr_rev(i,j,k) = 0.0_dp
460 tpc_wev(i,j,k) = 0.0_dp
461 tnc_wev(i,j,k) = 0.0_dp
471 tnccn_act(i,j,k,l,m) = 1.0
478 if (mpirank==mpiroot)
write (*,*)
'creating microphysics lookup tables ... '
479 if (mpirank==mpiroot)
write (*,
'(a, f5.2, a, f5.2, a, f5.2, a, f5.2)') &
480 ' using: mu_c_o=',mu_c_o,
' mu_i=',mu_i,
' mu_r=',mu_r,
' mu_g=',mu_g
485 if (mpirank==mpiroot)
write(*,*)
' calling table_ccnAct routine'
486 call table_ccnact(errmsg, errflg)
487 if (.not. errflg==0)
return
491 if (mpirank==mpiroot)
write(*,*)
' creating qc collision eff tables'
496 if (mpirank==mpiroot)
write(*,*)
' creating rain evap table'
500 if (mpirank==mpiroot)
write(*,*)
' creating ice converting to snow table'
504 if (mpirank==mpiroot) print
'("Calculating TEMPO tables part 1 took ",f10.3," seconds.")', etime-stime
506 end if precomputed_tables_1
516 xam_g = am_g(idx_bg1)
521 if (mpirank==mpiroot) print
'("Calling radar_init took ",f10.3," seconds.")', etime-stime
523 if_not_iiwarm:
if (.not. iiwarm)
then
525 precomputed_tables_2:
if (.not.precomputed_tables)
then
530 if (mpirank==mpiroot)
write(*,*)
' creating rain collecting graupel table'
532 if (dimnrhg == nrhg)
then
533 call qr_acr_qg_par(dimnrhg, qr_acr_qg_hailaware_file)
534 using_hail_aware_table = .true.
536 call qr_acr_qg_par(dimnrhg, qr_acr_qg_file)
539 if (mpirank==mpiroot) print
'("Computing rain collecting graupel table took ",f10.3," seconds.")', etime-stime
542 if (mpirank==mpiroot)
write (*,*)
' creating rain collecting snow table'
546 if (mpirank==mpiroot) print
'("Computing rain collecting snow table took ",f10.3," seconds.")', etime-stime
549 if (mpirank==mpiroot)
write(*,*)
' creating freezing of water drops table'
551 call freezeh2o_par(threads)
553 if (mpirank==mpiroot) print
'("Computing freezing of water drops table took ",f10.3," seconds.")', etime-stime
556 if (mpirank==mpiroot) print
'("Calculating TEMPO tables part 2 took ",f10.3," seconds.")', etime-stime
558 end if precomputed_tables_2
562 if (mpirank==mpiroot)
write(*,*)
' ... DONE microphysical lookup tables'
573 subroutine tempo_3d_to_1d_driver(qv, qc, qr, qi, qs, qg, qb, ni, nr, nc, ng, &
574 nwfa, nifa, nwfa2d, nifa2d, &
576 p, w, dz, dt_in, dt_inner, &
577 sedi_semi, decfl, lsm, &
581 GRAUPELNC, GRAUPELNCV, SR, &
582 refl_10cm, diagflag, do_radar_ref, &
584 vt_dbz_wt, first_time_step, &
585 re_cloud, re_ice, re_snow, &
586 has_reqc, has_reqi, has_reqs, &
587 aero_ind_fdb, rand_perturb_on, &
589 rand_pert, spp_prt_list, spp_var_list, &
590 spp_stddev_cutoff, n_var_spp, &
591 ids,ide, jds,jde, kds,kde, & ! domain dims
592 ims,ime, jms,jme, kms,kme, & ! memory dims
593 its,ite, jts,jte, kts,kte, & ! tile dims
594 fullradar_diag, istep, nsteps, &
596 ! Extended diagnostics, array pointers
597 ! only associated if ext_diag flag is .true.
601 prw_vcde, tpri_inu, tpri_ide_d, &
602 tpri_ide_s, tprs_ide, tprs_sde_d, &
603 tprs_sde_s, tprg_gde_d, &
604 tprg_gde_s, tpri_iha, tpri_wfz, &
605 tpri_rfz, tprg_rfz, tprs_scw, tprg_scw, &
606 tprg_rcs, tprs_rcs, &
607 tprr_rci, tprg_rcg, &
608 tprw_vcd_c, tprw_vcd_e, tprr_sml, &
609 tprr_gml, tprr_rcg, &
610 tprr_rcs, tprv_rev, tten3, qvten3, &
611 qrten3, qsten3, qgten3, qiten3, niten3, &
612 nrten3, ncten3, qcten3, &
616 integer,
intent(in):: ids,ide, jds,jde, kds,kde, &
617 ims,ime, jms,jme, kms,kme, &
618 its,ite, jts,jte, kts,kte
619 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
intent(inout):: &
620 qv, qc, qr, qi, qs, qg, ni, nr
621 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
optional,
intent(inout):: &
623 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
optional,
intent(in):: &
625 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
optional,
intent(inout):: &
626 nc, nwfa, nifa, qb, ng
627 real(wp),
dimension(ims:ime, jms:jme),
optional,
intent(in):: nwfa2d, nifa2d
628 integer,
dimension(ims:ime, jms:jme),
intent(in):: lsm
629 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
optional,
intent(inout):: &
630 re_cloud, re_ice, re_snow
631 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
intent(inout):: pfils, pflls
632 integer,
intent(in) :: rand_perturb_on, kme_stoch, n_var_spp
633 real(wp),
dimension(:,:),
optional,
intent(in) :: rand_pert
634 real(wp),
dimension(:),
optional,
intent(in) :: spp_prt_list
635 real(wp),
dimension(:),
intent(in),
optional :: spp_stddev_cutoff
636 character(len=10),
optional,
dimension(:),
intent(in) :: spp_var_list
637 integer,
intent(in):: has_reqc, has_reqi, has_reqs
639 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
intent(in):: &
641 real(wp),
dimension(ims:ime, jms:jme),
intent(inout):: &
643 real(wp),
dimension(ims:ime, jms:jme),
optional,
intent(inout):: &
646 GRAUPELNC, GRAUPELNCV
647 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
intent(inout):: &
649 real(wp),
dimension(ims:ime, jms:jme),
intent(inout):: &
651 real(wp),
dimension(ims:ime, kms:kme, jms:jme),
optional,
intent(inout):: &
653 logical,
intent(in) :: first_time_step
654 real(wp),
intent(in):: dt_in, dt_inner
655 logical,
intent(in) :: sedi_semi
656 integer,
intent(in) :: decfl
658 integer,
intent (in) :: istep, nsteps
659 logical,
intent (in) :: fullradar_diag
661 logical,
intent (in) :: ext_diag
662 logical,
optional,
intent(in):: aero_ind_fdb
663 real(wp),
optional,
dimension(:,:,:),
intent(inout):: &
666 prw_vcde, tpri_inu, tpri_ide_d, &
667 tpri_ide_s, tprs_ide, &
668 tprs_sde_d, tprs_sde_s, tprg_gde_d, &
669 tprg_gde_s, tpri_iha, tpri_wfz, &
670 tpri_rfz, tprg_rfz, tprs_scw, tprg_scw, &
671 tprg_rcs, tprs_rcs, &
672 tprr_rci, tprg_rcg, &
673 tprw_vcd_c, tprw_vcd_e, tprr_sml, &
674 tprr_gml, tprr_rcg, &
675 tprr_rcs, tprv_rev, tten3, qvten3, &
676 qrten3, qsten3, qgten3, qiten3, niten3, &
677 nrten3, ncten3, qcten3
680 real(wp),
dimension(kts:kte):: &
681 qv1d, qc1d, qi1d, qr1d, qs1d, qg1d, qb1d, &
682 ni1d, nr1d, nc1d, ng1d, nwfa1d, nifa1d, &
683 t1d, p1d, w1d, dz1d, rho, dBZ, pfil1, pfll1
685 real(wp),
dimension(:),
allocatable:: &
688 prw_vcde1, tpri_inu1, tpri_ide1_d, &
689 tpri_ide1_s, tprs_ide1, &
690 tprs_sde1_d, tprs_sde1_s, tprg_gde1_d, &
691 tprg_gde1_s, tpri_iha1, tpri_wfz1, &
692 tpri_rfz1, tprg_rfz1, tprs_scw1, tprg_scw1,&
693 tprg_rcs1, tprs_rcs1, &
694 tprr_rci1, tprg_rcg1, &
695 tprw_vcd1_c, tprw_vcd1_e, tprr_sml1, &
696 tprr_gml1, tprr_rcg1, &
697 tprr_rcs1, tprv_rev1, tten1, qvten1, &
698 qrten1, qsten1, qgten1, qiten1, niten1, &
699 nrten1, ncten1, qcten1
701 real(wp),
dimension(kts:kte):: re_qc1d, re_qi1d, re_qs1d
703 real(wp),
dimension(its:ite, jts:jte):: pcp_ra, pcp_sn, pcp_gr, pcp_ic
704 real(wp) :: dt, pptrain, pptsnow, pptgraul, pptice
705 real(wp) :: qc_max, qr_max, qs_max, qi_max, qg_max, ni_max, nr_max
706 real(wp) :: ygra1, zans1
707 real(dp) :: lamg, lam_exp, lamr, N0_min, N0_exp
709 real(wp) :: rand1, rand2, rand3, rand_pert_max
711 integer:: imax_qc,imax_qr,imax_qi,imax_qs,imax_qg,imax_ni,imax_nr
712 integer:: jmax_qc,jmax_qr,jmax_qi,jmax_qs,jmax_qg,jmax_ni,jmax_nr
713 integer:: kmax_qc,kmax_qr,kmax_qi,kmax_qs,kmax_qg,kmax_ni,kmax_nr
714 integer:: i_start, j_start, i_end, j_end
715 logical,
optional,
intent(in) :: diagflag
716 integer,
optional,
intent(in) :: do_radar_ref
717 logical :: melti = .false.
721 character(len=*),
optional,
intent( out) :: errmsg
722 integer,
optional,
intent( out) :: errflg
725 if (
present(errmsg)) errmsg =
''
726 if (
present(errflg)) errflg = 0
729 test_only_once:
if (first_time_step .and. istep==1)
then
732 if ( (
present(tt) .and. (
present(th) .or.
present(pii))) .or. &
733 (.not.
present(tt) .and. .not.(
present(th) .and.
present(pii))) )
then
734 if (
present(errmsg) .and.
present(errflg))
then
735 write(errmsg,
'(a)')
'Logic error in tempo_3d_to_1d_driver: provide either tt or th+pii'
739 write(*,
'(a)')
'Logic error in tempo_3d_to_1d_driver: provide either tt or th+pii'
744 if (configs%aerosol_aware .and. (.not.
present(nc) .or. &
745 .not.
present(nwfa) .or. &
746 .not.
present(nifa) .or. &
747 .not.
present(nwfa2d) .or. &
748 .not.
present(nifa2d) ))
then
749 if (
present(errmsg) .and.
present(errflg))
then
750 write(errmsg,
'(*(a))')
'Logic error in tempo_3d_to_1d_driver: provide nc, nwfa, nifa, nwfa2d', &
751 ' and nifa2d for aerosol-aware version of TEMPO microphysics'
755 write(*,
'(*(a))')
'Logic error in tempo_3d_to_1d_driver: provide nc, nwfa, nifa, nwfa2d', &
756 ' and nifa2d for aerosol-aware version of TEMPO microphysics'
759 else if (merra2_aerosol_aware .and. (.not.
present(nc) .or. &
760 .not.
present(nwfa) .or. &
761 .not.
present(nifa) ))
then
762 if (
present(errmsg) .and.
present(errflg))
then
763 write(errmsg,
'(*(a))')
'Logic error in tempo_3d_to_1d_driver: provide nc, nwfa, and nifa', &
764 ' for merra2 aerosol-aware version of TEMPO microphysics'
768 write(*,
'(*(a))')
'Logic error in tempo_3d_to_1d_driver: provide nc, nwfa, and nifa', &
769 ' for merra2 aerosol-aware version of TEMPO microphysics'
772 else if (.not.configs%aerosol_aware .and. .not.merra2_aerosol_aware .and. &
773 (
present(nwfa) .or.
present(nifa) .or.
present(nwfa2d) .or.
present(nifa2d)))
then
774 write(*,*)
'WARNING, nc/nwfa/nifa/nwfa2d/nifa2d present but is_aerosol_aware/merra2_aerosol_aware are FALSE'
776 end if test_only_once
782 allocate_extended_diagnostics:
if (ext_diag)
then
783 allocate (prw_vcdc1(kts:kte))
784 allocate (prw_vcde1(kts:kte))
785 allocate (tpri_inu1(kts:kte))
786 allocate (tpri_ide1_d(kts:kte))
787 allocate (tpri_ide1_s(kts:kte))
788 allocate (tprs_ide1(kts:kte))
789 allocate (tprs_sde1_d(kts:kte))
790 allocate (tprs_sde1_s(kts:kte))
791 allocate (tprg_gde1_d(kts:kte))
792 allocate (tprg_gde1_s(kts:kte))
793 allocate (tpri_iha1(kts:kte))
794 allocate (tpri_wfz1(kts:kte))
795 allocate (tpri_rfz1(kts:kte))
796 allocate (tprg_rfz1(kts:kte))
797 allocate (tprs_scw1(kts:kte))
798 allocate (tprg_scw1(kts:kte))
799 allocate (tprg_rcs1(kts:kte))
800 allocate (tprs_rcs1(kts:kte))
801 allocate (tprr_rci1(kts:kte))
802 allocate (tprg_rcg1(kts:kte))
803 allocate (tprw_vcd1_c(kts:kte))
804 allocate (tprw_vcd1_e(kts:kte))
805 allocate (tprr_sml1(kts:kte))
806 allocate (tprr_gml1(kts:kte))
807 allocate (tprr_rcg1(kts:kte))
808 allocate (tprr_rcs1(kts:kte))
809 allocate (tprv_rev1(kts:kte))
810 allocate (tten1(kts:kte))
811 allocate (qvten1(kts:kte))
812 allocate (qrten1(kts:kte))
813 allocate (qsten1(kts:kte))
814 allocate (qgten1(kts:kte))
815 allocate (qiten1(kts:kte))
816 allocate (niten1(kts:kte))
817 allocate (nrten1(kts:kte))
818 allocate (ncten1(kts:kte))
819 allocate (qcten1(kts:kte))
821 allocate (prw_vcdc1(0))
822 allocate (prw_vcde1(0))
823 allocate (tpri_inu1(0))
824 allocate (tpri_ide1_d(0))
825 allocate (tpri_ide1_s(0))
826 allocate (tprs_ide1(0))
827 allocate (tprs_sde1_d(0))
828 allocate (tprs_sde1_s(0))
829 allocate (tprg_gde1_d(0))
830 allocate (tprg_gde1_s(0))
831 allocate (tpri_iha1(0))
832 allocate (tpri_wfz1(0))
833 allocate (tpri_rfz1(0))
834 allocate (tprg_rfz1(0))
835 allocate (tprs_scw1(0))
836 allocate (tprg_scw1(0))
837 allocate (tprg_rcs1(0))
838 allocate (tprs_rcs1(0))
839 allocate (tprr_rci1(0))
840 allocate (tprg_rcg1(0))
841 allocate (tprw_vcd1_c(0))
842 allocate (tprw_vcd1_e(0))
843 allocate (tprr_sml1(0))
844 allocate (tprr_gml1(0))
845 allocate (tprr_rcg1(0))
846 allocate (tprr_rcs1(0))
847 allocate (tprv_rev1(0))
858 end if allocate_extended_diagnostics
887 ndt = max(nint(dt_in/dt_inner),1)
889 if(dt_in .le. dt_inner) dt= dt_in
894 if (rand_perturb_on .ne. 0)
then
896 select case (spp_var_list(k))
898 rand_pert_max = spp_prt_list(k)*spp_stddev_cutoff(k)
934 j_loop:
do j = j_start, j_end
935 i_loop:
do i = i_start, i_end
953 if (rand_perturb_on .ne. 0)
then
954 if (mod(rand_perturb_on,2) .ne. 0) rand1 = rand_pert(i,1)
955 m = rshift(abs(rand_perturb_on),1)
956 if (mod(m,2) .ne. 0) rand2 = rand_pert(i,1)*2.
957 m = rshift(abs(rand_perturb_on),2)
958 if (mod(m,2) .ne. 0) rand3 = 0.25*(rand_pert(i,1)+rand_pert_max)
959 m = rshift(abs(rand_perturb_on),3)
968 IF (
PRESENT (snowncv) )
THEN
971 IF (
PRESENT (icencv) )
THEN
974 IF (
PRESENT (graupelncv) )
THEN
980 if (
present(tt))
then
983 t1d(k) = th(i,k,j)*pii(i,k,j)
996 rho(k) = roverrv * p1d(k) / (r * t1d(k) * (qv1d(k)+roverrv))
1002 initialize_extended_diagnostics:
if (ext_diag)
then
1040 endif initialize_extended_diagnostics
1043 if (configs%aerosol_aware .or. merra2_aerosol_aware)
then
1046 nwfa1d(k) = nwfa(i,k,j)
1047 nifa1d(k) = nifa(i,k,j)
1052 nc1d(k) = nt_c_l / rho(k)
1054 nc1d(k) = nt_c_o / rho(k)
1056 nwfa1d(k) = nwfa_default
1057 nifa1d(k) = nifa_default
1062 if ((
present(ng)) .and. (
present(qb)))
then
1063 configs%hail_aware = .true.
1071 if (qg1d(k) > r1)
then
1072 ygra1 = log10(max(1.e-9, qg1d(k)*rho(k)))
1073 zans1 = 3.4 + 2.0/7.0*(ygra1+8.0)
1075 n0_exp = max(gonv_min, min(10.0**(zans1), gonv_max))
1076 lam_exp = (n0_exp*am_g(idx_bg1)*cgg(1,1) / (rho(k)*qg1d(k)))**oge1
1077 lamg = lam_exp * (cgg(3,1)*ogg2*ogg1)**obmg
1078 ng1d(k) = cgg(2,1) * ogg3*rho(k) * qg1d(k) * lamg**bm_g / am_g(idx_bg1)
1079 ng1d(k) = max(r2, (ng1d(k)/rho(k)))
1080 qb1d(k) = qg1d(k) / rho_g(idx_bg1)
1089 call mp_tempo_main(qv1d=qv1d, qc1d=qc1d, qi1d=qi1d, qr1d=qr1d, qs1d=qs1d, qg1d=qg1d, qb1d=qb1d, &
1090 ni1d=ni1d, nr1d=nr1d, nc1d=nc1d, ng1d=ng1d, nwfa1d=nwfa1d, nifa1d=nifa1d, t1d=t1d, p1d=p1d, &
1091 w1d=w1d, dzq=dz1d, pptrain=pptrain, pptsnow=pptsnow, pptgraul=pptgraul, pptice=pptice, &
1092 rand1=rand1, rand2=rand3, rand3=rand3, &
1093 ext_diag=ext_diag, sedi_semi=sedi_semi, decfl=decfl, &
1094 prw_vcdc1=prw_vcdc1, &
1095 prw_vcde1=prw_vcde1, &
1096 tpri_inu1=tpri_inu1, tpri_ide1_d=tpri_ide1_d, tpri_ide1_s=tpri_ide1_s, tprs_ide1=tprs_ide1, &
1097 tprs_sde1_d=tprs_sde1_d, tprs_sde1_s=tprs_sde1_s, &
1098 tprg_gde1_d=tprg_gde1_d, tprg_gde1_s=tprg_gde1_s, tpri_iha1=tpri_iha1, tpri_wfz1=tpri_wfz1, &
1099 tpri_rfz1=tpri_rfz1, tprg_rfz1=tprg_rfz1, tprs_scw1=tprs_scw1, tprg_scw1=tprg_scw1, &
1100 tprg_rcs1=tprg_rcs1, tprs_rcs1=tprs_rcs1, tprr_rci1=tprr_rci1, &
1101 tprg_rcg1=tprg_rcg1, tprw_vcd1_c=tprw_vcd1_c, &
1102 tprw_vcd1_e=tprw_vcd1_e, tprr_sml1=tprr_sml1, tprr_gml1=tprr_gml1, tprr_rcg1=tprr_rcg1, &
1103 tprr_rcs1=tprr_rcs1, tprv_rev1=tprv_rev1, &
1104 tten1=tten1, qvten1=qvten1, qrten1=qrten1, qsten1=qsten1, &
1105 qgten1=qgten1, qiten1=qiten1, niten1=niten1, nrten1=nrten1, ncten1=ncten1, qcten1=qcten1, &
1106 pfil1=pfil1, pfll1=pfll1, lsml=lsml, &
1107 kts=kts, kte=kte, dt=dt, ii=i, jj=j, configs=configs)
1110 pcp_ra(i,j) = pcp_ra(i,j) + pptrain
1111 pcp_sn(i,j) = pcp_sn(i,j) + pptsnow
1112 pcp_gr(i,j) = pcp_gr(i,j) + pptgraul
1113 pcp_ic(i,j) = pcp_ic(i,j) + pptice
1114 rainncv(i,j) = pptrain + pptsnow + pptgraul + pptice
1115 rainnc(i,j) = rainnc(i,j) + pptrain + pptsnow + pptgraul + pptice
1116 IF (
PRESENT(snowncv) .AND.
PRESENT(snownc) )
THEN
1118 IF ( .NOT.
PRESENT(icencv) .OR. .NOT.
PRESENT(icenc) )
THEN
1119 snowncv(i,j) = pptsnow + pptice
1120 snownc(i,j) = snownc(i,j) + pptsnow + pptice
1122 snowncv(i,j) = pptsnow
1123 snownc(i,j) = snownc(i,j) + pptsnow
1127 IF (
PRESENT(icencv) .AND.
PRESENT(icenc) )
THEN
1128 icencv(i,j) = pptice
1129 icenc(i,j) = icenc(i,j) + pptice
1131 IF (
PRESENT(graupelncv) .AND.
PRESENT(graupelnc) )
THEN
1132 graupelncv(i,j) = pptgraul
1133 graupelnc(i,j) = graupelnc(i,j) + pptgraul
1135 sr(i,j) = (pptsnow + pptgraul + pptice)/(rainncv(i,j)+1.e-12)
1142 if (configs%aerosol_aware)
then
1143 if (
present (aero_ind_fdb) )
then
1144 if ( .not. aero_ind_fdb)
then
1145 nwfa1d(kts) = nwfa1d(kts) + nwfa2d(i,j)*dt
1146 nifa1d(kts) = nifa1d(kts) + nifa2d(i,j)*dt
1149 nwfa1d(kts) = nwfa1d(kts) + nwfa2d(i,j)*dt
1150 nifa1d(kts) = nifa1d(kts) + nifa2d(i,j)*dt
1155 nwfa(i,k,j) = nwfa1d(k)
1156 nifa(i,k,j) = nifa1d(k)
1160 if (merra2_aerosol_aware)
then
1163 nwfa(i,k,j) = nwfa1d(k)
1164 nifa(i,k,j) = nifa1d(k)
1168 if ((
present(ng)) .and. (
present(qb)))
then
1184 pfils(i,k,j) = pfils(i,k,j) + pfil1(k)
1185 pflls(i,k,j) = pflls(i,k,j) + pfll1(k)
1186 if (
present(tt))
then
1189 th(i,k,j) = t1d(k)/pii(i,k,j)
1192 if (qc1d(k) .gt. qc_max)
then
1197 elseif (qc1d(k) .lt. 0.0)
then
1198 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative qc ', qc1d(k), &
1201 if (qr1d(k) .gt. qr_max)
then
1206 elseif (qr1d(k) .lt. 0.0)
then
1207 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative qr ', qr1d(k), &
1210 if (nr1d(k) .gt. nr_max)
then
1215 elseif (nr1d(k) .lt. 0.0)
then
1216 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative nr ', nr1d(k), &
1219 if (qs1d(k) .gt. qs_max)
then
1224 elseif (qs1d(k) .lt. 0.0)
then
1225 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative qs ', qs1d(k), &
1228 if (qi1d(k) .gt. qi_max)
then
1233 elseif (qi1d(k) .lt. 0.0)
then
1234 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative qi ', qi1d(k), &
1237 if (qg1d(k) .gt. qg_max)
then
1242 elseif (qg1d(k) .lt. 0.0)
then
1243 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative qg ', qg1d(k), &
1246 if (ni1d(k) .gt. ni_max)
then
1251 elseif (ni1d(k) .lt. 0.0)
then
1252 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative ni ', ni1d(k), &
1255 if (qv1d(k) .lt. 0.0)
then
1256 write(*,
'(a,e16.7,a,3i8)')
'WARNING, negative qv ', qv1d(k), &
1258 if (k.lt.kte-2 .and. k.gt.kts+1)
then
1259 write(*,*)
' below and above are: ', qv(i,k-1,j), qv(i,k+1,j)
1260 qv(i,k,j) = max(1.e-7, 0.5*(qv(i,k-1,j) + qv(i,k+1,j)))
1267 assign_extended_diagnostics:
if (ext_diag)
then
1272 prw_vcdc(i,k,j) = prw_vcdc(i,k,j) + prw_vcdc1(k)
1273 prw_vcde(i,k,j) = prw_vcde(i,k,j) + prw_vcde1(k)
1274 tpri_inu(i,k,j) = tpri_inu(i,k,j) + tpri_inu1(k)
1275 tpri_ide_d(i,k,j) = tpri_ide_d(i,k,j) + tpri_ide1_d(k)
1276 tpri_ide_s(i,k,j) = tpri_ide_s(i,k,j) + tpri_ide1_s(k)
1277 tprs_ide(i,k,j) = tprs_ide(i,k,j) + tprs_ide1(k)
1278 tprs_sde_s(i,k,j) = tprs_sde_s(i,k,j) + tprs_sde1_s(k)
1279 tprs_sde_d(i,k,j) = tprs_sde_d(i,k,j) + tprs_sde1_d(k)
1280 tprg_gde_d(i,k,j) = tprg_gde_d(i,k,j) + tprg_gde1_d(k)
1281 tprg_gde_s(i,k,j) = tprg_gde_s(i,k,j) + tprg_gde1_s(k)
1282 tpri_iha(i,k,j) = tpri_iha(i,k,j) + tpri_iha1(k)
1283 tpri_wfz(i,k,j) = tpri_wfz(i,k,j) + tpri_wfz1(k)
1284 tpri_rfz(i,k,j) = tpri_rfz(i,k,j) + tpri_rfz1(k)
1285 tprg_rfz(i,k,j) = tprg_rfz(i,k,j) + tprg_rfz1(k)
1286 tprs_scw(i,k,j) = tprs_scw(i,k,j) + tprs_scw1(k)
1287 tprg_scw(i,k,j) = tprg_scw(i,k,j) + tprg_scw1(k)
1288 tprg_rcs(i,k,j) = tprg_rcs(i,k,j) + tprg_rcs1(k)
1289 tprs_rcs(i,k,j) = tprs_rcs(i,k,j) + tprs_rcs1(k)
1290 tprr_rci(i,k,j) = tprr_rci(i,k,j) + tprr_rci1(k)
1291 tprg_rcg(i,k,j) = tprg_rcg(i,k,j) + tprg_rcg1(k)
1292 tprw_vcd_c(i,k,j) = tprw_vcd_c(i,k,j) + tprw_vcd1_c(k)
1293 tprw_vcd_e(i,k,j) = tprw_vcd_e(i,k,j) + tprw_vcd1_e(k)
1294 tprr_sml(i,k,j) = tprr_sml(i,k,j) + tprr_sml1(k)
1295 tprr_gml(i,k,j) = tprr_gml(i,k,j) + tprr_gml1(k)
1296 tprr_rcg(i,k,j) = tprr_rcg(i,k,j) + tprr_rcg1(k)
1297 tprr_rcs(i,k,j) = tprr_rcs(i,k,j) + tprr_rcs1(k)
1298 tprv_rev(i,k,j) = tprv_rev(i,k,j) + tprv_rev1(k)
1299 tten3(i,k,j) = tten3(i,k,j) + tten1(k)
1300 qvten3(i,k,j) = qvten3(i,k,j) + qvten1(k)
1301 qrten3(i,k,j) = qrten3(i,k,j) + qrten1(k)
1302 qsten3(i,k,j) = qsten3(i,k,j) + qsten1(k)
1303 qgten3(i,k,j) = qgten3(i,k,j) + qgten1(k)
1304 qiten3(i,k,j) = qiten3(i,k,j) + qiten1(k)
1305 niten3(i,k,j) = niten3(i,k,j) + niten1(k)
1306 nrten3(i,k,j) = nrten3(i,k,j) + nrten1(k)
1307 ncten3(i,k,j) = ncten3(i,k,j) + ncten1(k)
1308 qcten3(i,k,j) = qcten3(i,k,j) + qcten1(k)
1311 endif assign_extended_diagnostics
1313 if (ndt>1 .and. it==ndt)
then
1315 sr(i,j) = (pcp_sn(i,j) + pcp_gr(i,j) + pcp_ic(i,j))/(rainnc(i,j)+1.e-12)
1316 rainncv(i,j) = rainnc(i,j)
1317 IF (
PRESENT (snowncv) )
THEN
1318 snowncv(i,j) = snownc(i,j)
1320 IF (
PRESENT (icencv) )
THEN
1321 icencv(i,j) = icenc(i,j)
1323 IF (
PRESENT (graupelncv) )
THEN
1324 graupelncv(i,j) = graupelnc(i,j)
1330 last_step_only:
IF ((ndt>1 .and. it==ndt) .or. &
1331 (nsteps>1 .and. istep==nsteps) .or. &
1332 (nsteps==1 .and. ndt==1))
THEN
1338 diagflag_present:
IF (
PRESENT (diagflag) )
THEN
1339 if (diagflag .and. do_radar_ref == 1)
then
1342 if (fullradar_diag)
then
1348 if (
present(vt_dbz_wt))
then
1349 call calc_refl10cm (qv1d=qv1d, qc1d=qc1d, qr1d=qr1d, nr1d=nr1d, qs1d=qs1d, qg1d=qg1d, &
1350 ng1d=ng1d, qb1d=qb1d, t1d=t1d, p1d=p1d, dbz=dbz, rand1=rand1, kts=kts, kte=kte, ii=i, jj=j, &
1351 melti=melti, vt_dbz=vt_dbz_wt(i,:,j), &
1352 first_time_step=first_time_step, configs=configs)
1354 call calc_refl10cm (qv1d=qv1d, qc1d=qc1d, qr1d=qr1d, nr1d=nr1d, qs1d=qs1d, qg1d=qg1d, &
1355 ng1d=ng1d, qb1d=qb1d, t1d=t1d, p1d=p1d, dbz=dbz, rand1=rand1, kts=kts, kte=kte, ii=i, jj=j, &
1356 melti=melti, configs=configs)
1359 refl_10cm(i,k,j) = max(-35., dbz(k))
1362 ENDIF diagflag_present
1364 IF (has_reqc.ne.0 .and. has_reqi.ne.0 .and. has_reqs.ne.0)
THEN
1366 re_qc1d(k) = re_qc_min
1367 re_qi1d(k) = re_qi_min
1368 re_qs1d(k) = re_qs_min
1371 call calc_effectrad (t1d=t1d, p1d=p1d, qv1d=qv1d, qc1d=qc1d, &
1372 nc1d=nc1d, qi1d=qi1d, ni1d=ni1d, qs1d=qs1d, &
1373 re_qc1d=re_qc1d, re_qi1d=re_qi1d, re_qs1d=re_qs1d, &
1374 kts=kts, kte=kte, lsml=lsml, configs=configs)
1376 re_cloud(i,k,j) = max(re_qc_min, min(re_qc1d(k), re_qc_max))
1377 re_ice(i,k,j) = max(re_qi_min, min(re_qi1d(k), re_qi_max))
1378 re_snow(i,k,j) = max(re_qs_min, min(re_qs1d(k), re_qs_max))
1381 ENDIF last_step_only
1398 do j = j_start, j_end
1400 do i = i_start, i_end
1401 pfils(i,k,j) = pfils(i,k,j)/dt_in
1402 pflls(i,k,j) = pflls(i,k,j)/dt_in
1411 deallocate_extended_diagnostics:
if (ext_diag)
then
1412 deallocate (prw_vcdc1)
1413 deallocate (prw_vcde1)
1414 deallocate (tpri_inu1)
1415 deallocate (tpri_ide1_d)
1416 deallocate (tpri_ide1_s)
1417 deallocate (tprs_ide1)
1418 deallocate (tprs_sde1_d)
1419 deallocate (tprs_sde1_s)
1420 deallocate (tprg_gde1_d)
1421 deallocate (tprg_gde1_s)
1422 deallocate (tpri_iha1)
1423 deallocate (tpri_wfz1)
1424 deallocate (tpri_rfz1)
1425 deallocate (tprg_rfz1)
1426 deallocate (tprs_scw1)
1427 deallocate (tprg_scw1)
1428 deallocate (tprg_rcs1)
1429 deallocate (tprs_rcs1)
1430 deallocate (tprr_rci1)
1431 deallocate (tprg_rcg1)
1432 deallocate (tprw_vcd1_c)
1433 deallocate (tprw_vcd1_e)
1434 deallocate (tprr_sml1)
1435 deallocate (tprr_gml1)
1436 deallocate (tprr_rcg1)
1437 deallocate (tprr_rcs1)
1438 deallocate (tprv_rev1)
1449 end if deallocate_extended_diagnostics