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