CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
sfc_diag.f
1
3
5 module sfc_diag
6 contains
7
13 subroutine sfc_diag_run (im,xlat_d,xlon_d, &
14 & lsm,lsm_ruc,grav,cp,eps,epsm1,con_rocp, &
15 & con_karman, &
16 & shflx,cdq,wind, &
17 & usfco,vsfco,use_oceanuv, &
18 & zf,ps,u1,v1,t1,q1,prslki,evap,fm,fh,fm10,fh2, &
19 & ust,tskin,qsurf,thsfc_loc,diag_flux,diag_log, &
20 & use_lake_model,iopt_lake,iopt_lake_clm, &
21 & lake_t2m,lake_q2m,use_lake2m, &
22 & f10m,u10m,v10m,t2m,q2m,dpt2m,errmsg,errflg &
23 & )
24!
25 use machine , only : kind_phys, kind_dbl_prec
26 use funcphys, only : fpvs
27 use physcons, only : con_t0c
28 implicit none
29!
30 integer, intent(in) :: im, lsm, lsm_ruc, iopt_lake, iopt_lake_clm
31 logical, intent(in) :: use_lake2m
32 logical, intent(in) :: use_oceanuv
33 logical, intent(in) :: thsfc_loc ! Flag for reference pot. temp.
34 logical, intent(in) :: diag_flux ! Flag for flux method in 2-m diagnostics
35 logical, intent(in) :: diag_log ! Flag for 2-m log diagnostics under stable conditions
36 real(kind=kind_phys), intent(in) :: grav,cp,eps,epsm1,con_rocp
37 real(kind=kind_phys), intent(in) :: con_karman
38 real(kind=kind_phys), dimension(:), intent( in) :: &
39 & zf, ps, u1, v1, t1, q1, ust, tskin, &
40 & usfco, vsfco, &
41 & qsurf, prslki, evap, fm, fh, fm10, fh2, &
42 & shflx, cdq, wind, xlat_d, xlon_d
43 real(kind=kind_phys), dimension(:), intent(out) :: &
44 & f10m, u10m, v10m, t2m, q2m, dpt2m
45 real(kind=kind_phys), dimension(:), intent(in), optional :: &
46 & lake_t2m, lake_q2m
47 integer, dimension(:), intent(in) :: use_lake_model
48 character(len=*), intent(out) :: errmsg
49 integer, intent(out) :: errflg
50!
51! locals
52!
53 real (kind_phys), parameter :: zero = 0._kind_dbl_prec
54 real (kind_phys), parameter :: one = 1._kind_dbl_prec
55 real (kind_phys), parameter :: qmin=1.0e-8
56 integer :: k,i
57
58 logical :: debug_print
59 real(kind=kind_phys) :: q1c, qv, tem, qv1, th2m, x2m, rho
60 real(kind=kind_phys) :: dt, dq, qsfcmr, qsfcprox, ff, fac, dz1
61 real(kind=kind_phys) :: t2_alt, q2_alt
62 real(kind=kind_phys) :: thcon, cqs, chs, chs2, cqs2
63 real(kind=kind_phys) :: testptlat, testptlon
64 logical :: have_2m
65!
66 real(kind=kind_phys) :: fhi, qss, wrk
67! real(kind=kind_phys) sig2k, fhi, qss
68!
69! real, parameter :: g=grav
70!
71 ! Initialize CCPP error handling variables
72 errmsg = ''
73 errflg = 0
74
75 !--
76 testptlat = 35.3_kind_phys
77 testptlon = 273.0_kind_phys
78 !--
79 debug_print = .false.
80!
81! estimate sigma ** k at 2 m
82!
83! sig2k = 1. - 4. * g * 2. / (cp * 280.)
84!
85! initialize variables. all units are supposedly m.k.s. unless specified
86! ps is in pascals
87!
88!!
89
90 do i = 1, im
91 f10m(i) = fm10(i) / fm(i)
92 if (use_oceanuv) then
93 u10m(i) = usfco(i)+f10m(i) * (u1(i)-usfco(i))
94 v10m(i) = vsfco(i)+f10m(i) * (v1(i)-vsfco(i))
95 else
96 u10m(i) = f10m(i) * u1(i)
97 v10m(i) = f10m(i) * v1(i)
98 endif
99
100 have_2m = use_lake_model(i)>0 .and. use_lake2m .and. &
101 & iopt_lake==iopt_lake_clm
102 if(have_2m) then
103 t2m(i) = lake_t2m(i)
104 q2m(i) = lake_q2m(i)
105 endif
106 fhi = fh2(i) / fh(i)
107 wrk = 1.0 - fhi
108
109 if(lsm /= lsm_ruc) then
110 !-- original method
111 if(have_2m) then
112 ! already have 2m T & Q from lake
113 else
114 if(thsfc_loc) then ! Use local potential temperature
115 t2m(i)=tskin(i)*wrk + t1(i)*prslki(i)*fhi - (grav+grav)/cp
116 else ! Use potential temperature referenced to 1000 hPa
117 t2m(i) = tskin(i)*wrk + t1(i)*fhi - (grav+grav)/cp
118 endif
119 if(evap(i) >= 0.) then ! for evaporation>0, use inferred qsurf to deduce q2m
120 q2m(i) = qsurf(i)*wrk + max(qmin,q1(i))*fhi
121 else ! for dew formation, use saturated q at tskin
122 qss = fpvs(tskin(i))
123 qss = eps * qss/(ps(i) + epsm1 * qss)
124 q2m(i) = qss*wrk + max(qmin,q1(i))*fhi
125 endif
126 endif
127 qss = fpvs(t2m(i))
128 qss = eps * qss / (ps(i) + epsm1 * qss)
129 q2m(i) = min(q2m(i),qss)
130 else
131 !RUC lsm
132 thcon = (1.e5_kind_phys/ps(i))**con_rocp
133 !-- make sure 1st level q is not higher than saturated value
134 qss = fpvs(t1(i))
135 qss = eps * qss / (ps(i) + epsm1 * qss)
136 q1c = min(q1(i),qss) ! lev 1 spec. humidity
137
138 qv1 = q1c / (one - q1c) ! lev 1 mixing ratio
139 qsfcmr = qsurf(i)/(one - qsurf(i)) ! surface mixing ratio
140 chs = cdq(i) * wind(i)
141 cqs = chs
142 chs2 = ust(i)*con_karman/fh2(i)
143 cqs2 = chs2
144 qsfcprox = max(qmin,qv1 + evap(i)/cqs) ! surface mix. ratio computed from the flux
145
146 ruc_have_2m: if(.not.have_2m) then
147 if(diag_flux) then
148 !-- flux method
149 th2m = tskin(i)*thcon - shflx(i)/chs2
150 t2m(i) = th2m/thcon
151 x2m = max(qmin,qsfcprox - evap(i)/cqs2) ! mix. ratio
152 q2m(i) = x2m/(one + x2m) ! spec. humidity
153 else
154 t2m(i) = tskin(i)*wrk + t1(i)*fhi - (grav+grav)/cp
155 q2m(i) = qsurf(i)*wrk + max(qmin,q1c)*fhi
156 endif ! flux method
157
158 if(diag_log) then
159 !-- Alternative logarithmic diagnostics:
160 dt = t1(i) - tskin(i)
161 dq = qv1 - qsfcmr
162 dz1= zf(i) ! level of atm. forcing
163 IF (dt > zero) THEN
164 ff = min(max(one-dt/10._kind_phys,0.01_kind_phys), one)
165 !for now, set zt = 0.05
166 fac = log((2._kind_phys +.05_kind_phys)/(0.05_kind_phys + &
167 & ff))/log((dz1 + .05_kind_phys)/(0.05_kind_phys + ff))
168 t2_alt = tskin(i) + fac * dt
169 ELSE
170 !no alternatives (yet) for unstable conditions
171 t2_alt = t2m(i)
172 ENDIF
173
174 IF (dq > zero) THEN
175 ff = min(max(one-dq/0.003_kind_phys,0.01_kind_phys), one)
176 !-- for now, set zt = 0.05
177 fac = log((2._kind_phys +.05_kind_phys)/(0.05_kind_phys + &
178 & ff))/log((dz1 + .05_kind_phys)/(0.05_kind_phys + ff))
179 q2_alt = qsfcmr + fac * dq ! mix. ratio
180 q2_alt = q2_alt/(one + q2_alt) ! spec. humidity
181 ELSE
182 !no alternatives (yet) for unstable conditions
183 q2_alt = q2m(i)
184 ENDIF
185 !-- Note: use of alternative diagnostics will make
186 ! it cooler and drier with stable stratification
187 t2m(i) = t2_alt
188 q2m(i) = q2_alt
189 endif ! log method for stable regime
190
191 !-- check that T2m values lie in the range between tskin and t1
192 x2m = max(min(tskin(i),t1(i)) , t2m(i))
193 t2m(i) = min(max(tskin(i),t1(i)) , x2m)
194 !-- check that Q2m values lie in the range between qsurf and q1
195 x2m = max(min(qsurf(i),q1c) , q2m(i))
196 q2m(i) = min(max(qsurf(i),q1c) , x2m)
197
198 !-- make sure q2m is not oversaturated
199 qss = fpvs(t2m(i))
200 qss = eps * qss/(ps(i) + epsm1 * qss)
201 q2m(i) = min(q2m(i),qss)
202
203 if(diag_flux) then
204 !-- from WRF, applied in HRRR - Jimy Dudhia
205 ! Limit Q2m diagnostic to no more than 5 percent higher than lowest level value
206 ! This prevents unrealistic values when QFX is not mostly surface
207 ! flux because calculation is based on surface flux only.
208 ! Problems occurred in transition periods and weak winds and vegetation source
209 q2m(i) = min(q2m(i),1.05_kind_dbl_prec*q1c) ! works if qsurf > q1c, evaporation
210 endif
211 endif ruc_have_2m
212
213
214 !-- Compute dew point, using vapor pressure
215 qv = max(qmin,(q2m(i)/(1.-q2m(i))))
216 tem = max(ps(i) * qv/( eps - epsm1 *qv), qmin)
217 dpt2m(i) = 243.5_kind_dbl_prec/( ( 17.67_kind_dbl_prec / &
218 & log(tem/611.2_kind_dbl_prec) ) - one) + con_t0c
219 dpt2m(i) = min(dpt2m(i),t2m(i))
220
221
222 if (debug_print) then
223 !-- diagnostics for a test point with known lat/lon
224 if (abs(xlat_d(i)-testptlat).lt.0.2 .and. &
225 & abs(xlon_d(i)-testptlon).lt.0.2)then
226 print 100,'(ruc_lsm_diag) i=',i, &
227 & ' lat,lon=',xlat_d(i),xlon_d(i),'zf ',zf(i), &
228 & 'tskin ',tskin(i),'t2m ',t2m(i),'t1',t1(i),'shflx',shflx(i),&
229 & 'qsurf ',qsurf(i),'qsfcprox ',qsfcprox,'q2m ',q2m(i), &
230 & 'q1 ',q1(i),'evap ',evap(i),'dpt2m ',dpt2m(i), &
231 & 'chs2 ',chs2,'cqs2 ',cqs2,'cqs ',cqs,'cdq',cdq(i)
232 endif
233 endif
234 100 format (";;; ",a,i4,a,2f14.7/(4(a10,'='es11.4)))
235 endif ! RUC LSM
236 enddo
237
238 return
239 end subroutine sfc_diag_run
241
242 end module sfc_diag
subroutine shflx(nsoil, smc, smcmax, dt, yy, zz1, zsoil, zbot, psisat, bexp, df1, ice, quartz, csoil, vegtyp, shdfac, lheatstrg, stc, t1, tbot, sh2o, ssoil)
This subroutine updates the temperature state of the soil column based on the thermal diffusion equat...
Definition sflx.f:3303
subroutine sfc_diag_run(im, xlat_d, xlon_d, lsm, lsm_ruc, grav, cp, eps, epsm1, con_rocp, con_karman, shflx, cdq, wind, usfco, vsfco, use_oceanuv, zf, ps, u1, v1, t1, q1, prslki, evap, fm, fh, fm10, fh2, ust, tskin, qsurf, thsfc_loc, diag_flux, diag_log, use_lake_model, iopt_lake, iopt_lake_clm, lake_t2m, lake_q2m, use_lake2m, f10m, u10m, v10m, t2m, q2m, dpt2m, errmsg, errflg)
Definition sfc_diag.f:24
This module contain the RUC land surface model driver.
Definition lsm_ruc.F90:5
This module contains the land surface diagnose calcualtion.
Definition sfc_diag.f:5