CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
mp_tempo.F90
1
3
4
7module mp_tempo
8
9 use mpi_f08
10 use machine, only : kind_phys
11
12 use module_mp_tempo_params
13 use module_mp_tempo_utils, only : make_icenumber, make_rainnumber, make_dropletnumber
14 use module_mp_tempo, only : tempo_init, tempo_3d_to_1d_driver, tempo_finalize
15
16 implicit none
17
18 public :: mp_tempo_init, mp_tempo_run, mp_tempo_finalize
19
20 private
21
22 integer, parameter :: ext_ndiag3d = 37
23
24 contains
25
30 subroutine mp_tempo_init(ncol, nlev, con_pi, con_t0c, con_rv, &
31 con_cp, con_rgas, con_boltz, con_amd, &
32 con_amw, con_avgd, con_hvap, con_hfus, &
33 con_g, con_rd, con_eps, &
34 restart, imp_physics, &
35 imp_physics_tempo, convert_dry_rho, &
36 spechum, qc, qr, qi, qs, qg, ni, nr, &
37 chw, vh, &
38 is_aerosol_aware, merra2_aerosol_aware, &
39 is_hail_aware, &
40 nc, nwfa2d, nifa2d, &
41 nwfa, nifa, tgrs, prsl, phil, area, &
42 aerfld, mpicomm, mpirank, mpiroot, &
43 threads, ext_diag, diag3d, &
44 is_initialized, errmsg, errflg)
45
46 implicit none
47
48 ! Interface variables
49 integer, intent(in ) :: ncol
50 integer, intent(in ) :: nlev
51 real(kind_phys), intent(in ) :: con_pi, con_t0c, con_rv, con_cp, con_rgas, &
52 con_boltz, con_amd, con_amw, con_avgd, &
53 con_hvap, con_hfus, con_g, con_rd, con_eps
54 logical, intent(in ) :: restart
55 logical, intent(inout) :: is_initialized
56 integer, intent(in ) :: imp_physics
57 integer, intent(in ) :: imp_physics_tempo
58 ! Hydrometeors
59 logical, intent(in ) :: convert_dry_rho
60 real(kind_phys), intent(inout) :: spechum(:,:)
61 real(kind_phys), intent(inout) :: qc(:,:)
62 real(kind_phys), intent(inout) :: qr(:,:)
63 real(kind_phys), intent(inout) :: qi(:,:)
64 real(kind_phys), intent(inout) :: qs(:,:)
65 real(kind_phys), intent(inout) :: qg(:,:)
66 real(kind_phys), intent(inout) :: ni(:,:)
67 real(kind_phys), intent(inout) :: nr(:,:)
68 real(kind_phys), intent(inout), optional :: chw(:,:), vh(:,:)
69 ! Aerosols
70 logical, intent(in ) :: is_aerosol_aware
71 logical, intent(in ) :: merra2_aerosol_aware
72 logical, intent(in ) :: is_hail_aware
73 real(kind_phys), intent(inout), optional :: nc(:,:)
74 real(kind_phys), intent(inout), optional :: nwfa(:,:)
75 real(kind_phys), intent(inout), optional :: nifa(:,:)
76 real(kind_phys), intent(inout), optional :: nwfa2d(:)
77 real(kind_phys), intent(inout), optional :: nifa2d(:)
78 real(kind_phys), intent(in) :: aerfld(:,:,:)
79 ! State variables
80 real(kind_phys), intent(in ) :: tgrs(:,:)
81 real(kind_phys), intent(in ) :: prsl(:,:)
82 real(kind_phys), intent(in ) :: phil(:,:)
83 real(kind_phys), intent(in ) :: area(:)
84 ! MPI information
85 type(mpi_comm), intent(in ) :: mpicomm
86 integer, intent(in ) :: mpirank
87 integer, intent(in ) :: mpiroot
88 ! Threading/blocking information
89 integer, intent(in ) :: threads
90 ! Extended diagnostics
91 logical, intent(in ) :: ext_diag
92 real(kind_phys), intent(in ), optional :: diag3d(:,:,:)
93 ! CCPP error handling
94 character(len=*), intent( out) :: errmsg
95 integer, intent( out) :: errflg
96
97 !
98 real(kind_phys) :: qv(1:ncol,1:nlev) ! kg kg-1 (water vapor mixing ratio)
99 real(kind_phys) :: hgt(1:ncol,1:nlev) ! m
100 real(kind_phys) :: rho(1:ncol,1:nlev) ! kg m-3
101 real(kind_phys) :: orho(1:ncol,1:nlev) ! m3 kg-1
102 real(kind_phys) :: nc_local(1:ncol,1:nlev) ! needed because nc is only allocated if is_aerosol_aware is true
103 !
104 real (kind=kind_phys) :: h_01, z1, niin3, niccn3
105 integer :: i, k
106
107 ! Initialize the CCPP error handling variables
108 errmsg = ''
109 errflg = 0
110
111 if (is_initialized) return
112
113 ! Consistency checks
114 if (imp_physics/=imp_physics_tempo) then
115 write(errmsg,'(*(a))') "Logic error: namelist choice of microphysics is different from TEMPO MP"
116 errflg = 1
117 return
118 end if
119
120 if (ext_diag) then
121 if (size(diag3d,dim=3) /= ext_ndiag3d) then
122 write(errmsg,'(*(a))') "Logic error: number of diagnostic 3d arrays from model does not match requirements"
123 errflg = 1
124 return
125 end if
126 end if
127
128 if (is_aerosol_aware .and. merra2_aerosol_aware) then
129 write(errmsg,'(*(a))') "Logic error: Only one TEMPO aerosol option can be true, either is_aerosol_aware or merra2_aerosol_aware)"
130 errflg = 1
131 return
132 end if
133
134 ! Call TEMPO init (also sets initial default values of physical constants)
135 if (mpirank==mpiroot) write(*,*) 'Calling tempo_init() with is_aerosol_aware = ', is_aerosol_aware
136
137 call tempo_init(is_aerosol_aware_in=is_aerosol_aware, &
138 merra2_aerosol_aware_in=merra2_aerosol_aware, &
139 is_hail_aware_in=is_hail_aware, &
140 mpicomm=mpicomm, mpirank=mpirank, mpiroot=mpiroot, &
141 threads=threads, errmsg=errmsg, errflg=errflg)
142 if (errflg /= 0) return
143
144 ! For restart runs, the init is done here
145 if (restart) then
146 is_initialized = .true.
147 return
148 end if
149
150 ! Set local TEMPO MP module constants from host model and overwrite derived constants calculated in module_mp_tempo_params/mp_tempo_params_init()
151 pi = con_pi
152 lvap0 = con_hvap
153 lfus = con_hfus
154 lsub = lvap0 + lfus
155 olfus = 1./lfus
156
157 rv = con_rv
158 r = con_rd
159 roverrv = con_eps
160 cp2 = con_cp
161 t_0 = con_t0c
162 r_uni = con_rgas
163 k_b = con_boltz
164 n_avo = con_avgd
165
166 orv = 1.0 / rv
167 am_r = pi * rho_w2 / 6.0
168 am_i = pi * rho_i / 6.0
169 am_g = (/pi*rho_g(1)/6.0, &
170 pi*rho_g(2)/6.0, &
171 pi*rho_g(3)/6.0, &
172 pi*rho_g(4)/6.0, &
173 pi*rho_g(5)/6.0, &
174 pi*rho_g(6)/6.0, &
175 pi*rho_g(7)/6.0, &
176 pi*rho_g(8)/6.0, &
177 pi*rho_g(9)/6.0/)
178
179 m_w = con_amw*1.0e-3 !module_mp_tempo expects kg/mol
180 m_a = con_amd*1.0e-3 !module_mp_tempo expects kg/mol
181 ma_w = m_w/n_avo
182
183 ar_volume = 4.0 / 3.0 * pi * (2.5e-6)**3
184
185 ! Geopotential height in m2 s-2 to height in m
186 hgt = phil/con_g
187
188 ! Ensure non-negative mass mixing ratios of all water variables
189 where(spechum<0) spechum = 1.0e-10 ! COMMENT, gthompsn, spechum should *never* be identically zero.
190 where(qc<0) qc = 0.0
191 where(qr<0) qr = 0.0
192 where(qi<0) qi = 0.0
193 where(qs<0) qs = 0.0
194 where(qg<0) qg = 0.0
195
199 if (merra2_aerosol_aware) then
200 call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
201 end if
202
203
204 qv = spechum/(1.0_kind_phys-spechum)
205
206 if (convert_dry_rho) then
207 qc = qc/(1.0_kind_phys-spechum)
208 qr = qr/(1.0_kind_phys-spechum)
209 qi = qi/(1.0_kind_phys-spechum)
210 qs = qs/(1.0_kind_phys-spechum)
211 qg = qg/(1.0_kind_phys-spechum)
212
213 ni = ni/(1.0_kind_phys-spechum)
214 nr = nr/(1.0_kind_phys-spechum)
215 if (is_hail_aware) then
216 chw = chw/(1.0_kind_phys-spechum)
217 vh = vh/(1.0_kind_phys-spechum)
218 endif
219 if (is_aerosol_aware .or. merra2_aerosol_aware) then
220 nc = nc/(1.0_kind_phys-spechum)
221 nwfa = nwfa/(1.0_kind_phys-spechum)
222 nifa = nifa/(1.0_kind_phys-spechum)
223 end if
224 end if
225
226 ! Density of moist air in kg m-3 and inverse density of air
227 rho = con_eps*prsl/(con_rd*tgrs*(qv+con_eps))
228 orho = 1.0/rho
229
230 ! Ensure we have 1st guess ice number where mass non-zero but no number.
231 where(qi .LE. 0.0) ni=0.0
232 where(qi .GT. 0 .and. ni .LE. 0.0) ni = make_icenumber(qi*rho, tgrs) * orho
233 where(qi .EQ. 0.0 .and. ni .GT. 0.0) ni=0.0
234
235 ! Ensure we have 1st guess rain number where mass non-zero but no number.
236 where(qr .LE. 0.0) nr=0.0
237 where(qr .GT. 0 .and. nr .LE. 0.0) nr = make_rainnumber(qr*rho, tgrs) * orho
238 where(qr .EQ. 0.0 .and. nr .GT. 0.0) nr=0.0
239
240 if (is_hail_aware) then
241 where(qg .LE. 0.0) chw=0.0
242 where(qg .LE. 0.0) vh=0.0
243 endif
244
245 !..Check for existing aerosol data, both CCN and IN aerosols. If missing
246 !.. fill in just a basic vertical profile, somewhat boundary-layer following.
247 if (is_aerosol_aware) then
248
249 ! Potential cloud condensation nuclei (CCN)
250 if (maxval(nwfa) .lt. eps) then
251 if (mpirank==mpiroot) write(*,*) ' Apparently there are no initial CCN aerosols.'
252 do i = 1, ncol
253 if (hgt(i,1).le.1000.0) then
254 h_01 = 0.8
255 elseif (hgt(i,1).ge.2500.0) then
256 h_01 = 0.01
257 else
258 h_01 = 0.8*cos(hgt(i,1)*0.001 - 1.0)
259 endif
260 niccn3 = -1.0*alog(naccn1/naccn0)/h_01
261 nwfa(i,1) = naccn1+naccn0*exp(-((hgt(i,2)-hgt(i,1))/1000.)*niccn3)
262 z1 = hgt(i,2)-hgt(i,1)
263 nwfa2d(i) = nwfa(i,1) * 0.000196 * (50./z1)
264 do k = 2, nlev
265 nwfa(i,k) = naccn1+naccn0*exp(-((hgt(i,k)-hgt(i,1))/1000.)*niccn3)
266 enddo
267 enddo
268 else
269 if (mpirank==mpiroot) write(*,*) ' Apparently initial CCN aerosols are present.'
270 if (maxval(nwfa2d) .lt. eps) then
271 !+---+-----------------------------------------------------------------+
272 !..Scale the lowest level aerosol data into an emissions rate. This is
273 !.. very far from ideal, but need higher emissions where larger amount
274 !.. of (climo) existing and lesser emissions where there exists fewer to
275 !.. begin as a first-order simplistic approach. Later, proper connection to
276 !.. emission inventory would be better, but, for now, scale like this:
277 !.. where: Nwfa=50 per cc, emit 0.875E4 aerosols per second per grid box unit
278 !.. that was tested as ~(20kmx20kmx50m = 2.E10 m**-3)
279 !+---+-----------------------------------------------------------------+
280 if (mpirank==mpiroot) write(*,*) ' Apparently there are no initial CCN aerosol surface emission rates.'
281 do i = 1, ncol
282 z1 = hgt(i,2)-hgt(i,1)
283 nwfa2d(i) = nwfa(i,1) * 0.000196 * (50./z1)
284 enddo
285 else
286 if (mpirank==mpiroot) write(*,*) ' Apparently initial CCN aerosol surface emission rates are present.'
287 endif
288 endif
289
290 ! Potential ice nuclei (IN)
291 if (maxval(nifa) .lt. eps) then
292 if (mpirank==mpiroot) write(*,*) ' Apparently there are no initial IN aerosols.'
293 do i = 1, ncol
294 if (hgt(i,1).le.1000.0) then
295 h_01 = 0.8
296 elseif (hgt(i,1).ge.2500.0) then
297 h_01 = 0.01
298 else
299 h_01 = 0.8*cos(hgt(i,1)*0.001 - 1.0)
300 endif
301 niin3 = -1.0*alog(nain1/nain0)/h_01
302 nifa(i,1) = nain1+nain0*exp(-((hgt(i,2)-hgt(i,1))/1000.)*niin3)
303 nifa2d(i) = 0.
304 do k = 2, nlev
305 nifa(i,k) = nain1+nain0*exp(-((hgt(i,k)-hgt(i,1))/1000.)*niin3)
306 enddo
307 enddo
308 else
309 if (mpirank==mpiroot) write(*,*) ' Apparently initial IN aerosols are present.'
310 if (maxval(nifa2d) .lt. eps) then
311 if (mpirank==mpiroot) write(*,*) ' Apparently there are no initial IN aerosol surface emission rates, set to zero.'
312 ! calculate IN surface flux here, right now just set to zero
313 nifa2d = 0.
314 else
315 if (mpirank==mpiroot) write(*,*) ' Apparently initial IN aerosol surface emission rates are present.'
316 endif
317 endif
318
319 ! Ensure we have 1st guess cloud droplet number where mass non-zero but no number.
320 where(qc .LE. 0.0) nc=0.0
321 where(qc .GT. 0 .and. nc .LE. 0.0) nc = make_dropletnumber(qc*rho, nwfa*rho) * orho
322 where(qc .EQ. 0.0 .and. nc .GT. 0.0) nc = 0.0
323
324 ! Ensure non-negative aerosol number concentrations.
325 where(nwfa .LE. 0.0) nwfa = 1.1e6
326 where(nifa .LE. 0.0) nifa = nain1*0.01
327
328 ! Copy to local array for calculating cloud effective radii below
329 nc_local = nc
330
331 else if (merra2_aerosol_aware) then
332
333 ! Ensure we have 1st guess cloud droplet number where mass non-zero but no number.
334 where(qc .LE. 0.0) nc=0.0
335 where(qc .GT. 0 .and. nc .LE. 0.0) nc = make_dropletnumber(qc*rho, nwfa*rho) * orho
336 where(qc .EQ. 0.0 .and. nc .GT. 0.0) nc = 0.0
337
338 else
339
340 ! Constant droplet concentration for single moment cloud water as in
341 ! module_mp_thompson.F90, only needed for effective radii calculation
342 nc_local = nt_c_l/rho
343
344 end if
345
346 if (convert_dry_rho) then
347 !qc = qc/(1.0_kind_phys+qv)
348 !qr = qr/(1.0_kind_phys+qv)
349 !qi = qi/(1.0_kind_phys+qv)
350 !qs = qs/(1.0_kind_phys+qv)
351 !qg = qg/(1.0_kind_phys+qv)
352
353 ni = ni/(1.0_kind_phys+qv)
354 nr = nr/(1.0_kind_phys+qv)
355 if (is_hail_aware) then
356 chw = chw/(1.0_kind_phys+qv)
357 vh = vh/(1.0_kind_phys+qv)
358 endif
359 if (is_aerosol_aware .or. merra2_aerosol_aware) then
360 nc = nc/(1.0_kind_phys+qv)
361 nwfa = nwfa/(1.0_kind_phys+qv)
362 nifa = nifa/(1.0_kind_phys+qv)
363 end if
364 end if
365
366 is_initialized = .true.
367
368 end subroutine mp_tempo_init
369
370
377 subroutine mp_tempo_run(ncol, nlev, con_g, con_rd, &
378 con_eps, convert_dry_rho, &
379 spechum, qc, qr, qi, qs, qg, ni, nr, &
380 chw, vh, &
381 is_aerosol_aware, is_hail_aware, &
382 merra2_aerosol_aware, nc, nwfa, nifa,&
383 nwfa2d, nifa2d, aero_ind_fdb, &
384 tgrs, prsl, phii, omega, &
385 sedi_semi, decfl, islmsk, dtp, &
386 dt_inner, &
387 first_time_step, istep, nsteps, &
388 prcp, rain, graupel, ice, snow, sr, &
389 refl_10cm, fullradar_diag, &
390 max_hail_diam_sfc, &
391 do_radar_ref, aerfld, &
392 mpicomm, mpirank, mpiroot, blkno, &
393 ext_diag, diag3d, reset_diag3d, &
394 spp_wts_mp, spp_mp, n_var_spp, &
395 spp_prt_list, spp_var_list, &
396 spp_stddev_cutoff, &
397 cplchm, pfi_lsan, pfl_lsan, &
398 is_initialized, errmsg, errflg)
399
400 implicit none
401
402 ! Interface variables
403 logical, intent(inout) :: is_initialized
404 ! Dimensions and constants
405 integer, intent(in ) :: ncol
406 integer, intent(in ) :: nlev
407 real(kind_phys), intent(in ) :: con_g
408 real(kind_phys), intent(in ) :: con_rd
409 real(kind_phys), intent(in ) :: con_eps
410 ! Hydrometeors
411 logical, intent(in ) :: convert_dry_rho
412 real(kind_phys), intent(inout) :: spechum(:,:)
413 real(kind_phys), intent(inout) :: qc(:,:)
414 real(kind_phys), intent(inout) :: qr(:,:)
415 real(kind_phys), intent(inout) :: qi(:,:)
416 real(kind_phys), intent(inout) :: qs(:,:)
417 real(kind_phys), intent(inout) :: qg(:,:)
418 real(kind_phys), intent(inout) :: ni(:,:)
419 real(kind_phys), intent(inout) :: nr(:,:)
420 real(kind_phys), optional, intent(inout) :: chw(:,:), vh(:,:)
421 ! Aerosols
422 logical, intent(in) :: is_aerosol_aware, fullradar_diag
423 logical, intent(in) :: merra2_aerosol_aware, is_hail_aware
424 real(kind_phys), optional, intent(inout) :: nc(:,:)
425 real(kind_phys), optional, intent(inout) :: nwfa(:,:)
426 real(kind_phys), optional, intent(inout) :: nifa(:,:)
427 real(kind_phys), optional, intent(in ) :: nwfa2d(:)
428 real(kind_phys), optional, intent(in ) :: nifa2d(:)
429 real(kind_phys), intent(in) :: aerfld(:,:,:)
430 logical, optional, intent(in ) :: aero_ind_fdb
431 ! State variables and timestep information
432 real(kind_phys), intent(inout) :: tgrs(:,:)
433 real(kind_phys), intent(in ) :: prsl(:,:)
434 real(kind_phys), intent(in ) :: phii(:,:)
435 real(kind_phys), intent(in ) :: omega(:,:)
436 integer, intent(in ) :: islmsk(:)
437 real(kind_phys), intent(in ) :: dtp
438 logical, intent(in ) :: first_time_step
439 integer, intent(in ) :: istep, nsteps
440 real, intent(in ) :: dt_inner
441 ! Precip/rain/snow/graupel fall amounts and fraction of frozen precip
442 real(kind_phys), intent(inout) :: prcp(:)
443 real(kind_phys), intent(inout), optional :: rain(:)
444 real(kind_phys), intent(inout), optional :: graupel(:)
445 real(kind_phys), intent(inout), optional :: ice(:)
446 real(kind_phys), intent(inout), optional :: snow(:)
447 real(kind_phys), intent( out) :: sr(:)
448 ! Radar reflectivity
449 real(kind_phys), intent(inout) :: refl_10cm(:,:)
450 real(kind_phys), intent(inout) :: max_hail_diam_sfc(:)
451 logical, intent(in ) :: do_radar_ref
452 logical, intent(in) :: sedi_semi
453 integer, intent(in) :: decfl
454 ! MPI and block information
455 integer, intent(in) :: blkno
456 type(mpi_comm), intent(in) :: mpicomm
457 integer, intent(in) :: mpirank
458 integer, intent(in) :: mpiroot
459 ! Extended diagnostic output
460 logical, intent(in) :: ext_diag
461 real(kind_phys), target, intent(inout), optional :: diag3d(:,:,:)
462 logical, intent(in) :: reset_diag3d
463
464 ! CCPP error handling
465 character(len=*), intent( out) :: errmsg
466 integer, intent( out) :: errflg
467
468 ! SPP
469 integer, intent(in) :: spp_mp
470 integer, intent(in) :: n_var_spp
471 real(kind_phys), intent(in), optional :: spp_wts_mp(:,:)
472 real(kind_phys), intent(in), optional :: spp_prt_list(:)
473 character(len=10), intent(in), optional :: spp_var_list(:)
474 real(kind_phys), intent(in), optional :: spp_stddev_cutoff(:)
475
476 logical, intent (in) :: cplchm
477 ! ice and liquid water 3d precipitation fluxes - only allocated if cplchm is .true.
478 real(kind=kind_phys), intent(inout), dimension(:,:), optional :: pfi_lsan
479 real(kind=kind_phys), intent(inout), dimension(:,:), optional :: pfl_lsan
480
481 ! Local variables
482
483 ! Reduced time step if subcycling is used
484 real(kind_phys) :: dtstep
485 integer :: ndt
486 ! Air density
487 real(kind_phys) :: rho(1:ncol,1:nlev)
488 ! Water vapor mixing ratio (instead of specific humidity)
489 real(kind_phys) :: qv(1:ncol,1:nlev)
490 ! Vertical velocity and level width
491 real(kind_phys) :: w(1:ncol,1:nlev)
492 real(kind_phys) :: dz(1:ncol,1:nlev)
493 ! Rain/snow/graupel fall amounts
494 real(kind_phys) :: rain_mp(1:ncol) ! mm, dummy, not used
495 real(kind_phys) :: graupel_mp(1:ncol) ! mm, dummy, not used
496 real(kind_phys) :: ice_mp(1:ncol) ! mm, dummy, not used
497 real(kind_phys) :: snow_mp(1:ncol) ! mm, dummy, not used
498 real(kind_phys) :: delta_rain_mp(1:ncol) ! mm
499 real(kind_phys) :: delta_graupel_mp(1:ncol) ! mm
500 real(kind_phys) :: delta_ice_mp(1:ncol) ! mm
501 real(kind_phys) :: delta_snow_mp(1:ncol) ! mm
502
503 real(kind_phys) :: pfils(1:ncol,1:nlev,1)
504 real(kind_phys) :: pflls(1:ncol,1:nlev,1)
505 ! Radar reflectivity
506 logical :: diagflag ! must be true if do_radar_ref is true, not used otherwise
507 integer :: do_radar_ref_mp ! integer instead of logical do_radar_ref
508 ! Effective cloud radii - turned off in CCPP (taken care off in radiation)
509 logical, parameter :: do_effective_radii = .false.
510 integer, parameter :: has_reqc = 0
511 integer, parameter :: has_reqi = 0
512 integer, parameter :: has_reqs = 0
513 integer, parameter :: kme_stoch = 1
514 integer :: spp_mp_opt
515 ! Dimensions used in mp_gt_driver
516 integer :: ids,ide, jds,jde, kds,kde, &
517 ims,ime, jms,jme, kms,kme, &
518 its,ite, jts,jte, kts,kte
519 ! Pointer arrays for extended diagnostics
520 !real(kind_phys), dimension(:,:,:), pointer :: vts1 => null()
521 !real(kind_phys), dimension(:,:,:), pointer :: txri => null()
522 !real(kind_phys), dimension(:,:,:), pointer :: txrc => null()
523 real(kind_phys), dimension(:,:,:), pointer :: prw_vcdc => null()
524 real(kind_phys), dimension(:,:,:), pointer :: prw_vcde => null()
525 real(kind_phys), dimension(:,:,:), pointer :: tpri_inu => null()
526 real(kind_phys), dimension(:,:,:), pointer :: tpri_ide_d => null()
527 real(kind_phys), dimension(:,:,:), pointer :: tpri_ide_s => null()
528 real(kind_phys), dimension(:,:,:), pointer :: tprs_ide => null()
529 real(kind_phys), dimension(:,:,:), pointer :: tprs_sde_d => null()
530 real(kind_phys), dimension(:,:,:), pointer :: tprs_sde_s => null()
531 real(kind_phys), dimension(:,:,:), pointer :: tprg_gde_d => null()
532 real(kind_phys), dimension(:,:,:), pointer :: tprg_gde_s => null()
533 real(kind_phys), dimension(:,:,:), pointer :: tpri_iha => null()
534 real(kind_phys), dimension(:,:,:), pointer :: tpri_wfz => null()
535 real(kind_phys), dimension(:,:,:), pointer :: tpri_rfz => null()
536 real(kind_phys), dimension(:,:,:), pointer :: tprg_rfz => null()
537 real(kind_phys), dimension(:,:,:), pointer :: tprs_scw => null()
538 real(kind_phys), dimension(:,:,:), pointer :: tprg_scw => null()
539 real(kind_phys), dimension(:,:,:), pointer :: tprg_rcs => null()
540 real(kind_phys), dimension(:,:,:), pointer :: tprs_rcs => null()
541 real(kind_phys), dimension(:,:,:), pointer :: tprr_rci => null()
542 real(kind_phys), dimension(:,:,:), pointer :: tprg_rcg => null()
543 real(kind_phys), dimension(:,:,:), pointer :: tprw_vcd_c => null()
544 real(kind_phys), dimension(:,:,:), pointer :: tprw_vcd_e => null()
545 real(kind_phys), dimension(:,:,:), pointer :: tprr_sml => null()
546 real(kind_phys), dimension(:,:,:), pointer :: tprr_gml => null()
547 real(kind_phys), dimension(:,:,:), pointer :: tprr_rcg => null()
548 real(kind_phys), dimension(:,:,:), pointer :: tprr_rcs => null()
549 real(kind_phys), dimension(:,:,:), pointer :: tprv_rev => null()
550 real(kind_phys), dimension(:,:,:), pointer :: tten3 => null()
551 real(kind_phys), dimension(:,:,:), pointer :: qvten3 => null()
552 real(kind_phys), dimension(:,:,:), pointer :: qrten3 => null()
553 real(kind_phys), dimension(:,:,:), pointer :: qsten3 => null()
554 real(kind_phys), dimension(:,:,:), pointer :: qgten3 => null()
555 real(kind_phys), dimension(:,:,:), pointer :: qiten3 => null()
556 real(kind_phys), dimension(:,:,:), pointer :: niten3 => null()
557 real(kind_phys), dimension(:,:,:), pointer :: nrten3 => null()
558 real(kind_phys), dimension(:,:,:), pointer :: ncten3 => null()
559 real(kind_phys), dimension(:,:,:), pointer :: qcten3 => null()
560
561 ! Initialize the CCPP error handling variables
562 errmsg = ''
563 errflg = 0
564
565 if (is_hail_aware .and. sedi_semi) then
566 write(errmsg, fmt='((a))') 'Cannot use hail-aware TEMPO with sedi_semi... plese set sedi_semi=.false.'
567 errflg = 1
568 return
569 endif
570
571 if (first_time_step .and. istep==1 .and. blkno==1) then
572 ! Check initialization state
573 if (.not.is_initialized) then
574 write(errmsg, fmt='((a))') 'mp_tempo_run called before mp_tempo_init'
575 errflg = 1
576 return
577 end if
578 ! Check forr optional arguments of aerosol-aware microphysics
579 if (is_aerosol_aware .and. .not. (present(nc) .and. &
580 present(nwfa) .and. &
581 present(nifa) .and. &
582 present(nwfa2d) .and. &
583 present(nifa2d) )) then
584 write(errmsg,fmt='(*(a))') 'Logic error in mp_tempo_run:', &
585 ' aerosol-aware microphysics require all of the', &
586 ' following optional arguments:', &
587 ' nc, nwfa, nifa, nwfa2d, nifa2d'
588 errflg = 1
589 return
590 else if (merra2_aerosol_aware .and. .not. (present(nc) .and. &
591 present(nwfa) .and. &
592 present(nifa) )) then
593 write(errmsg,fmt='(*(a))') 'Logic error in mp_tempo_run:', &
594 ' merra2 aerosol-aware microphysics require the', &
595 ' following optional arguments: nc, nwfa, nifa'
596 errflg = 1
597 return
598 end if
599 ! Consistency cheecks - subcycling and inner loop at the same time are not supported
600 if (nsteps>1 .and. dt_inner < dtp) then
601 write(errmsg,'(*(a))') "Logic error: Subcycling and inner loop cannot be used at the same time"
602 errflg = 1
603 return
604 else if (mpirank==mpiroot .and. nsteps>1) then
605 write(*,'(a,i0,a,a,f6.2,a)') 'TEMPO MP is using ', nsteps, ' substep(s) per time step with an ', &
606 'effective time step of ', dtp/real(nsteps, kind=kind_phys), ' seconds'
607 else if (mpirank==mpiroot .and. dt_inner < dtp) then
608 ndt = max(nint(dtp/dt_inner),1)
609 write(*,'(a,i0,a,a,f6.2,a)') 'TEMPO MP is using ', ndt, ' inner loops per time step with an ', &
610 'effective time step of ', dtp/real(ndt, kind=kind_phys), ' seconds'
611 end if
612 end if
613
614 ! Set stochastic physics selection to apply all perturbations
615 if ( spp_mp==7 ) then
616 spp_mp_opt=7
617 else
618 spp_mp_opt=0
619 endif
620
621 ! Set reduced time step if subcycling is used
622 if (nsteps>1) then
623 dtstep = dtp/real(nsteps, kind=kind_phys)
624 else
625 dtstep = dtp
626 end if
627 if (merra2_aerosol_aware) then
628 call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
629 end if
630
634
635 ! DH* - do this only if istep == 1? Would be ok if it was
636 ! guaranteed that nothing else in the same subcycle group
637 ! was using these arrays, but it is somewhat dangerous.
638 qv = spechum/(1.0_kind_phys-spechum)
639
640 if (convert_dry_rho) then
641 qc = qc/(1.0_kind_phys-spechum)
642 qr = qr/(1.0_kind_phys-spechum)
643 qi = qi/(1.0_kind_phys-spechum)
644 qs = qs/(1.0_kind_phys-spechum)
645 qg = qg/(1.0_kind_phys-spechum)
646
647 ni = ni/(1.0_kind_phys-spechum)
648 nr = nr/(1.0_kind_phys-spechum)
649 if (is_hail_aware) then
650 chw = chw/(1.0_kind_phys-spechum)
651 vh = vh/(1.0_kind_phys-spechum)
652 endif
653 if (is_aerosol_aware .or. merra2_aerosol_aware) then
654 nc = nc/(1.0_kind_phys-spechum)
655 nwfa = nwfa/(1.0_kind_phys-spechum)
656 nifa = nifa/(1.0_kind_phys-spechum)
657 end if
658 end if
659 ! *DH
660
662 rho = con_eps*prsl/(con_rd*tgrs*(qv+con_eps))
663
665 w = -omega/(rho*con_g)
666
668 dz = (phii(:,2:nlev+1) - phii(:,1:nlev)) / con_g
669
670 ! Accumulated values inside Thompson scheme, not used;
671 ! only use delta and add to inout variables (different units)
672 rain_mp = 0
673 graupel_mp = 0
674 ice_mp = 0
675 snow_mp = 0
676 delta_rain_mp = 0
677 delta_graupel_mp = 0
678 delta_ice_mp = 0
679 delta_snow_mp = 0
680
681 ! Flags for calculating radar reflectivity; diagflag is redundant
682 if (do_radar_ref) then
683 diagflag = .true.
684 do_radar_ref_mp = 1
685 else
686 diagflag = .false.
687 do_radar_ref_mp = 0
688 end if
689
690 ! Set internal dimensions
691 ids = 1
692 ims = 1
693 its = 1
694 ide = ncol
695 ime = ncol
696 ite = ncol
697 jds = 1
698 jms = 1
699 jts = 1
700 jde = 1
701 jme = 1
702 jte = 1
703 kds = 1
704 kms = 1
705 kts = 1
706 kde = nlev
707 kme = nlev
708 kte = nlev
709 if(cplchm) then
710 pfi_lsan = 0.0
711 pfl_lsan = 0.0
712 end if
713
714 ! Set pointers for extended diagnostics
715 set_extended_diagnostic_pointers: if (ext_diag) then
716 if (reset_diag3d) then
717 diag3d = 0.0
718 end if
719 !vts1 => diag3d(:,:,X:X)
720 !txri => diag3d(:,:,X:X)
721 !txrc => diag3d(:,:,X:X)
722 prw_vcdc => diag3d(:,:,1:1)
723 prw_vcde => diag3d(:,:,2:2)
724 tpri_inu => diag3d(:,:,3:3)
725 tpri_ide_d => diag3d(:,:,4:4)
726 tpri_ide_s => diag3d(:,:,5:5)
727 tprs_ide => diag3d(:,:,6:6)
728 tprs_sde_d => diag3d(:,:,7:7)
729 tprs_sde_s => diag3d(:,:,8:8)
730 tprg_gde_d => diag3d(:,:,9:9)
731 tprg_gde_s => diag3d(:,:,10:10)
732 tpri_iha => diag3d(:,:,11:11)
733 tpri_wfz => diag3d(:,:,12:12)
734 tpri_rfz => diag3d(:,:,13:13)
735 tprg_rfz => diag3d(:,:,14:14)
736 tprs_scw => diag3d(:,:,15:15)
737 tprg_scw => diag3d(:,:,16:16)
738 tprg_rcs => diag3d(:,:,17:17)
739 tprs_rcs => diag3d(:,:,18:18)
740 tprr_rci => diag3d(:,:,19:19)
741 tprg_rcg => diag3d(:,:,20:20)
742 tprw_vcd_c => diag3d(:,:,21:21)
743 tprw_vcd_e => diag3d(:,:,22:22)
744 tprr_sml => diag3d(:,:,23:23)
745 tprr_gml => diag3d(:,:,24:24)
746 tprr_rcg => diag3d(:,:,25:25)
747 tprr_rcs => diag3d(:,:,26:26)
748 tprv_rev => diag3d(:,:,27:27)
749 tten3 => diag3d(:,:,28:28)
750 qvten3 => diag3d(:,:,29:29)
751 qrten3 => diag3d(:,:,30:30)
752 qsten3 => diag3d(:,:,31:31)
753 qgten3 => diag3d(:,:,32:32)
754 qiten3 => diag3d(:,:,33:33)
755 niten3 => diag3d(:,:,34:34)
756 nrten3 => diag3d(:,:,35:35)
757 ncten3 => diag3d(:,:,36:36)
758 qcten3 => diag3d(:,:,37:37)
759 else
760 allocate(prw_vcdc(0,0,0))
761 allocate(prw_vcde(0,0,0))
762 allocate(tpri_inu(0,0,0))
763 allocate(tpri_ide_d(0,0,0))
764 allocate(tpri_ide_s(0,0,0))
765 allocate(tprs_ide(0,0,0))
766 allocate(tprs_sde_d(0,0,0))
767 allocate(tprs_sde_s(0,0,0))
768 allocate(tprg_gde_d(0,0,0))
769 allocate(tprg_gde_s(0,0,0))
770 allocate(tpri_iha(0,0,0))
771 allocate(tpri_wfz(0,0,0))
772 allocate(tpri_rfz(0,0,0))
773 allocate(tprg_rfz(0,0,0))
774 allocate(tprs_scw(0,0,0))
775 allocate(tprg_scw(0,0,0))
776 allocate(tprg_rcs(0,0,0))
777 allocate(tprs_rcs(0,0,0))
778 allocate(tprr_rci(0,0,0))
779 allocate(tprg_rcg(0,0,0))
780 allocate(tprw_vcd_c(0,0,0))
781 allocate(tprw_vcd_e(0,0,0))
782 allocate(tprr_sml(0,0,0))
783 allocate(tprr_gml(0,0,0))
784 allocate(tprr_rcg(0,0,0))
785 allocate(tprr_rcs(0,0,0))
786 allocate(tprv_rev(0,0,0))
787 allocate(tten3(0,0,0))
788 allocate(qvten3(0,0,0))
789 allocate(qrten3(0,0,0))
790 allocate(qsten3(0,0,0))
791 allocate(qgten3(0,0,0))
792 allocate(qiten3(0,0,0))
793 allocate(niten3(0,0,0))
794 allocate(nrten3(0,0,0))
795 allocate(ncten3(0,0,0))
796 allocate(qcten3(0,0,0))
797 end if set_extended_diagnostic_pointers
799 if (is_aerosol_aware) then
800 if (is_hail_aware) then
801 call tempo_3d_to_1d_driver(qv=qv, qc=qc, qr=qr, qi=qi, qs=qs, qg=qg, qb=vh, ni=ni, nr=nr, &
802 nc=nc, ng=chw, nwfa=nwfa, nifa=nifa, nwfa2d=nwfa2d, nifa2d=nifa2d, &
803 tt=tgrs, p=prsl, w=w, dz=dz, dt_in=dtstep, dt_inner=dt_inner, &
804 sedi_semi=sedi_semi, decfl=decfl, lsm=islmsk, &
805 rainnc=rain_mp, rainncv=delta_rain_mp, &
806 snownc=snow_mp, snowncv=delta_snow_mp, &
807 icenc=ice_mp, icencv=delta_ice_mp, &
808 graupelnc=graupel_mp, graupelncv=delta_graupel_mp, sr=sr, &
809 refl_10cm=refl_10cm, &
810 diagflag=diagflag, do_radar_ref=do_radar_ref_mp, &
811 max_hail_diam_sfc=max_hail_diam_sfc, &
812 has_reqc=has_reqc, has_reqi=has_reqi, has_reqs=has_reqs, &
813 aero_ind_fdb=aero_ind_fdb, rand_perturb_on=spp_mp_opt, &
814 kme_stoch=kme_stoch, &
815 rand_pert=spp_wts_mp, spp_var_list=spp_var_list, &
816 spp_prt_list=spp_prt_list, n_var_spp=n_var_spp, &
817 spp_stddev_cutoff=spp_stddev_cutoff, &
818 ids=ids, ide=ide, jds=jds, jde=jde, kds=kds, kde=kde, &
819 ims=ims, ime=ime, jms=jms, jme=jme, kms=kms, kme=kme, &
820 its=its, ite=ite, jts=jts, jte=jte, kts=kts, kte=kte, &
821 fullradar_diag=fullradar_diag, istep=istep, nsteps=nsteps, &
822 first_time_step=first_time_step, errmsg=errmsg, errflg=errflg, &
823 ! Extended diagnostics
824 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
825 else
826
827! write(errmsg,'(*(a))') "TEMPO aerosol-aware UNTESTED -- DO NOT USE"
828! errflg = 1
829! return
830 call tempo_3d_to_1d_driver(qv=qv, qc=qc, qr=qr, qi=qi, qs=qs, qg=qg, ni=ni, nr=nr, &
831 nc=nc, nwfa=nwfa, nifa=nifa, nwfa2d=nwfa2d, nifa2d=nifa2d, &
832 tt=tgrs, p=prsl, w=w, dz=dz, dt_in=dtstep, dt_inner=dt_inner, &
833 sedi_semi=sedi_semi, decfl=decfl, lsm=islmsk, &
834 rainnc=rain_mp, rainncv=delta_rain_mp, &
835 snownc=snow_mp, snowncv=delta_snow_mp, &
836 icenc=ice_mp, icencv=delta_ice_mp, &
837 graupelnc=graupel_mp, graupelncv=delta_graupel_mp, sr=sr, &
838 refl_10cm=refl_10cm, &
839 diagflag=diagflag, do_radar_ref=do_radar_ref_mp, &
840 max_hail_diam_sfc=max_hail_diam_sfc, &
841 has_reqc=has_reqc, has_reqi=has_reqi, has_reqs=has_reqs, &
842 aero_ind_fdb=aero_ind_fdb, rand_perturb_on=spp_mp_opt, &
843 kme_stoch=kme_stoch, &
844 rand_pert=spp_wts_mp, spp_var_list=spp_var_list, &
845 spp_prt_list=spp_prt_list, n_var_spp=n_var_spp, &
846 spp_stddev_cutoff=spp_stddev_cutoff, &
847 ids=ids, ide=ide, jds=jds, jde=jde, kds=kds, kde=kde, &
848 ims=ims, ime=ime, jms=jms, jme=jme, kms=kms, kme=kme, &
849 its=its, ite=ite, jts=jts, jte=jte, kts=kts, kte=kte, &
850 fullradar_diag=fullradar_diag, istep=istep, nsteps=nsteps, &
851 first_time_step=first_time_step, errmsg=errmsg, errflg=errflg, &
852 ! Extended diagnostics
853 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
854 ! ! vts1=vts1, txri=txri, txrc=txrc, &
855 ! prw_vcdc=prw_vcdc, &
856 ! prw_vcde=prw_vcde, tpri_inu=tpri_inu, tpri_ide_d=tpri_ide_d, &
857 ! tpri_ide_s=tpri_ide_s, tprs_ide=tprs_ide, &
858 ! tprs_sde_d=tprs_sde_d, &
859 ! tprs_sde_s=tprs_sde_s, tprg_gde_d=tprg_gde_d, &
860 ! tprg_gde_s=tprg_gde_s, tpri_iha=tpri_iha, &
861 ! tpri_wfz=tpri_wfz, tpri_rfz=tpri_rfz, tprg_rfz=tprg_rfz, &
862 ! tprs_scw=tprs_scw, tprg_scw=tprg_scw, tprg_rcs=tprg_rcs, &
863 ! tprs_rcs=tprs_rcs, &
864 ! tprr_rci=tprr_rci, tprg_rcg=tprg_rcg, tprw_vcd_c=tprw_vcd_c, &
865 ! tprw_vcd_e=tprw_vcd_e, tprr_sml=tprr_sml, tprr_gml=tprr_gml, &
866 ! tprr_rcg=tprr_rcg, tprr_rcs=tprr_rcs, &
867 ! tprv_rev=tprv_rev, tten3=tten3, &
868 ! qvten3=qvten3, qrten3=qrten3, qsten3=qsten3, qgten3=qgten3, &
869 ! qiten3=qiten3, niten3=niten3, nrten3=nrten3, ncten3=ncten3, &
870 ! qcten3=qcten3,
871 endif
872 else if (merra2_aerosol_aware) then
873 write(errmsg,'(*(a))') "TEMPO aerosol-aware with MERRA2 UNTESTED -- DO NOT USE"
874 errflg = 1
875 return
876 call tempo_3d_to_1d_driver(qv=qv, qc=qc, qr=qr, qi=qi, qs=qs, qg=qg, ni=ni, nr=nr, &
877 nc=nc, nwfa=nwfa, nifa=nifa, &
878 tt=tgrs, p=prsl, w=w, dz=dz, dt_in=dtstep, dt_inner=dt_inner, &
879 sedi_semi=sedi_semi, decfl=decfl, lsm=islmsk, &
880 rainnc=rain_mp, rainncv=delta_rain_mp, &
881 snownc=snow_mp, snowncv=delta_snow_mp, &
882 icenc=ice_mp, icencv=delta_ice_mp, &
883 graupelnc=graupel_mp, graupelncv=delta_graupel_mp, sr=sr, &
884 refl_10cm=refl_10cm, &
885 diagflag=diagflag, do_radar_ref=do_radar_ref_mp, &
886 max_hail_diam_sfc=max_hail_diam_sfc, &
887 has_reqc=has_reqc, has_reqi=has_reqi, has_reqs=has_reqs, &
888 aero_ind_fdb=aero_ind_fdb, rand_perturb_on=spp_mp_opt, &
889 kme_stoch=kme_stoch, &
890 rand_pert=spp_wts_mp, spp_var_list=spp_var_list, &
891 spp_prt_list=spp_prt_list, n_var_spp=n_var_spp, &
892 spp_stddev_cutoff=spp_stddev_cutoff, &
893 ids=ids, ide=ide, jds=jds, jde=jde, kds=kds, kde=kde, &
894 ims=ims, ime=ime, jms=jms, jme=jme, kms=kms, kme=kme, &
895 its=its, ite=ite, jts=jts, jte=jte, kts=kts, kte=kte, &
896 fullradar_diag=fullradar_diag, istep=istep, nsteps=nsteps, &
897 first_time_step=first_time_step, errmsg=errmsg, errflg=errflg, &
898 ! Extended diagnostics
899 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
900 ! ! vts1=vts1, txri=txri, txrc=txrc, &
901 ! prw_vcdc=prw_vcdc, &
902 ! prw_vcde=prw_vcde, tpri_inu=tpri_inu, tpri_ide_d=tpri_ide_d, &
903 ! tpri_ide_s=tpri_ide_s, tprs_ide=tprs_ide, &
904 ! tprs_sde_d=tprs_sde_d, &
905 ! tprs_sde_s=tprs_sde_s, tprg_gde_d=tprg_gde_d, &
906 ! tprg_gde_s=tprg_gde_s, tpri_iha=tpri_iha, &
907 ! tpri_wfz=tpri_wfz, tpri_rfz=tpri_rfz, tprg_rfz=tprg_rfz, &
908 ! tprs_scw=tprs_scw, tprg_scw=tprg_scw, tprg_rcs=tprg_rcs, &
909 ! tprs_rcs=tprs_rcs, &
910 ! tprr_rci=tprr_rci, tprg_rcg=tprg_rcg, tprw_vcd_c=tprw_vcd_c, &
911 ! tprw_vcd_e=tprw_vcd_e, tprr_sml=tprr_sml, tprr_gml=tprr_gml, &
912 ! tprr_rcg=tprr_rcg, tprr_rcs=tprr_rcs, &
913 ! tprv_rev=tprv_rev, tten3=tten3, &
914 ! qvten3=qvten3, qrten3=qrten3, qsten3=qsten3, qgten3=qgten3, &
915 ! qiten3=qiten3, niten3=niten3, nrten3=nrten3, ncten3=ncten3, &
916 ! qcten3=qcten3,
917 else
918 call tempo_3d_to_1d_driver(qv=qv, qc=qc, qr=qr, qi=qi, qs=qs, qg=qg, ni=ni, nr=nr, &
919 tt=tgrs, p=prsl, w=w, dz=dz, dt_in=dtstep, dt_inner=dt_inner, &
920 sedi_semi=sedi_semi, decfl=decfl, lsm=islmsk, &
921 rainnc=rain_mp, rainncv=delta_rain_mp, &
922 snownc=snow_mp, snowncv=delta_snow_mp, &
923 icenc=ice_mp, icencv=delta_ice_mp, &
924 graupelnc=graupel_mp, graupelncv=delta_graupel_mp, sr=sr, &
925 refl_10cm=refl_10cm, &
926 diagflag=diagflag, do_radar_ref=do_radar_ref_mp, &
927 max_hail_diam_sfc=max_hail_diam_sfc, &
928 has_reqc=has_reqc, has_reqi=has_reqi, has_reqs=has_reqs, &
929 rand_perturb_on=spp_mp_opt, kme_stoch=kme_stoch, &
930 rand_pert=spp_wts_mp, spp_var_list=spp_var_list, &
931 spp_prt_list=spp_prt_list, n_var_spp=n_var_spp, &
932 spp_stddev_cutoff=spp_stddev_cutoff, &
933 ids=ids, ide=ide, jds=jds, jde=jde, kds=kds, kde=kde, &
934 ims=ims, ime=ime, jms=jms, jme=jme, kms=kms, kme=kme, &
935 its=its, ite=ite, jts=jts, jte=jte, kts=kts, kte=kte, &
936 fullradar_diag=fullradar_diag, istep=istep, nsteps=nsteps, &
937 first_time_step=first_time_step, errmsg=errmsg, errflg=errflg, &
938 ! Extended diagnostics
939 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
940 !! vts1=vts1, txri=txri, txrc=txrc, &
941 ! prw_vcdc=prw_vcdc, &
942 ! prw_vcde=prw_vcde, tpri_inu=tpri_inu, tpri_ide_d=tpri_ide_d, &
943 ! tpri_ide_s=tpri_ide_s, tprs_ide=tprs_ide, &
944 ! tprs_sde_d=tprs_sde_d, &
945 ! tprs_sde_s=tprs_sde_s, tprg_gde_d=tprg_gde_d, &
946 ! tprg_gde_s=tprg_gde_s, tpri_iha=tpri_iha, &
947 ! tpri_wfz=tpri_wfz, tpri_rfz=tpri_rfz, tprg_rfz=tprg_rfz, &
948 ! tprs_scw=tprs_scw, tprg_scw=tprg_scw, tprg_rcs=tprg_rcs, &
949 ! tprs_rcs=tprs_rcs, &
950 ! tprr_rci=tprr_rci, tprg_rcg=tprg_rcg, tprw_vcd_c=tprw_vcd_c, &
951 ! tprw_vcd_e=tprw_vcd_e, tprr_sml=tprr_sml, tprr_gml=tprr_gml, &
952 ! tprr_rcg=tprr_rcg, tprr_rcs=tprr_rcs, &
953 ! tprv_rev=tprv_rev, tten3=tten3, &
954 ! qvten3=qvten3, qrten3=qrten3, qsten3=qsten3, qgten3=qgten3, &
955 ! qiten3=qiten3, niten3=niten3, nrten3=nrten3, ncten3=ncten3, &
956 ! qcten3=qcten3)
957 end if
958 if (errflg/=0) return
959
960 ! DH* - do this only if istep == nsteps? Would be ok if it was
961 ! guaranteed that nothing else in the same subcycle group
962 ! was using these arrays, but it is somewhat dangerous.
963
965 spechum = qv/(1.0_kind_phys+qv)
966
967 if (convert_dry_rho) then
968 qc = qc/(1.0_kind_phys+qv)
969 qr = qr/(1.0_kind_phys+qv)
970 qi = qi/(1.0_kind_phys+qv)
971 qs = qs/(1.0_kind_phys+qv)
972 qg = qg/(1.0_kind_phys+qv)
973
974 ni = ni/(1.0_kind_phys+qv)
975 nr = nr/(1.0_kind_phys+qv)
976 if (is_hail_aware) then
977 chw = chw/(1.0_kind_phys+qv)
978 vh = vh/(1.0_kind_phys+qv)
979 endif
980 if (is_aerosol_aware .or. merra2_aerosol_aware) then
981 nc = nc/(1.0_kind_phys+qv)
982 nwfa = nwfa/(1.0_kind_phys+qv)
983 nifa = nifa/(1.0_kind_phys+qv)
984 end if
985 end if
986 ! *DH
987
989 ! "rain" in Thompson MP refers to precipitation (total of liquid rainfall+snow+graupel+ice)
990 prcp = prcp + max(0.0, delta_rain_mp/1000.0_kind_phys)
991 graupel = graupel + max(0.0, delta_graupel_mp/1000.0_kind_phys)
992 ice = ice + max(0.0, delta_ice_mp/1000.0_kind_phys)
993 snow = snow + max(0.0, delta_snow_mp/1000.0_kind_phys)
994 rain = rain + max(0.0, (delta_rain_mp - (delta_graupel_mp + delta_ice_mp + delta_snow_mp))/1000.0_kind_phys)
995
996 ! Recompute sr at last subcycling step
997 if (nsteps>1 .and. istep == nsteps) then
998 ! Unlike inside mp_gt_driver, rain does not contain frozen precip
999 sr = (snow + graupel + ice)/(rain + snow + graupel + ice +1.e-12)
1000 end if
1001
1002 ! output instantaneous ice/snow and rain water 3d precipitation fluxes
1003 if(cplchm) then
1004 pfi_lsan(:,:) = pfils(:,:,1)
1005 pfl_lsan(:,:) = pflls(:,:,1)
1006 end if
1007
1008 ! DH* Not really needed because they go out of scope ...
1009 ! But having them in here seems to cause problems with Intel?
1010 ! It looked like this is also nullifying the pointers passed
1011 ! from the CCPP caps.
1012 !unset_extended_diagnostic_pointers: if (ext_diag) then
1013 ! !vts1 => null()
1014 ! !txri => null()
1015 ! !txrc => null()
1016 ! prw_vcdc => null()
1017 ! prw_vcde => null()
1018 ! tpri_inu => null()
1019 ! tpri_ide_d => null()
1020 ! tpri_ide_s => null()
1021 ! tprs_ide => null()
1022 ! tprs_sde_d => null()
1023 ! tprs_sde_s => null()
1024 ! tprg_gde_d => null()
1025 ! tprg_gde_s => null()
1026 ! tpri_iha => null()
1027 ! tpri_wfz => null()
1028 ! tpri_rfz => null()
1029 ! tprg_rfz => null()
1030 ! tprs_scw => null()
1031 ! tprg_scw => null()
1032 ! tprg_rcs => null()
1033 ! tprs_rcs => null()
1034 ! tprr_rci => null()
1035 ! tprg_rcg => null()
1036 ! tprw_vcd_c => null()
1037 ! tprw_vcd_e => null()
1038 ! tprr_sml => null()
1039 ! tprr_gml => null()
1040 ! tprr_rcg => null()
1041 ! tprr_rcs => null()
1042 ! tprv_rev => null()
1043 ! tten3 => null()
1044 ! qvten3 => null()
1045 ! qrten3 => null()
1046 ! qsten3 => null()
1047 ! qgten3 => null()
1048 ! qiten3 => null()
1049 ! niten3 => null()
1050 ! nrten3 => null()
1051 ! ncten3 => null()
1052 ! qcten3 => null()
1053 !end if unset_extended_diagnostic_pointers
1054 ! *DH
1055
1056 end subroutine mp_tempo_run
1058
1062 subroutine mp_tempo_finalize(is_initialized, errmsg, errflg)
1063
1064 implicit none
1065 logical, intent(inout) :: is_initialized
1066 character(len=*), intent( out) :: errmsg
1067 integer, intent( out) :: errflg
1068
1069 ! Initialize the CCPP error handling variables
1070 errmsg = ''
1071 errflg = 0
1072
1073 if (.not.is_initialized) return
1074
1075 call tempo_finalize()
1076
1077 is_initialized = .false.
1078
1079 end subroutine mp_tempo_finalize
1080
1081 subroutine get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
1082 ! To calculate nifa and nwfa from bins of aerosols.
1083 ! In GOCART and MERRA2, aerosols are given as mixing ratio (kg/kg). To
1084 ! convert from kg/kg to #/kg, the "unit mass" (mass of one particle)
1085 ! within the mass bins is calculated. A lognormal size distribution
1086 ! within aerosol bins is used to find the size based upon the median
1087 ! mass. NIFA is mainly summarized over five dust bins and NWFA over the
1088 ! other 10 bins. The parameters besides each bins are carefully tuned
1089 ! for a good performance of the scheme.
1090 !
1091 ! The fields for the last index of the aerfld array
1092 ! are specified as below.
1093 ! 1: dust bin 1, 0.1 to 1.0 micrometers
1094 ! 2: dust bin 2, 1.0 to 1.8 micrometers
1095 ! 3: dust bin 3, 1.8 to 3.0 micrometers
1096 ! 4: dust bin 4, 3.0 to 6.0 micrometers
1097 ! 5: dust bin 5, 6.0 to 10.0 micrometers
1098 ! 6: sea salt bin 1, 0.03 to 0.1 micrometers
1099 ! 7: sea salt bin 2, 0.1 to 0.5 micrometers
1100 ! 8: sea salt bin 3, 0.5 to 1.5 micrometers
1101 ! 9: sea salt bin 4, 1.5 to 5.0 micrometers
1102 ! 10: sea salt bin 5, 5.0 to 10.0 micrometers
1103 ! 11: Sulfate, 0.35 (mean) micrometers
1104 ! 15: water-friendly organic carbon, 0.35 (mean) micrometers
1105 !
1106 ! Bin densities are as follows:
1107 ! 1: dust bin 1: 2500 kg/m2
1108 ! 2-5: dust bin 2-5: 2650 kg/m2
1109 ! 6-10: sea salt bins 6-10: 2200 kg/m2
1110 ! 11: sulfate: 1700 kg/m2
1111 ! 15: organic carbon: 1800 kg/m2
1112
1113 implicit none
1114 integer, intent(in)::ncol, nlev
1115 real (kind=kind_phys), dimension(:,:,:), intent(in) :: aerfld
1116 real (kind=kind_phys), dimension(:,:), intent(out ):: nifa, nwfa
1117
1118 nifa=(aerfld(:,:,1)/4.0737762+aerfld(:,:,2)/30.459203+aerfld(:,:,3)/153.45048+ &
1119 aerfld(:,:,4)/1011.5142+ aerfld(:,:,5)/5683.3501)*1.e15
1120
1121 nwfa=((aerfld(:,:,6)/0.0045435214+aerfld(:,:,7)/0.2907854+aerfld(:,:,8)/12.91224+ &
1122 aerfld(:,:,9)/206.2216+ aerfld(:,:,10)/4326.23)*9.+aerfld(:,:,11)/0.3053104*5+ &
1123 aerfld(:,:,15)/0.3232698*8)*1.e15
1124 end subroutine get_niwfa
1125
1126end module mp_tempo
subroutine, public mp_tempo_run(ncol, nlev, con_g, con_rd, con_eps, convert_dry_rho, spechum, qc, qr, qi, qs, qg, ni, nr, chw, vh, is_aerosol_aware, is_hail_aware, merra2_aerosol_aware, nc, nwfa, nifa, nwfa2d, nifa2d, aero_ind_fdb, tgrs, prsl, phii, omega, sedi_semi, decfl, islmsk, dtp, dt_inner, first_time_step, istep, nsteps, prcp, rain, graupel, ice, snow, sr, refl_10cm, fullradar_diag, max_hail_diam_sfc, do_radar_ref, aerfld, mpicomm, mpirank, mpiroot, blkno, ext_diag, diag3d, reset_diag3d, spp_wts_mp, spp_mp, n_var_spp, spp_prt_list, spp_var_list, spp_stddev_cutoff, cplchm, pfi_lsan, pfl_lsan, is_initialized, errmsg, errflg)
Definition mp_tempo.F90:399
subroutine tempo_finalize()