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, &
38 is_aerosol_aware, merra2_aerosol_aware, &
41 nwfa, nifa, tgrs, prsl, phil, area, &
42 aerfld, mpicomm, mpirank, mpiroot, &
43 threads, ext_diag, diag3d, &
44 is_initialized, errmsg, errflg)
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
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(:,:)
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(:,:,:)
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(:)
85 type(mpi_comm),
intent(in ) :: mpicomm
86 integer,
intent(in ) :: mpirank
87 integer,
intent(in ) :: mpiroot
89 integer,
intent(in ) :: threads
91 logical,
intent(in ) :: ext_diag
92 real(kind_phys),
intent(in ),
optional :: diag3d(:,:,:)
94 character(len=*),
intent( out) :: errmsg
95 integer,
intent( out) :: errflg
98 real(kind_phys) :: qv(1:ncol,1:nlev)
99 real(kind_phys) :: hgt(1:ncol,1:nlev)
100 real(kind_phys) :: rho(1:ncol,1:nlev)
101 real(kind_phys) :: orho(1:ncol,1:nlev)
102 real(kind_phys) :: nc_local(1:ncol,1:nlev)
104 real (kind=kind_phys) :: h_01, z1, niin3, niccn3
111 if (is_initialized)
return
114 if (imp_physics/=imp_physics_tempo)
then
115 write(errmsg,
'(*(a))')
"Logic error: namelist choice of microphysics is different from TEMPO MP"
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"
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)"
135 if (mpirank==mpiroot)
write(*,*)
'Calling tempo_init() with is_aerosol_aware = ', is_aerosol_aware
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
146 is_initialized = .true.
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, &
183 ar_volume = 4.0 / 3.0 * pi * (2.5e-6)**3
189 where(spechum<0) spechum = 1.0e-10
199 if (merra2_aerosol_aware)
then
200 call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
204 qv = spechum/(1.0_kind_phys-spechum)
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)
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)
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)
227 rho = con_eps*prsl/(con_rd*tgrs*(qv+con_eps))
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
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
240 if (is_hail_aware)
then
241 where(qg .LE. 0.0) chw=0.0
242 where(qg .LE. 0.0) vh=0.0
247 if (is_aerosol_aware)
then
250 if (maxval(nwfa) .lt. eps)
then
251 if (mpirank==mpiroot)
write(*,*)
' Apparently there are no initial CCN aerosols.'
253 if (hgt(i,1).le.1000.0)
then
255 elseif (hgt(i,1).ge.2500.0)
then
258 h_01 = 0.8*cos(hgt(i,1)*0.001 - 1.0)
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)
265 nwfa(i,k) = naccn1+naccn0*exp(-((hgt(i,k)-hgt(i,1))/1000.)*niccn3)
269 if (mpirank==mpiroot)
write(*,*)
' Apparently initial CCN aerosols are present.'
270 if (maxval(nwfa2d) .lt. eps)
then
280 if (mpirank==mpiroot)
write(*,*)
' Apparently there are no initial CCN aerosol surface emission rates.'
282 z1 = hgt(i,2)-hgt(i,1)
283 nwfa2d(i) = nwfa(i,1) * 0.000196 * (50./z1)
286 if (mpirank==mpiroot)
write(*,*)
' Apparently initial CCN aerosol surface emission rates are present.'
291 if (maxval(nifa) .lt. eps)
then
292 if (mpirank==mpiroot)
write(*,*)
' Apparently there are no initial IN aerosols.'
294 if (hgt(i,1).le.1000.0)
then
296 elseif (hgt(i,1).ge.2500.0)
then
299 h_01 = 0.8*cos(hgt(i,1)*0.001 - 1.0)
301 niin3 = -1.0*alog(nain1/nain0)/h_01
302 nifa(i,1) = nain1+nain0*exp(-((hgt(i,2)-hgt(i,1))/1000.)*niin3)
305 nifa(i,k) = nain1+nain0*exp(-((hgt(i,k)-hgt(i,1))/1000.)*niin3)
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.'
315 if (mpirank==mpiroot)
write(*,*)
' Apparently initial IN aerosol surface emission rates are present.'
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
325 where(nwfa .LE. 0.0) nwfa = 1.1e6
326 where(nifa .LE. 0.0) nifa = nain1*0.01
331 else if (merra2_aerosol_aware)
then
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
342 nc_local = nt_c_l/rho
346 if (convert_dry_rho)
then
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)
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)
366 is_initialized = .true.
378 con_eps, convert_dry_rho, &
379 spechum, qc, qr, qi, qs, qg, ni, nr, &
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, &
387 first_time_step, istep, nsteps, &
388 prcp, rain, graupel, ice, snow, sr, &
389 refl_10cm, fullradar_diag, &
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, &
397 cplchm, pfi_lsan, pfl_lsan, &
398 is_initialized, errmsg, errflg)
403 logical,
intent(inout) :: is_initialized
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
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(:,:)
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
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
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(:)
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
455 integer,
intent(in) :: blkno
456 type(mpi_comm),
intent(in) :: mpicomm
457 integer,
intent(in) :: mpirank
458 integer,
intent(in) :: mpiroot
460 logical,
intent(in) :: ext_diag
461 real(kind_phys),
target,
intent(inout),
optional :: diag3d(:,:,:)
462 logical,
intent(in) :: reset_diag3d
465 character(len=*),
intent( out) :: errmsg
466 integer,
intent( out) :: errflg
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(:)
476 logical,
intent (in) :: cplchm
478 real(kind=kind_phys),
intent(inout),
dimension(:,:),
optional :: pfi_lsan
479 real(kind=kind_phys),
intent(inout),
dimension(:,:),
optional :: pfl_lsan
484 real(kind_phys) :: dtstep
487 real(kind_phys) :: rho(1:ncol,1:nlev)
489 real(kind_phys) :: qv(1:ncol,1:nlev)
491 real(kind_phys) :: w(1:ncol,1:nlev)
492 real(kind_phys) :: dz(1:ncol,1:nlev)
494 real(kind_phys) :: rain_mp(1:ncol)
495 real(kind_phys) :: graupel_mp(1:ncol)
496 real(kind_phys) :: ice_mp(1:ncol)
497 real(kind_phys) :: snow_mp(1:ncol)
498 real(kind_phys) :: delta_rain_mp(1:ncol)
499 real(kind_phys) :: delta_graupel_mp(1:ncol)
500 real(kind_phys) :: delta_ice_mp(1:ncol)
501 real(kind_phys) :: delta_snow_mp(1:ncol)
503 real(kind_phys) :: pfils(1:ncol,1:nlev,1)
504 real(kind_phys) :: pflls(1:ncol,1:nlev,1)
507 integer :: do_radar_ref_mp
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
516 integer :: ids,ide, jds,jde, kds,kde, &
517 ims,ime, jms,jme, kms,kme, &
518 its,ite, jts,jte, kts,kte
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()
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.'
571 if (first_time_step .and. istep==1 .and. blkno==1)
then
573 if (.not.is_initialized)
then
574 write(errmsg, fmt=
'((a))')
'mp_tempo_run called before mp_tempo_init'
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'
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'
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"
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'
615 if ( spp_mp==7 )
then
623 dtstep = dtp/real(nsteps, kind=kind_phys)
627 if (merra2_aerosol_aware)
then
628 call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
638 qv = spechum/(1.0_kind_phys-spechum)
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)
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)
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)
662 rho = con_eps*prsl/(con_rd*tgrs*(qv+con_eps))
665 w = -omega/(rho*con_g)
668 dz = (phii(:,2:nlev+1) - phii(:,1:nlev)) / con_g
682 if (do_radar_ref)
then
715 set_extended_diagnostic_pointers:
if (ext_diag)
then
716 if (reset_diag3d)
then
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)
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, &
824 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
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, &
853 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
872 else if (merra2_aerosol_aware)
then
873 write(errmsg,
'(*(a))')
"TEMPO aerosol-aware with MERRA2 UNTESTED -- DO NOT USE"
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, &
899 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
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, &
939 ext_diag=ext_diag, pfils=pfils, pflls=pflls)
958 if (errflg/=0)
return
965 spechum = qv/(1.0_kind_phys+qv)
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)
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)
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)
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)
997 if (nsteps>1 .and. istep == nsteps)
then
999 sr = (snow + graupel + ice)/(rain + snow + graupel + ice +1.e-12)
1004 pfi_lsan(:,:) = pfils(:,:,1)
1005 pfl_lsan(:,:) = pflls(:,:,1)