Line data Source code
1 : module zm_conv_convtran
2 :
3 : use ccpp_kinds, only: kind_phys
4 :
5 : implicit none
6 :
7 : save
8 : private ! Make default type private to the module
9 : !
10 : ! PUBLIC: interfaces
11 : !
12 : public zm_conv_convtran_init
13 : public zm_conv_convtran_run ! convective transport
14 :
15 :
16 : contains
17 :
18 : !===============================================================================
19 : !> \section arg_table_zm_conv_convtran_init Argument Table
20 : !! \htmlinclude zm_conv_convtran_init.html
21 : !!
22 0 : subroutine zm_conv_convtran_init(qprops, ncnst, doconvtran, errmsg, errflg)
23 :
24 : use ccpp_constituent_prop_mod, only: ccpp_constituent_prop_ptr_t
25 :
26 : implicit none
27 :
28 : type(ccpp_constituent_prop_ptr_t), intent(in) :: qprops(:)
29 : integer, intent(in) :: ncnst ! number of tracers to transport
30 :
31 : logical, intent(out) :: doconvtran(:) ! flag for doing convective transport (ncnst)
32 : character(len=512), intent(out) :: errmsg
33 : integer, intent(out) :: errflg
34 :
35 :
36 : integer :: q_index
37 :
38 0 : errmsg = ''
39 0 : errflg = 0
40 :
41 : ! Only convectively transport constituents that are water species
42 0 : do q_index=1,ncnst
43 0 : call qprops(q_index)%is_water_species(doconvtran(q_index), errflg, errmsg)
44 0 : if (errflg /= 0) return
45 : end do
46 :
47 : end subroutine zm_conv_convtran_init
48 : !===============================================================================
49 : !> \section arg_table_zm_conv_convtran_run Argument Table
50 : !! \htmlinclude zm_conv_convtran_run.html
51 : !!
52 123840 : subroutine zm_conv_convtran_run(ncol, pver, &
53 123840 : doconvtran,q ,ncnst ,mu ,md , &
54 123840 : du ,eu ,ed ,dp ,dsubcld , &
55 247680 : jt ,mx ,ideep ,il1g ,il2g , &
56 247680 : nstep ,fracis ,dqdt ,dpdry ,const_metadata, &
57 123840 : scheme_name, errmsg, errflg)
58 :
59 : !-----------------------------------------------------------------------
60 : !
61 : ! Purpose:
62 : ! Convective transport of trace species
63 : !
64 : ! Mixing ratios may be with respect to either dry or moist air
65 : !
66 : ! Method:
67 : ! <Describe the algorithm(s) used in the routine.>
68 : ! <Also include any applicable external references.>
69 : !
70 : ! Author: P. Rasch
71 : !
72 : !-----------------------------------------------------------------------
73 : use ccpp_constituent_prop_mod, only: ccpp_constituent_prop_ptr_t
74 :
75 :
76 : implicit none
77 : !-----------------------------------------------------------------------
78 : !
79 : ! Input arguments
80 : !
81 : integer, intent(in) :: ncol
82 : integer, intent(in) :: pver
83 : integer, intent(in) :: ncnst ! number of tracers to transport
84 : logical, intent(in) :: doconvtran(:) ! flag for doing convective transport (ncnst)
85 : real(kind_phys), intent(in) :: q(:,:,:) ! Tracer array including moisture (ncol,pver,ncnst)
86 : real(kind_phys), intent(in) :: mu(:,:) ! Mass flux up (ncol,pver)
87 : real(kind_phys), intent(in) :: md(:,:) ! Mass flux down (ncol,pver)
88 : real(kind_phys), intent(in) :: du(:,:) ! Mass detraining from updraft (ncol,pver)
89 : real(kind_phys), intent(in) :: eu(:,:) ! Mass entraining from updraft (ncol,pver)
90 : real(kind_phys), intent(in) :: ed(:,:) ! Mass entraining from downdraft (ncol,pver)
91 : real(kind_phys), intent(in) :: dp(:,:) ! Delta pressure between interfaces (ncol,pver)
92 : real(kind_phys), intent(in) :: dsubcld(:) ! Delta pressure from cloud base to sfc (ncol)
93 : real(kind_phys), intent(in) :: fracis(:,:,:) ! fraction of tracer that is insoluble (ncol,pver,ncnst)
94 :
95 : integer, intent(in) :: jt(:) ! Index of cloud top for each column (ncol)
96 : integer, intent(in) :: mx(:) ! Index of cloud top for each column (ncol)
97 : integer, intent(in) :: ideep(:) ! Gathering array (ncol)
98 : integer, intent(in) :: il1g ! Gathered min lon indices over which to operate
99 : integer, intent(in) :: il2g ! Gathered max lon indices over which to operate
100 : integer, intent(in) :: nstep ! Time step index
101 :
102 : real(kind_phys), intent(in) :: dpdry(:,:) ! Delta pressure between interfaces (ncol,pver)
103 :
104 :
105 : type(ccpp_constituent_prop_ptr_t), intent(in) :: const_metadata(:)
106 : character(len=40), intent(out) :: scheme_name
107 : character(len=512), intent(out) :: errmsg
108 : integer, intent(out) :: errflg
109 :
110 : ! input/output
111 :
112 : real(kind_phys), intent(out) :: dqdt(:,:,:) ! Tracer tendency array (ncol,pver,ncnst)
113 :
114 : !--------------------------Local Variables------------------------------
115 :
116 : integer i ! Work index
117 : integer k ! Work index
118 : integer kbm ! Highest altitude index of cloud base
119 : integer kk ! Work index
120 : integer kkp1 ! Work index
121 : integer km1 ! Work index
122 : integer kp1 ! Work index
123 : integer ktm ! Highest altitude index of cloud top
124 : integer m ! Work index
125 :
126 : logical :: is_dry
127 :
128 : real(kind_phys) cabv ! Mix ratio of constituent above
129 : real(kind_phys) cbel ! Mix ratio of constituent below
130 : real(kind_phys) cdifr ! Normalized diff between cabv and cbel
131 247680 : real(kind_phys) chat(ncol,pver) ! Mix ratio in env at interfaces
132 247680 : real(kind_phys) cond(ncol,pver) ! Mix ratio in downdraft at interfaces
133 247680 : real(kind_phys) const(ncol,pver) ! Gathered tracer array
134 247680 : real(kind_phys) fisg(ncol,pver) ! gathered insoluble fraction of tracer
135 247680 : real(kind_phys) conu(ncol,pver) ! Mix ratio in updraft at interfaces
136 247680 : real(kind_phys) dcondt(ncol,pver) ! Gathered tend array
137 : real(kind_phys) small ! A small number
138 : real(kind_phys) mbsth ! Threshold for mass fluxes
139 : real(kind_phys) mupdudp ! A work variable
140 : real(kind_phys) minc ! A work variable
141 : real(kind_phys) maxc ! A work variable
142 : real(kind_phys) fluxin ! A work variable
143 : real(kind_phys) fluxout ! A work variable
144 : real(kind_phys) netflux ! A work variable
145 :
146 247680 : real(kind_phys) dutmp(ncol,pver) ! Mass detraining from updraft
147 247680 : real(kind_phys) eutmp(ncol,pver) ! Mass entraining from updraft
148 247680 : real(kind_phys) edtmp(ncol,pver) ! Mass entraining from downdraft
149 123840 : real(kind_phys) dptmp(ncol,pver) ! Delta pressure between interfaces
150 : real(kind_phys) total(ncol)
151 : real(kind_phys) negadt,qtmp
152 :
153 : character(len=256) :: standard_name
154 :
155 : !-----------------------------------------------------------------------
156 : !
157 123840 : scheme_name = "zm_conv_convtran_run"
158 123840 : errmsg = ''
159 123840 : errflg = 0
160 :
161 123840 : small = 1.e-36_kind_phys
162 : ! mbsth is the threshold below which we treat the mass fluxes as zero (in mb/s)
163 123840 : mbsth = 1.e-15_kind_phys
164 :
165 : ! Find the highest level top and bottom levels of convection
166 123840 : ktm = pver
167 123840 : kbm = pver
168 464242 : do i = il1g, il2g
169 340402 : ktm = min(ktm,jt(i))
170 464242 : kbm = min(kbm,mx(i))
171 : end do
172 :
173 : ! Loop ever each constituent
174 7889875200 : dqdt(:,:,:) = 0._kind_phys
175 5201280 : do m = 1, ncnst
176 :
177 5077440 : call const_metadata(m)%standard_name(standard_name)
178 5077440 : if (standard_name == 'water_vapor_mixing_ratio_wrt_moist_air_and_condensed_water') then
179 : cycle
180 : end if
181 :
182 5077440 : if (doconvtran(m)) then
183 :
184 1297224 : call const_metadata(m)%is_dry(is_dry, errflg, errmsg)
185 1297224 : if (is_dry) then
186 60824016 : do k = 1,pver
187 226276875 : do i =il1g,il2g
188 165452859 : dptmp(i,k) = dpdry(i,k)
189 165452859 : dutmp(i,k) = du(i,k)*dp(i,k)/dpdry(i,k)
190 165452859 : eutmp(i,k) = eu(i,k)*dp(i,k)/dpdry(i,k)
191 225629811 : edtmp(i,k) = ed(i,k)*dp(i,k)/dpdry(i,k)
192 : end do
193 : end do
194 : else
195 61115040 : do k = 1,pver
196 227277210 : do i =il1g,il2g
197 166162170 : dptmp(i,k) = dp(i,k)
198 166162170 : dutmp(i,k) = du(i,k)
199 166162170 : eutmp(i,k) = eu(i,k)
200 226627050 : edtmp(i,k) = ed(i,k)
201 : end do
202 : end do
203 : endif
204 :
205 : ! Gather up the constituent and set tend to zero
206 121939056 : do k = 1,pver
207 453554085 : do i =il1g,il2g
208 331615029 : const(i,k) = q(ideep(i),k,m)
209 452256861 : fisg(i,k) = fracis(ideep(i),k,m)
210 : end do
211 : end do
212 :
213 : ! From now on work only with gathered data
214 :
215 : ! Interpolate environment tracer values to interfaces
216 121939056 : do k = 1,pver
217 120641832 : km1 = max(1,k-1)
218 453554085 : do i = il1g, il2g
219 331615029 : minc = min(const(i,km1),const(i,k))
220 331615029 : maxc = max(const(i,km1),const(i,k))
221 331615029 : if (minc < 0) then
222 : cdifr = 0._kind_phys
223 : else
224 331615029 : cdifr = abs(const(i,k)-const(i,km1))/max(maxc,small)
225 : endif
226 :
227 : ! If the two layers differ significantly use a geometric averaging
228 : ! procedure
229 331615029 : if (cdifr > 1.E-6_kind_phys) then
230 239709643 : cabv = max(const(i,km1),maxc*1.e-12_kind_phys)
231 239709643 : cbel = max(const(i,k),maxc*1.e-12_kind_phys)
232 239709643 : chat(i,k) = log(cabv/cbel)/(cabv-cbel)*cabv*cbel
233 :
234 : else ! Small diff, so just arithmetic mean
235 91905386 : chat(i,k) = 0.5_kind_phys* (const(i,k)+const(i,km1))
236 : end if
237 :
238 : ! Provisional up and down draft values
239 331615029 : conu(i,k) = chat(i,k)
240 331615029 : cond(i,k) = chat(i,k)
241 :
242 : ! provisional tends
243 452256861 : dcondt(i,k) = 0._kind_phys
244 :
245 : end do
246 : end do
247 :
248 : ! Do levels adjacent to top and bottom
249 4862977 : k = 2
250 4862977 : km1 = 1
251 4862977 : kk = pver
252 4862977 : do i = il1g,il2g
253 3565753 : mupdudp = mu(i,kk) + dutmp(i,kk)*dptmp(i,kk)
254 3565753 : if (mupdudp > mbsth) then
255 0 : conu(i,kk) = (+eutmp(i,kk)*fisg(i,kk)*const(i,kk)*dptmp(i,kk))/mupdudp
256 : endif
257 4862977 : if (md(i,k) < -mbsth) then
258 0 : cond(i,k) = (-edtmp(i,km1)*fisg(i,km1)*const(i,km1)*dptmp(i,km1))/md(i,k)
259 : endif
260 : end do
261 :
262 : ! Updraft from bottom to top
263 120641832 : do kk = pver-1,1,-1
264 119344608 : kkp1 = min(pver,kk+1)
265 448691108 : do i = il1g,il2g
266 328049276 : mupdudp = mu(i,kk) + dutmp(i,kk)*dptmp(i,kk)
267 447393884 : if (mupdudp > mbsth) then
268 56686832 : conu(i,kk) = ( mu(i,kkp1)*conu(i,kkp1)+eutmp(i,kk)*fisg(i,kk)* &
269 113373664 : const(i,kk)*dptmp(i,kk) )/mupdudp
270 : endif
271 : end do
272 : end do
273 :
274 : ! Downdraft from top to bottom
275 119344608 : do k = 3,pver
276 118047384 : km1 = max(1,k-1)
277 443828131 : do i = il1g,il2g
278 442530907 : if (md(i,k) < -mbsth) then
279 87364302 : cond(i,k) = ( md(i,km1)*cond(i,km1)-edtmp(i,km1)*fisg(i,km1)*const(i,km1) &
280 87364302 : *dptmp(i,km1) )/md(i,k)
281 : endif
282 : end do
283 : end do
284 :
285 :
286 12970142 : do k = ktm,pver
287 11672918 : km1 = max(1,k-1)
288 11672918 : kp1 = min(pver,k+1)
289 97893358 : do i = il1g,il2g
290 :
291 : ! limit fluxes outside convection to mass in appropriate layer
292 : ! these limiters are probably only safe for positive definite quantitities
293 : ! it assumes that mu and md already satify a courant number limit of 1
294 339692864 : fluxin = mu(i,kp1)*conu(i,kp1)+ mu(i,k)*min(chat(i,k),const(i,km1)) &
295 424616080 : -(md(i,k) *cond(i,k) + md(i,kp1)*min(chat(i,kp1),const(i,kp1)))
296 : fluxout = mu(i,k)*conu(i,k) + mu(i,kp1)*min(chat(i,kp1),const(i,k)) &
297 84923216 : -(md(i,kp1)*cond(i,kp1) + md(i,k)*min(chat(i,k),const(i,k)))
298 :
299 84923216 : netflux = fluxin - fluxout
300 84923216 : if (abs(netflux) < max(fluxin,fluxout)*1.e-12_kind_phys) then
301 5346420 : netflux = 0._kind_phys
302 : endif
303 96596134 : dcondt(i,k) = netflux/dptmp(i,k)
304 : end do
305 : end do
306 : ! %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
307 : !
308 5010604 : do k = kbm,pver
309 3713380 : km1 = max(1,k-1)
310 26957653 : do i = il1g,il2g
311 25660429 : if (k == mx(i)) then
312 :
313 3565753 : fluxin = mu(i,k)*min(chat(i,k),const(i,km1)) - md(i,k)*cond(i,k)
314 3565753 : fluxout = mu(i,k)*conu(i,k) - md(i,k)*min(chat(i,k),const(i,k))
315 :
316 3565753 : netflux = fluxin - fluxout
317 3565753 : if (abs(netflux) < max(fluxin,fluxout)*1.e-12_kind_phys) then
318 0 : netflux = 0._kind_phys
319 : endif
320 3565753 : dcondt(i,k) = netflux/dptmp(i,k)
321 18381296 : else if (k > mx(i)) then
322 15674324 : dcondt(i,k) = 0._kind_phys
323 : end if
324 : end do
325 : end do
326 :
327 : ! Scatter tendency back to full array
328 121939056 : do k = 1,pver
329 120641832 : kp1 = min(pver,k+1)
330 453554085 : do i = il1g,il2g
331 452256861 : dqdt(ideep(i),k,m) = dcondt(i,k)
332 : end do
333 : end do
334 :
335 : end if ! for doconvtran
336 :
337 : end do
338 :
339 123840 : return
340 : end subroutine zm_conv_convtran_run
341 :
342 :
343 : end module zm_conv_convtran
|