CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
module_mp_tempo.F90
1! 3D TEMPO Driver for CCPP
2!=================================================================================================================
4
5 use machine, only: wp => kind_phys, sp => kind_sngl_prec, dp => kind_dbl_prec
6 use module_mp_tempo_params
7 use module_mp_tempo_utils, only : create_bins, table_efrw, table_efsw, table_dropevap, &
8 table_ccnact, qi_aut_qs, qr_acr_qg_par, qr_acr_qs_par, freezeh2o_par, calc_refl10cm, calc_effectrad
9 use module_mp_tempo_main, only : mp_tempo_main
11
12 implicit none
13
14contains
15 !=================================================================================================================
16 ! This subroutine handles initialzation of the microphysics scheme including building of lookup tables,
17 ! allocating arrays for the microphysics scheme, and defining gamma function variables.
18 subroutine tempo_init(is_aerosol_aware_in, merra2_aerosol_aware_in, is_hail_aware_in, &
19 mpicomm, mpirank, mpiroot, threads, errmsg, errflg)
20
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
29
30 integer :: i, j, k, l, m, n
31 logical :: micro_init
32 real(wp) :: stime, etime
33 logical, parameter :: precomputed_tables = .false.
34
35 ! Initialize physical constants
36 call mp_tempo_params_init()
37
38 ! Set module variable is_aerosol_aware/merra2_aerosol_aware
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
43 else
44 configs%hail_aware = .false.
45 endif
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'
49 errflg = 1
50 return
51 end if
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'
57 else
58 write (*,'(a)') 'Using non-aerosol-aware version of TEMPO microphysics'
59 end if
60 end if
61
62 micro_init = .false.
63
64 if (configs%hail_aware) then
65 dimnrhg = nrhg
66 else
67 av_g(idx_bg1) = av_g_old
68 bv_g(idx_bg1) = bv_g_old
69 dimnrhg = nrhg1
70 endif
71
72 if (mpirank==mpiroot) then
73 write (*,*) 'Hail-aware option is: ', configs%hail_aware
74 write (*,*) 'Hail-aware option dimNRHG is: ', dimnrhg
75 endif
76
77 ! Allocate space for lookup tables (J. Michalakes 2009Jun08).
78 if (.not. allocated(tcg_racg)) then
79 allocate(tcg_racg(ntb_g1,ntb_g,dimnrhg,ntb_r1,ntb_r))
80 micro_init = .true.
81 endif
82
83 ! Rain-graupel (including array above tcg_racg)
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))
88
89 ! Rain-snow
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))
102
103 ! Cloud water freezing
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))
106
107 ! Rain freezing
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))
112
113 ! Ice growth and conversion to snow
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))
117
118 ! Collision efficiencies
119 if (.not. allocated(t_efrw)) allocate(t_efrw(nbr,nbc))
120 if (.not. allocated(t_efsw)) allocate(t_efsw(nbs,nbc))
121
122 ! Cloud water
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))
126
127 ! CCN
128 if (.not. allocated(tnccn_act)) allocate(tnccn_act(ntb_arc,ntb_arw,ntb_art,ntb_arr,ntb_ark))
129
130 !=================================================================================================================
131 if_micro_init: if (micro_init) then
132
136 !.. disp=SQRT((mu+2)/(mu+1) - 1) so mu varies from 15 for Maritime
137 !.. to 2 for really dirty air. This not used in 2-moment cloud water
138 !.. scheme and nu_c used instead and varies from 2 to 15 (integer-only).
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.))
141
143 sc3 = sc**(1./3.)
144
146 d0i = (xm0i/am_i)**(1./bm_i)
147 xm0s = am_s * d0s**bm_s
148 xm0g = am_g(nrhg) * d0g**bm_g
149
152 do n = 1, 15
153 cce(1,n) = n + 1.
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)
165 enddo
166
167 cie(1) = mu_i + 1.
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.
171 cie(5) = mu_i + 2.
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))
181 oig1 = 1.0 / cig(1)
182 oig2 = 1.0 / cig(2)
183 obmi = 1.0 / bm_i
184
185 cre(1) = bm_r + 1.
186 cre(2) = mu_r + 1.
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.
194 cre(10) = mu_r + 2.
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.
198
199 do n = 1, 13
200 crg(n) = gamma(cre(n))
201 enddo
202
203 obmr = 1.0 / bm_r
204 ore1 = 1.0 / cre(1)
205 org1 = 1.0 / crg(1)
206 org2 = 1.0 / crg(2)
207 org3 = 1.0 / crg(3)
208
209 cse(1) = bm_s + 1.
210 cse(2) = bm_s + 2.
211 cse(3) = bm_s*2.
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.
221 cse(13) = bv_s + 2.
222 cse(14) = bm_s + bv_s
223 cse(15) = mu_s + 1.
224 cse(16) = 1.0 + (1.0 + bv_s)/2.
225
226 if (original_thompson) then
227 cse(17) = cse(16) + mu_s + 1.
228 cse(18) = bv_s + mu_s + 3.
229 do n = 1, 18
230 csg(n) = gamma(cse(n))
231 enddo
232 else
233 cse(17) = bm_s + bv_s + 2.
234 do n = 1, 17
235 csg(n) = gamma(cse(n))
236 enddo
237 endif
238
239 oams = 1.0 / am_s
240 obms = 1.0 / bm_s
241 ocms = oams**obms
242
243 cge(1,:) = bm_g + 1.
244 cge(2,:) = mu_g + 1.
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.
249
250 do m = 1, nrhg
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. ! not used
255 cge(9,m) = mu_g + bv_g(m) + 3.
256 cge(11,m) = 0.5*(bv_g(m) + 5. + 2.*mu_g)
257 enddo
258
259 do m = 1, nrhg
260 do n = 1, 12
261 cgg(n,m) = gamma(cge(n,m))
262 enddo
263 enddo
264
265 oamg = 1.0 / am_g
266 obmg = 1.0 / bm_g
267
268 do m = 1, nrhg
269 oamg(m) = 1.0 / am_g(m)
270 ocmg(m) = oamg(m)**obmg
271 enddo
272
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)
277
278 !=================================================================================================================
279 ! Simplify various rate eqns the best we can now.
280
281 ! Rain collecting cloud water and cloud ice
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)
285
286 ! Graupel collecting cloud water
287 ! t1_qg_qc = PI*.25*av_g * cgg(9)
288
289 ! Snow collecting cloud water
290 t1_qs_qc = pi * 0.25 * av_s
291
292 ! Snow collecting cloud ice
293 t1_qs_qi = pi * 0.25 * av_s
294
295 ! Evaporation of rain; ignore depositional growth of rain.
296 t1_qr_ev = 0.78 * crg(10)
297 t2_qr_ev = 0.308 * sc3 * sqrt(av_r) * crg(11)
298
299 ! Sublimation/depositional growth of snow
300 t1_qs_sd = 0.86
301 t2_qs_sd = 0.28 * sc3 * sqrt(av_s)
302
303 ! Melting of snow
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)
306
307 ! Sublimation/depositional growth of graupel
308 t1_qg_sd = 0.86 * cgg(10,1)
309 ! t2_qg_sd = 0.28*Sc3*SQRT(av_g) * cgg(11)
310
311 ! Melting of graupel
312 t1_qg_me = pi * 4. * c_cube * olfus * 0.86 * cgg(10,1)
313 ! t2_qg_me = PI*4.*C_cube*olfus * 0.28*Sc3*SQRT(av_g) * cgg(11)
314
315
316 ! Constants for helping find lookup table indexes.
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)))
326
327 ! Create bins of cloud water (from minimum diameter to 100 microns).
328 dc(1) = d0c*1.0_dp
329 dtc(1) = d0c*1.0_dp
330 do n = 2, nbc
331 dc(n) = dc(n-1) + 1.0e-6_dp
332 dtc(n) = (dc(n) - dc(n-1))
333 enddo
334
335 ! Create bins of cloud ice (from min diameter up to 2x min snow size).
336 call create_bins(numbins=nbi, lowbin=d0i*1.0_dp, highbin=d0s*2.0_dp, &
337 bins=di, deltabins=dti)
338
339 ! Create bins of rain (from min diameter up to 5 mm).
340 call create_bins(numbins=nbr, lowbin=d0r*1.0_dp, highbin=0.005_dp, &
341 bins=dr, deltabins=dtr)
342
343 ! Create bins of snow (from min diameter up to 2 cm).
344 call create_bins(numbins=nbs, lowbin=d0s*1.0_dp, highbin=0.02_dp, &
345 bins=ds, deltabins=dts)
346
347 ! Create bins of graupel (from min diameter up to 5 cm).
348 call create_bins(numbins=nbg, lowbin=d0g*1.0_dp, highbin=0.05_dp, &
349 bins=dg, deltabins=dtg)
350
351 ! Create bins of cloud droplet number concentration (1 to 3000 per cc).
352 call create_bins(numbins=nbc, lowbin=1.0_dp, highbin=3000.0_dp, &
353 bins=t_nc)
354 t_nc = t_nc * 1.0e6_dp
355 nic1 = log(t_nc(nbc)/t_nc(1))
356
357 !=================================================================================================================
358 ! Create lookup tables for most costly calculations
359
360 ! Assign mpicomm to module variable
361 mpi_communicator = mpicomm
362
363 ! Standard tables are only written by master MPI task;
364 ! (physics init cannot be called by multiple threads,
365 ! hence no need to test for a specific thread number)
366 if (mpirank==mpiroot) then
367 thompson_table_writer = .true.
368 else
369 thompson_table_writer = .false.
370 end if
371
372 precomputed_tables_1: if (.not.precomputed_tables) then
373
374 call cpu_time(stime)
375
376 do m = 1, ntb_r
377 do k = 1, ntb_r1
378 do n = 1, dimnrhg
379 do j = 1, ntb_g
380 do i = 1, ntb_g1
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
386 enddo
387 enddo
388 enddo
389 enddo
390 enddo
391
392 do m = 1, ntb_r
393 do k = 1, ntb_r1
394 do j = 1, ntb_t
395 do i = 1, ntb_s
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
408 enddo
409 enddo
410 enddo
411 enddo
412
413 do m = 1, ntb_in
414 do k = 1, ntb_t1
415 do j = 1, ntb_r1
416 do i = 1, ntb_r
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
421 enddo
422 enddo
423 do j = 1, nbc
424 do i = 1, ntb_c
425 tpi_qcfz(i,j,k,m) = 0.0_dp
426 tni_qcfz(i,j,k,m) = 0.0_dp
427 enddo
428 enddo
429 enddo
430 enddo
431
432 do j = 1, ntb_i1
433 do i = 1, ntb_i
434 tps_iaus(i,j) = 0.0_dp
435 tni_iaus(i,j) = 0.0_dp
436 tpi_ide(i,j) = 0.0_dp
437 enddo
438 enddo
439
440 do j = 1, nbc
441 do i = 1, nbr
442 t_efrw(i,j) = 0.0
443 enddo
444 do i = 1, nbs
445 t_efsw(i,j) = 0.0
446 enddo
447 enddo
448
449 do k = 1, ntb_r
450 do j = 1, ntb_r1
451 do i = 1, nbr
452 tnr_rev(i,j,k) = 0.0_dp
453 enddo
454 enddo
455 enddo
456
457 do k = 1, nbc
458 do j = 1, ntb_c
459 do i = 1, nbc
460 tpc_wev(i,j,k) = 0.0_dp
461 tnc_wev(i,j,k) = 0.0_dp
462 enddo
463 enddo
464 enddo
465
466 do m = 1, ntb_ark
467 do l = 1, ntb_arr
468 do k = 1, ntb_art
469 do j = 1, ntb_arw
470 do i = 1, ntb_arc
471 tnccn_act(i,j,k,l,m) = 1.0
472 enddo
473 enddo
474 enddo
475 enddo
476 enddo
477
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
481
485 if (mpirank==mpiroot) write(*,*) ' calling table_ccnAct routine'
486 call table_ccnact(errmsg, errflg)
487 if (.not. errflg==0) return
488
491 if (mpirank==mpiroot) write(*,*) ' creating qc collision eff tables'
492 call table_efrw
493 call table_efsw
494
496 if (mpirank==mpiroot) write(*,*) ' creating rain evap table'
497 call table_dropevap
498
500 if (mpirank==mpiroot) write(*,*) ' creating ice converting to snow table'
501 call qi_aut_qs
502
503 call cpu_time(etime)
504 if (mpirank==mpiroot) print '("Calculating TEMPO tables part 1 took ",f10.3," seconds.")', etime-stime
505
506 end if precomputed_tables_1
507
509 call cpu_time(stime)
510 xam_r = am_r
511 xbm_r = bm_r
512 xmu_r = mu_r
513 xam_s = am_s
514 xbm_s = bm_s
515 xmu_s = mu_s
516 xam_g = am_g(idx_bg1)
517 xbm_g = bm_g
518 xmu_g = mu_g
519 call radar_init
520 call cpu_time(etime)
521 if (mpirank==mpiroot) print '("Calling radar_init took ",f10.3," seconds.")', etime-stime
522
523 if_not_iiwarm: if (.not. iiwarm) then
524
525 precomputed_tables_2: if (.not.precomputed_tables) then
526
527 call cpu_time(stime)
528
530 if (mpirank==mpiroot) write(*,*) ' creating rain collecting graupel table'
531 call cpu_time(stime)
532 if (dimnrhg == nrhg) then
533 call qr_acr_qg_par(dimnrhg, qr_acr_qg_hailaware_file)
534 using_hail_aware_table = .true.
535 else
536 call qr_acr_qg_par(dimnrhg, qr_acr_qg_file)
537 endif
538 call cpu_time(etime)
539 if (mpirank==mpiroot) print '("Computing rain collecting graupel table took ",f10.3," seconds.")', etime-stime
540
542 if (mpirank==mpiroot) write (*,*) ' creating rain collecting snow table'
543 call cpu_time(stime)
544 call qr_acr_qs_par
545 call cpu_time(etime)
546 if (mpirank==mpiroot) print '("Computing rain collecting snow table took ",f10.3," seconds.")', etime-stime
547
549 if (mpirank==mpiroot) write(*,*) ' creating freezing of water drops table'
550 call cpu_time(stime)
551 call freezeh2o_par(threads)
552 call cpu_time(etime)
553 if (mpirank==mpiroot) print '("Computing freezing of water drops table took ",f10.3," seconds.")', etime-stime
554
555 call cpu_time(etime)
556 if (mpirank==mpiroot) print '("Calculating TEMPO tables part 2 took ",f10.3," seconds.")', etime-stime
557
558 end if precomputed_tables_2
559
560 endif if_not_iiwarm
561
562 if (mpirank==mpiroot) write(*,*) ' ... DONE microphysical lookup tables'
563
564 endif if_micro_init
565
566 end subroutine tempo_init
567
568 !=================================================================================================================
569 ! This is a wrapper routine designed to transfer values from 3D to 1D.
570 ! Required microphysics variables are qv, qc, qr, nr, qi, ni, qs, qg
571 ! Optional microphysics variables are aerosol aware (nc, nwfa, nifa, nwfa2d, nifa2d), and hail aware (ng, qg)
572
573 subroutine tempo_3d_to_1d_driver(qv, qc, qr, qi, qs, qg, qb, ni, nr, nc, ng, &
574 nwfa, nifa, nwfa2d, nifa2d, &
575 tt, th, pii, &
576 p, w, dz, dt_in, dt_inner, &
577 sedi_semi, decfl, lsm, &
578 RAINNC, RAINNCV, &
579 SNOWNC, SNOWNCV, &
580 ICENC, ICENCV, &
581 GRAUPELNC, GRAUPELNCV, SR, &
582 refl_10cm, diagflag, do_radar_ref, &
583 max_hail_diam_sfc, &
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, &
588 kme_stoch, &
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, &
595 errmsg, errflg, &
596 ! Extended diagnostics, array pointers
597 ! only associated if ext_diag flag is .true.
598 ext_diag, &
599 !vts1, txri, txrc, &
600 prw_vcdc, &
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, &
613 pfils, pflls)
614
615 !..Subroutine arguments
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):: &
622 tt, th
623 real(wp), dimension(ims:ime, kms:kme, jms:jme), optional, intent(in):: &
624 pii
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
638
639 real(wp), dimension(ims:ime, kms:kme, jms:jme), intent(in):: &
640 p, w, dz
641 real(wp), dimension(ims:ime, jms:jme), intent(inout):: &
642 RAINNC, RAINNCV, SR
643 real(wp), dimension(ims:ime, jms:jme), optional, intent(inout):: &
644 SNOWNC, SNOWNCV, &
645 ICENC, ICENCV, &
646 GRAUPELNC, GRAUPELNCV
647 real(wp), dimension(ims:ime, kms:kme, jms:jme), intent(inout):: &
648 refl_10cm
649 real(wp), dimension(ims:ime, jms:jme), intent(inout):: &
650 max_hail_diam_sfc
651 real(wp), dimension(ims:ime, kms:kme, jms:jme), optional, intent(inout):: &
652 vt_dbz_wt
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
657 ! To support subcycling: current step and maximum number of steps
658 integer, intent (in) :: istep, nsteps
659 logical, intent (in) :: fullradar_diag
660 ! Extended diagnostics, array pointers only associated if ext_diag flag is .true.
661 logical, intent (in) :: ext_diag
662 logical, optional, intent(in):: aero_ind_fdb
663 real(wp), optional, dimension(:,:,:), intent(inout):: &
664 !vts1, txri, txrc, &
665 prw_vcdc, &
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
678
679 !..Local variables
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
684 !..Extended diagnostics, single column arrays
685 real(wp), dimension(:), allocatable:: &
686 !vtsk1, txri1, txrc1, &
687 prw_vcdc1, &
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
700
701 real(wp), dimension(kts:kte):: re_qc1d, re_qi1d, re_qs1d
702
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
708 integer:: lsml
709 real(wp) :: rand1, rand2, rand3, rand_pert_max
710 integer:: i, j, k, m
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.
718 integer :: ndt, it
719
720 ! CCPP error handling
721 character(len=*), optional, intent( out) :: errmsg
722 integer, optional, intent( out) :: errflg
723
724 ! CCPP
725 if (present(errmsg)) errmsg = ''
726 if (present(errflg)) errflg = 0
727
728 ! No need to test for every subcycling step
729 test_only_once: if (first_time_step .and. istep==1) then
730 ! Activate this code when removing the guard above
731
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'
736 errflg = 1
737 return
738 else
739 write(*,'(a)') 'Logic error in tempo_3d_to_1d_driver: provide either tt or th+pii'
740 stop
741 end if
742 end if
743
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'
752 errflg = 1
753 return
754 else
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'
757 stop
758 end if
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'
765 errflg = 1
766 return
767 else
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'
770 stop
771 end if
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'
775 end if
776 end if test_only_once
777
778 ! These must be alwyas allocated
779 !allocate (vtsk1(kts:kte))
780 !allocate (txri1(kts:kte))
781 !allocate (txrc1(kts:kte))
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))
820 else
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))
848 allocate (tten1(0))
849 allocate (qvten1(0))
850 allocate (qrten1(0))
851 allocate (qsten1(0))
852 allocate (qgten1(0))
853 allocate (qiten1(0))
854 allocate (niten1(0))
855 allocate (nrten1(0))
856 allocate (ncten1(0))
857 allocate (qcten1(0))
858 end if allocate_extended_diagnostics
859
860 !+---+
861 i_start = its
862 j_start = jts
863 i_end = ite
864 j_end = jte
865
866 !..For idealized testing by developer.
867 ! if ( (ide-ids+1).gt.4 .and. (jde-jds+1).lt.4 .and. &
868 ! ids.eq.its.and.ide.eq.ite.and.jds.eq.jts.and.jde.eq.jte) then
869 ! i_start = its + 2
870 ! i_end = ite
871 ! j_start = jts
872 ! j_end = jte
873 ! endif
874
875 ! dt = dt_in
876 rainnc(:,:) = 0.0
877 snownc(:,:) = 0.0
878 icenc(:,:) = 0.0
879 graupelnc(:,:) = 0.0
880 pcp_ra(:,:) = 0.0
881 pcp_sn(:,:) = 0.0
882 pcp_gr(:,:) = 0.0
883 pcp_ic(:,:) = 0.0
884 pfils(:,:,:) = 0.0
885 pflls(:,:,:) = 0.0
886 rand_pert_max = 0.0
887 ndt = max(nint(dt_in/dt_inner),1)
888 dt = dt_in/ndt
889 if(dt_in .le. dt_inner) dt= dt_in
890
891 !Get the Thompson MP SPP magnitude and standard deviation cutoff,
892 !then compute rand_pert_max
893
894 if (rand_perturb_on .ne. 0) then
895 do k =1,n_var_spp
896 select case (spp_var_list(k))
897 case('mp')
898 rand_pert_max = spp_prt_list(k)*spp_stddev_cutoff(k)
899 end select
900 enddo
901 endif
902
903 do it = 1, ndt
904
905 qc_max = 0.
906 qr_max = 0.
907 qs_max = 0.
908 qi_max = 0.
909 qg_max = 0
910 ni_max = 0.
911 nr_max = 0.
912 imax_qc = 0
913 imax_qr = 0
914 imax_qi = 0
915 imax_qs = 0
916 imax_qg = 0
917 imax_ni = 0
918 imax_nr = 0
919 jmax_qc = 0
920 jmax_qr = 0
921 jmax_qi = 0
922 jmax_qs = 0
923 jmax_qg = 0
924 jmax_ni = 0
925 jmax_nr = 0
926 kmax_qc = 0
927 kmax_qr = 0
928 kmax_qi = 0
929 kmax_qs = 0
930 kmax_qg = 0
931 kmax_ni = 0
932 kmax_nr = 0
933
934 j_loop: do j = j_start, j_end
935 i_loop: do i = i_start, i_end
936
937 !+---+-----------------------------------------------------------------+
938 !..Introduce stochastic parameter perturbations by creating as many scalar rand1, rand2, ...
939 !.. variables as needed to perturb different pieces of microphysics. gthompsn 21Mar2018
940 ! Setting spp_mp_opt to 1 gives graupel Y-intercept pertubations (2^0)
941 ! 2 gives cloud water distribution gamma shape parameter perturbations (2^1)
942 ! 4 gives CCN & IN activation perturbations (2^2)
943 ! 3 gives both 1+2
944 ! 5 gives both 1+4
945 ! 6 gives both 2+4
946 ! 7 gives all 1+2+4
947 ! For now (22Mar2018), standard deviation should be up to 0.75 and cut-off at 3.0
948 ! stddev in order to constrain the various perturbations from being too extreme.
949 !+---+-----------------------------------------------------------------+
950 rand1 = 0.0
951 rand2 = 0.0
952 rand3 = 0.0
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)
960 endif
961 !+---+-----------------------------------------------------------------+
962
963 pptrain = 0.
964 pptsnow = 0.
965 pptgraul = 0.
966 pptice = 0.
967 rainncv(i,j) = 0.
968 IF ( PRESENT (snowncv) ) THEN
969 snowncv(i,j) = 0.
970 ENDIF
971 IF ( PRESENT (icencv) ) THEN
972 icencv(i,j) = 0.
973 ENDIF
974 IF ( PRESENT (graupelncv) ) THEN
975 graupelncv(i,j) = 0.
976 ENDIF
977 sr(i,j) = 0.
978
979 do k = kts, kte
980 if (present(tt)) then
981 t1d(k) = tt(i,k,j)
982 else
983 t1d(k) = th(i,k,j)*pii(i,k,j)
984 end if
985 p1d(k) = p(i,k,j)
986 w1d(k) = w(i,k,j)
987 dz1d(k) = dz(i,k,j)
988 qv1d(k) = qv(i,k,j)
989 qc1d(k) = qc(i,k,j)
990 qi1d(k) = qi(i,k,j)
991 qr1d(k) = qr(i,k,j)
992 qs1d(k) = qs(i,k,j)
993 qg1d(k) = qg(i,k,j)
994 ni1d(k) = ni(i,k,j)
995 nr1d(k) = nr(i,k,j)
996 rho(k) = roverrv * p1d(k) / (r * t1d(k) * (qv1d(k)+roverrv))
997
998 ! These arrays are always allocated and must be initialized
999 !vtsk1(k) = 0.
1000 !txrc1(k) = 0.
1001 !txri1(k) = 0.
1002 initialize_extended_diagnostics: if (ext_diag) then
1003 prw_vcdc1(k) = 0.
1004 prw_vcde1(k) = 0.
1005 tpri_inu1(k) = 0.
1006 tpri_ide1_d(k) = 0.
1007 tpri_ide1_s(k) = 0.
1008 tprs_ide1(k) = 0.
1009 tprs_sde1_d(k) = 0.
1010 tprs_sde1_s(k) = 0.
1011 tprg_gde1_d(k) = 0.
1012 tprg_gde1_s(k) = 0.
1013 tpri_iha1(k) = 0.
1014 tpri_wfz1(k) = 0.
1015 tpri_rfz1(k) = 0.
1016 tprg_rfz1(k) = 0.
1017 tprs_scw1(k) = 0.
1018 tprg_scw1(k) = 0.
1019 tprg_rcs1(k) = 0.
1020 tprs_rcs1(k) = 0.
1021 tprr_rci1(k) = 0.
1022 tprg_rcg1(k) = 0.
1023 tprw_vcd1_c(k) = 0.
1024 tprw_vcd1_e(k) = 0.
1025 tprr_sml1(k) = 0.
1026 tprr_gml1(k) = 0.
1027 tprr_rcg1(k) = 0.
1028 tprr_rcs1(k) = 0.
1029 tprv_rev1(k) = 0.
1030 tten1(k) = 0.
1031 qvten1(k) = 0.
1032 qrten1(k) = 0.
1033 qsten1(k) = 0.
1034 qgten1(k) = 0.
1035 qiten1(k) = 0.
1036 niten1(k) = 0.
1037 nrten1(k) = 0.
1038 ncten1(k) = 0.
1039 qcten1(k) = 0.
1040 endif initialize_extended_diagnostics
1041 enddo
1042 lsml = lsm(i,j)
1043 if (configs%aerosol_aware .or. merra2_aerosol_aware) then
1044 do k = kts, kte
1045 nc1d(k) = nc(i,k,j)
1046 nwfa1d(k) = nwfa(i,k,j)
1047 nifa1d(k) = nifa(i,k,j)
1048 enddo
1049 else
1050 do k = kts, kte
1051 if(lsml == 1) then
1052 nc1d(k) = nt_c_l / rho(k)
1053 else
1054 nc1d(k) = nt_c_o / rho(k)
1055 endif
1056 nwfa1d(k) = nwfa_default
1057 nifa1d(k) = nifa_default
1058 enddo
1059 endif
1060
1061 ! ng and qb are optional hail-aware variables
1062 if ((present(ng)) .and. (present(qb))) then
1063 configs%hail_aware = .true.
1064 do k = kts, kte
1065 ng1d(k) = ng(i,k,j)
1066 qb1d(k) = qb(i,k,j)
1067 enddo
1068 else
1069 do k = kte, kts, -1
1070 ! This is the one-moment graupel formulation
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)
1074 ! zans1 = max(2.0, min(zans1, 6.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)
1081 else
1082 ng1d(k) = 0
1083 qb1d(k) = 0
1084 endif
1085 enddo
1086 endif
1087
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)
1108
1109
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
1117 ! Add ice to snow if separate ice not present
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
1121 ELSE
1122 snowncv(i,j) = pptsnow
1123 snownc(i,j) = snownc(i,j) + pptsnow
1124 ENDIF
1125 ENDIF
1126 ! Use separate ice if present (as in FV3)
1127 IF ( PRESENT(icencv) .AND. PRESENT(icenc) ) THEN
1128 icencv(i,j) = pptice
1129 icenc(i,j) = icenc(i,j) + pptice
1130 ENDIF
1131 IF ( PRESENT(graupelncv) .AND. PRESENT(graupelnc) ) THEN
1132 graupelncv(i,j) = pptgraul
1133 graupelnc(i,j) = graupelnc(i,j) + pptgraul
1134 ENDIF
1135 sr(i,j) = (pptsnow + pptgraul + pptice)/(rainncv(i,j)+1.e-12)
1136
1137
1138
1139 !..Reset lowest model level to initial state aerosols (fake sfc source).
1140 !.. Changed 13 May 2013 to fake emissions in which nwfa2d is aerosol
1141 !.. number tendency (number per kg per second).
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
1147 endif
1148 else
1149 nwfa1d(kts) = nwfa1d(kts) + nwfa2d(i,j)*dt
1150 nifa1d(kts) = nifa1d(kts) + nifa2d(i,j)*dt
1151 end if
1152
1153 do k = kts, kte
1154 nc(i,k,j) = nc1d(k)
1155 nwfa(i,k,j) = nwfa1d(k)
1156 nifa(i,k,j) = nifa1d(k)
1157 enddo
1158 endif
1159
1160 if (merra2_aerosol_aware) then
1161 do k = kts, kte
1162 nc(i,k,j) = nc1d(k)
1163 nwfa(i,k,j) = nwfa1d(k)
1164 nifa(i,k,j) = nifa1d(k)
1165 enddo
1166 endif
1167
1168 if ((present(ng)) .and. (present(qb))) then
1169 do k = kts, kte
1170 ng(i,k,j) = ng1d(k)
1171 qb(i,k,j) = qb1d(k)
1172 enddo
1173 endif
1174
1175 do k = kts, kte
1176 qv(i,k,j) = qv1d(k)
1177 qc(i,k,j) = qc1d(k)
1178 qi(i,k,j) = qi1d(k)
1179 qr(i,k,j) = qr1d(k)
1180 qs(i,k,j) = qs1d(k)
1181 qg(i,k,j) = qg1d(k)
1182 ni(i,k,j) = ni1d(k)
1183 nr(i,k,j) = nr1d(k)
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
1187 tt(i,k,j) = t1d(k)
1188 else
1189 th(i,k,j) = t1d(k)/pii(i,k,j)
1190 endif
1191
1192 if (qc1d(k) .gt. qc_max) then
1193 imax_qc = i
1194 jmax_qc = j
1195 kmax_qc = k
1196 qc_max = qc1d(k)
1197 elseif (qc1d(k) .lt. 0.0) then
1198 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative qc ', qc1d(k), &
1199 ' at i,j,k=', i,j,k
1200 endif
1201 if (qr1d(k) .gt. qr_max) then
1202 imax_qr = i
1203 jmax_qr = j
1204 kmax_qr = k
1205 qr_max = qr1d(k)
1206 elseif (qr1d(k) .lt. 0.0) then
1207 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative qr ', qr1d(k), &
1208 ' at i,j,k=', i,j,k
1209 endif
1210 if (nr1d(k) .gt. nr_max) then
1211 imax_nr = i
1212 jmax_nr = j
1213 kmax_nr = k
1214 nr_max = nr1d(k)
1215 elseif (nr1d(k) .lt. 0.0) then
1216 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative nr ', nr1d(k), &
1217 ' at i,j,k=', i,j,k
1218 endif
1219 if (qs1d(k) .gt. qs_max) then
1220 imax_qs = i
1221 jmax_qs = j
1222 kmax_qs = k
1223 qs_max = qs1d(k)
1224 elseif (qs1d(k) .lt. 0.0) then
1225 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative qs ', qs1d(k), &
1226 ' at i,j,k=', i,j,k
1227 endif
1228 if (qi1d(k) .gt. qi_max) then
1229 imax_qi = i
1230 jmax_qi = j
1231 kmax_qi = k
1232 qi_max = qi1d(k)
1233 elseif (qi1d(k) .lt. 0.0) then
1234 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative qi ', qi1d(k), &
1235 ' at i,j,k=', i,j,k
1236 endif
1237 if (qg1d(k) .gt. qg_max) then
1238 imax_qg = i
1239 jmax_qg = j
1240 kmax_qg = k
1241 qg_max = qg1d(k)
1242 elseif (qg1d(k) .lt. 0.0) then
1243 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative qg ', qg1d(k), &
1244 ' at i,j,k=', i,j,k
1245 endif
1246 if (ni1d(k) .gt. ni_max) then
1247 imax_ni = i
1248 jmax_ni = j
1249 kmax_ni = k
1250 ni_max = ni1d(k)
1251 elseif (ni1d(k) .lt. 0.0) then
1252 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative ni ', ni1d(k), &
1253 ' at i,j,k=', i,j,k
1254 endif
1255 if (qv1d(k) .lt. 0.0) then
1256 write(*,'(a,e16.7,a,3i8)') 'WARNING, negative qv ', qv1d(k), &
1257 ' at i,j,k=', i,j,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)))
1261 else
1262 qv(i,k,j) = 1.e-7
1263 endif
1264 endif
1265 enddo
1266
1267 assign_extended_diagnostics: if (ext_diag) then
1268 do k=kts,kte
1269 !vts1(i,k,j) = vtsk1(k)
1270 !txri(i,k,j) = txri(i,k,j) + txri1(k)
1271 !txrc(i,k,j) = txrc(i,k,j) + txrc1(k)
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)
1309
1310 enddo
1311 endif assign_extended_diagnostics
1312
1313 if (ndt>1 .and. it==ndt) then
1314
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)
1319 ENDIF
1320 IF ( PRESENT (icencv) ) THEN
1321 icencv(i,j) = icenc(i,j)
1322 ENDIF
1323 IF ( PRESENT (graupelncv) ) THEN
1324 graupelncv(i,j) = graupelnc(i,j)
1325 ENDIF
1326 endif
1327
1328 ! Diagnostic calculations only for last step
1329 ! if Thompson MP is called multiple times
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
1333
1334!! max_hail_diam_sfc(i,j) = hail_mass_99th_percentile(kts, kte, qg1d, t1d, p1d, qv1d)
1335
1337
1338 diagflag_present: IF ( PRESENT (diagflag) ) THEN
1339 if (diagflag .and. do_radar_ref == 1) then
1340 !
1341 ! Only set melti to true at the output times
1342 if (fullradar_diag) then
1343 melti=.true.
1344 else
1345 melti=.false.
1346 endif
1347 !
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)
1353 else
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)
1357 end if
1358 do k = kts, kte
1359 refl_10cm(i,k,j) = max(-35., dbz(k))
1360 enddo
1361 endif
1362 ENDIF diagflag_present
1363
1364 IF (has_reqc.ne.0 .and. has_reqi.ne.0 .and. has_reqs.ne.0) THEN
1365 do k = kts, kte
1366 re_qc1d(k) = re_qc_min
1367 re_qi1d(k) = re_qi_min
1368 re_qs1d(k) = re_qs_min
1369 enddo
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)
1375 do k = kts, kte
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))
1379 enddo
1380 ENDIF
1381 ENDIF last_step_only
1382
1383 enddo i_loop
1384 enddo j_loop
1385
1386 ! DEBUG - GT
1387 ! write(*,'(a,7(a,e13.6,1x,a,i3,a,i3,a,i3,a,1x))') 'MP-GT:', &
1388 ! 'qc: ', qc_max, '(', imax_qc, ',', jmax_qc, ',', kmax_qc, ')', &
1389 ! 'qr: ', qr_max, '(', imax_qr, ',', jmax_qr, ',', kmax_qr, ')', &
1390 ! 'qi: ', qi_max, '(', imax_qi, ',', jmax_qi, ',', kmax_qi, ')', &
1391 ! 'qs: ', qs_max, '(', imax_qs, ',', jmax_qs, ',', kmax_qs, ')', &
1392 ! 'qg: ', qg_max, '(', imax_qg, ',', jmax_qg, ',', kmax_qg, ')', &
1393 ! 'ni: ', ni_max, '(', imax_ni, ',', jmax_ni, ',', kmax_ni, ')', &
1394 ! 'nr: ', nr_max, '(', imax_nr, ',', jmax_nr, ',', kmax_nr, ')'
1395 ! END DEBUG - GT
1396 enddo ! end of nt loop
1397
1398 do j = j_start, j_end
1399 do k = kts, kte
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
1403 enddo
1404 enddo
1405 enddo
1406
1407 ! These are always allocated
1408 !deallocate (vtsk1)
1409 !deallocate (txri1)
1410 !deallocate (txrc1)
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)
1439 deallocate (tten1)
1440 deallocate (qvten1)
1441 deallocate (qrten1)
1442 deallocate (qsten1)
1443 deallocate (qgten1)
1444 deallocate (qiten1)
1445 deallocate (niten1)
1446 deallocate (nrten1)
1447 deallocate (ncten1)
1448 deallocate (qcten1)
1449 end if deallocate_extended_diagnostics
1450
1451 END SUBROUTINE tempo_3d_to_1d_driver
1453
1455 SUBROUTINE tempo_finalize()
1456
1457 IMPLICIT NONE
1458
1459 if (ALLOCATED(tcg_racg)) DEALLOCATE(tcg_racg)
1460 if (ALLOCATED(tmr_racg)) DEALLOCATE(tmr_racg)
1461 if (ALLOCATED(tcr_gacr)) DEALLOCATE(tcr_gacr)
1462 if (ALLOCATED(tnr_racg)) DEALLOCATE(tnr_racg)
1463 if (ALLOCATED(tnr_gacr)) DEALLOCATE(tnr_gacr)
1464
1465 if (ALLOCATED(tcs_racs1)) DEALLOCATE(tcs_racs1)
1466 if (ALLOCATED(tmr_racs1)) DEALLOCATE(tmr_racs1)
1467 if (ALLOCATED(tcs_racs2)) DEALLOCATE(tcs_racs2)
1468 if (ALLOCATED(tmr_racs2)) DEALLOCATE(tmr_racs2)
1469 if (ALLOCATED(tcr_sacr1)) DEALLOCATE(tcr_sacr1)
1470 if (ALLOCATED(tms_sacr1)) DEALLOCATE(tms_sacr1)
1471 if (ALLOCATED(tcr_sacr2)) DEALLOCATE(tcr_sacr2)
1472 if (ALLOCATED(tms_sacr2)) DEALLOCATE(tms_sacr2)
1473 if (ALLOCATED(tnr_racs1)) DEALLOCATE(tnr_racs1)
1474 if (ALLOCATED(tnr_racs2)) DEALLOCATE(tnr_racs2)
1475 if (ALLOCATED(tnr_sacr1)) DEALLOCATE(tnr_sacr1)
1476 if (ALLOCATED(tnr_sacr2)) DEALLOCATE(tnr_sacr2)
1477
1478 if (ALLOCATED(tpi_qcfz)) DEALLOCATE(tpi_qcfz)
1479 if (ALLOCATED(tni_qcfz)) DEALLOCATE(tni_qcfz)
1480
1481 if (ALLOCATED(tpi_qrfz)) DEALLOCATE(tpi_qrfz)
1482 if (ALLOCATED(tpg_qrfz)) DEALLOCATE(tpg_qrfz)
1483 if (ALLOCATED(tni_qrfz)) DEALLOCATE(tni_qrfz)
1484 if (ALLOCATED(tnr_qrfz)) DEALLOCATE(tnr_qrfz)
1485
1486 if (ALLOCATED(tps_iaus)) DEALLOCATE(tps_iaus)
1487 if (ALLOCATED(tni_iaus)) DEALLOCATE(tni_iaus)
1488 if (ALLOCATED(tpi_ide)) DEALLOCATE(tpi_ide)
1489
1490 if (ALLOCATED(t_efrw)) DEALLOCATE(t_efrw)
1491 if (ALLOCATED(t_efsw)) DEALLOCATE(t_efsw)
1492
1493 if (ALLOCATED(tnr_rev)) DEALLOCATE(tnr_rev)
1494 if (ALLOCATED(tpc_wev)) DEALLOCATE(tpc_wev)
1495 if (ALLOCATED(tnc_wev)) DEALLOCATE(tnc_wev)
1496
1497 if (ALLOCATED(tnccn_act)) DEALLOCATE(tnccn_act)
1498
1499 END SUBROUTINE tempo_finalize
1500
1501end module module_mp_tempo
1502 !+---+-----------------------------------------------------------------+
1503 !ctrlL
1504 !+---+-----------------------------------------------------------------+
1505 !+---+-----------------------------------------------------------------+
1506
subroutine tempo_finalize()
This module is more library code whereas the individual microphysics schemes contains specific detail...