Skip to content

Commit 61b419f

Browse files
authored
Merge pull request #1032 from grantfirl/ufs-dev-PR88
UFS-dev PR#88
2 parents c83cd17 + e1ccc99 commit 61b419f

File tree

3 files changed

+24
-23
lines changed

3 files changed

+24
-23
lines changed

physics/GFS_rrtmgp_cloud_mp.F90

Lines changed: 1 addition & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -875,17 +875,12 @@ subroutine cmp_reff_Thompson(nLev, nCol, i_cldliq, i_cldice, i_cldsnow, i_cldice
875875
qi_mp(iCol,iLay) = tracer(iCol,iLay,i_cldice) / (1.-q_lay(iCol,iLay))
876876
qs_mp(iCol,iLay) = tracer(iCol,iLay,i_cldsnow) / (1.-q_lay(iCol,iLay))
877877
ni_mp(iCol,iLay) = tracer(iCol,iLay,i_cldice_nc) / (1.-q_lay(iCol,iLay))
878-
if (ltaerosol) then
878+
if (ltaerosol .or. mraerosol) then
879879
nc_mp(iCol,iLay) = tracer(iCol,iLay,i_cldliq_nc) / (1.-q_lay(iCol,iLay))
880880
nwfa(iCol,iLay) = tracer(iCol,iLay,i_twa)
881881
if (qc_mp(iCol,iLay) > 1.e-12 .and. nc_mp(iCol,iLay) < 100.) then
882882
nc_mp(iCol,iLay) = make_DropletNumber(qc_mp(iCol,iLay)*rho, nwfa(iCol,iLay)*rho) * orho
883883
endif
884-
elseif (mraerosol) then
885-
nc_mp(iCol,iLay) = tracer(iCol,iLay,i_cldliq_nc) / (1.-q_lay(iCol,iLay))
886-
if (qc_mp(iCol,iLay) > 1.e-12 .and. nc_mp(iCol,iLay) < 100.) then
887-
nc_mp(iCol,iLay) = make_DropletNumber(qc_mp(iCol,iLay)*rho, nwfa(iCol,iLay)*rho) * orho
888-
endif
889884
else
890885
if (nint(lsmask(iCol)) == 1) then !land
891886
nc_mp(iCol,iLay) = nt_c_l*orho

physics/module_mp_thompson.F90

Lines changed: 12 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -3400,8 +3400,8 @@ subroutine mp_thompson (qv1d, qc1d, qi1d, qr1d, qs1d, qg1d, ni1d, &
34003400
tcond(k) = (5.69 + 0.0168*tempc)*1.0E-5 * 418.936
34013401
ocp(k) = 1./(Cp*(1.+0.887*qv(k)))
34023402
lvt2(k)=lvap(k)*lvap(k)*ocp(k)*oRv*otemp*otemp
3403-
3404-
nwfa(k) = MAX(11.1E6*rho(k), (nwfa1d(k) + nwfaten(k)*DT)*rho(k))
3403+
if (is_aerosol_aware) &
3404+
nwfa(k) = MAX(11.1E6*rho(k), (nwfa1d(k) + nwfaten(k)*DT)*rho(k))
34053405
enddo
34063406

34073407
do k = kts, kte
@@ -3654,7 +3654,8 @@ subroutine mp_thompson (qv1d, qc1d, qi1d, qr1d, qs1d, qg1d, ni1d, &
36543654
qvten(k) = qvten(k) - prw_vcd(k)
36553655
qcten(k) = qcten(k) + prw_vcd(k)
36563656
ncten(k) = ncten(k) + pnc_wcd(k)
3657-
nwfaten(k) = nwfaten(k) - pnc_wcd(k)
3657+
if (is_aerosol_aware) &
3658+
nwfaten(k) = nwfaten(k) - pnc_wcd(k)
36583659
tten(k) = tten(k) + lvap(k)*ocp(k)*prw_vcd(k)*(1-IFDRY)
36593660
rc(k) = MAX(R1, (qc1d(k) + DT*qcten(k))*rho(k))
36603661
if (rc(k).eq.R1) L_qc(k) = .false.
@@ -3743,7 +3744,8 @@ subroutine mp_thompson (qv1d, qc1d, qi1d, qr1d, qs1d, qg1d, ni1d, &
37433744
qrten(k) = qrten(k) - prv_rev(k)
37443745
qvten(k) = qvten(k) + prv_rev(k)
37453746
nrten(k) = nrten(k) - pnr_rev(k)
3746-
nwfaten(k) = nwfaten(k) + pnr_rev(k)
3747+
if (is_aerosol_aware) &
3748+
nwfaten(k) = nwfaten(k) + pnr_rev(k)
37473749
tten(k) = tten(k) - lvap(k)*ocp(k)*prv_rev(k)*(1-IFDRY)
37483750

37493751
rr(k) = MAX(R1, (qr1d(k) + DT*qrten(k))*rho(k))
@@ -4232,10 +4234,12 @@ subroutine mp_thompson (qv1d, qc1d, qi1d, qr1d, qs1d, qg1d, ni1d, &
42324234
qv1d(k) = MAX(1.E-10, qv1d(k) + qvten(k)*DT)
42334235
qc1d(k) = qc1d(k) + qcten(k)*DT
42344236
nc1d(k) = MAX(2./rho(k), MIN(nc1d(k) + ncten(k)*DT, Nt_c_max))
4235-
nwfa1d(k) = MAX(11.1E6, MIN(9999.E6, &
4236-
(nwfa1d(k)+nwfaten(k)*DT)))
4237-
nifa1d(k) = MAX(naIN1*0.01, MIN(9999.E6, &
4238-
(nifa1d(k)+nifaten(k)*DT)))
4237+
if (is_aerosol_aware) then
4238+
nwfa1d(k) = MAX(11.1E6, MIN(9999.E6, &
4239+
(nwfa1d(k)+nwfaten(k)*DT)))
4240+
nifa1d(k) = MAX(naIN1*0.01, MIN(9999.E6, &
4241+
(nifa1d(k)+nifaten(k)*DT)))
4242+
end if
42394243
if (qc1d(k) .le. R1) then
42404244
qc1d(k) = 0.0
42414245
nc1d(k) = 0.0

physics/mp_thompson.F90

Lines changed: 11 additions & 9 deletions
Original file line numberDiff line numberDiff line change
@@ -151,6 +151,10 @@ subroutine mp_thompson_init(ncol, nlev, con_g, con_rd, con_eps, &
151151
!> - Convert specific humidity to water vapor mixing ratio.
152152
!> - Also, hydrometeor variables are mass or number mixing ratio
153153
!> - either kg of species per kg of dry air, or per kg of (dry + vapor).
154+
if (merra2_aerosol_aware) then
155+
call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
156+
end if
157+
154158

155159
qv = spechum/(1.0_kind_phys-spechum)
156160

@@ -163,7 +167,7 @@ subroutine mp_thompson_init(ncol, nlev, con_g, con_rd, con_eps, &
163167

164168
ni = ni/(1.0_kind_phys-spechum)
165169
nr = nr/(1.0_kind_phys-spechum)
166-
if (is_aerosol_aware) then
170+
if (is_aerosol_aware .or. merra2_aerosol_aware) then
167171
nc = nc/(1.0_kind_phys-spechum)
168172
nwfa = nwfa/(1.0_kind_phys-spechum)
169173
nifa = nifa/(1.0_kind_phys-spechum)
@@ -208,8 +212,6 @@ subroutine mp_thompson_init(ncol, nlev, con_g, con_rd, con_eps, &
208212
nwfa(i,k) = naCCN1+naCCN0*exp(-((hgt(i,k)-hgt(i,1))/1000.)*niCCN3)
209213
enddo
210214
enddo
211-
else if (merra2_aerosol_aware) then
212-
call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
213215
else
214216
if (mpirank==mpiroot) write(*,*) ' Apparently initial CCN aerosols are present.'
215217
if (MAXVAL(nwfa2d) .lt. eps) then
@@ -555,6 +557,9 @@ subroutine mp_thompson_run(ncol, nlev, con_g, con_rd, &
555557
else
556558
dtstep = dtp
557559
end if
560+
if (merra2_aerosol_aware) then
561+
call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
562+
end if
558563

559564
!> - Convert specific humidity to water vapor mixing ratio.
560565
!> - Also, hydrometeor variables are mass or number mixing ratio
@@ -574,7 +579,7 @@ subroutine mp_thompson_run(ncol, nlev, con_g, con_rd, &
574579

575580
ni = ni/(1.0_kind_phys-spechum)
576581
nr = nr/(1.0_kind_phys-spechum)
577-
if (is_aerosol_aware) then
582+
if (is_aerosol_aware .or. merra2_aerosol_aware) then
578583
nc = nc/(1.0_kind_phys-spechum)
579584
nwfa = nwfa/(1.0_kind_phys-spechum)
580585
nifa = nifa/(1.0_kind_phys-spechum)
@@ -681,9 +686,6 @@ subroutine mp_thompson_run(ncol, nlev, con_g, con_rd, &
681686
ncten3 => diag3d(:,:,36:36)
682687
qcten3 => diag3d(:,:,37:37)
683688
end if set_extended_diagnostic_pointers
684-
if (merra2_aerosol_aware) then
685-
call get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
686-
end if
687689
!> - Call mp_gt_driver() with or without aerosols, with or without effective radii, ...
688690
if (is_aerosol_aware .or. merra2_aerosol_aware) then
689691
call mp_gt_driver(qv=qv, qc=qc, qr=qr, qi=qi, qs=qs, qg=qg, ni=ni, nr=nr, &
@@ -921,8 +923,8 @@ subroutine get_niwfa(aerfld, nifa, nwfa, ncol, nlev)
921923
aerfld(:,:,4)/1011.5142+ aerfld(:,:,5)/5683.3501)*1.e15
922924

923925
nwfa=((aerfld(:,:,6)/0.0045435214+aerfld(:,:,7)/0.2907854+aerfld(:,:,8)/12.91224+ &
924-
aerfld(:,:,9)/206.2216+ aerfld(:,:,10)/4326.23)*1.+aerfld(:,:,11)/0.3053104*5+ &
925-
aerfld(:,:,15)/0.3232698*1)*1.e15
926+
aerfld(:,:,9)/206.2216+ aerfld(:,:,10)/4326.23)*9.+aerfld(:,:,11)/0.3053104*5+ &
927+
aerfld(:,:,15)/0.3232698*8)*1.e15
926928
end subroutine get_niwfa
927929

928930
end module mp_thompson

0 commit comments

Comments
 (0)