1 | !======================================================================= |
---|
2 | ! CALLTHERM_INTERFACE |
---|
3 | !======================================================================= |
---|
4 | ! Main interface to the Martian thermal plume model |
---|
5 | ! This interface handles sub-timesteps for this model |
---|
6 | ! A call to this interface must be inserted in the main 'physics' routine |
---|
7 | ! NB: for information: |
---|
8 | ! In the Mars LMD-GCM, the thermal plume model is called after |
---|
9 | ! the vertical turbulent mixing scheme (Mellor and Yamada) |
---|
10 | ! and the surface layer scheme (Richardson-based surface layer + subgrid gustiness) |
---|
11 | ! Other routines called before the thermals model are |
---|
12 | ! radiative transfer and (orographic) gravity wave drag. |
---|
13 | ! ----------------------------------------------------------------------- |
---|
14 | ! Author : A. Colaitis 2011-01-05 (with updates 2011-2013) |
---|
15 | ! after C. Rio and F. Hourdin |
---|
16 | ! Institution : Laboratoire de Meteorologie Dynamique (LMD) Paris, France |
---|
17 | ! ----------------------------------------------------------------------- |
---|
18 | ! Corresponding author : A. Spiga aymeric.spiga_AT_upmc.fr |
---|
19 | ! ----------------------------------------------------------------------- |
---|
20 | ! ASSOCIATED FILES |
---|
21 | ! --> thermcell_main_mars.F90 |
---|
22 | ! --> thermcell_dqup.F90 |
---|
23 | ! --> comtherm_h.F90 |
---|
24 | ! ----------------------------------------------------------------------- |
---|
25 | ! Reference paper: |
---|
26 | ! A. Colaïtis, A. Spiga, F. Hourdin, C. Rio, F. Forget, and E. Millour. |
---|
27 | ! A thermal plume model for the Martian convective boundary layer. |
---|
28 | ! Journal of Geophysical Research (Planets), 118:1468-1487, July 2013. |
---|
29 | ! http://dx.doi.org/10.1002/jgre.20104 |
---|
30 | ! http://arxiv.org/abs/1306.6215 |
---|
31 | ! ----------------------------------------------------------------------- |
---|
32 | ! Reference paper for terrestrial plume model: |
---|
33 | ! C. Rio and F. Hourdin. |
---|
34 | ! A thermal plume model for the convective boundary layer : Representation of cumulus clouds. |
---|
35 | ! Journal of the Atmospheric Sciences, 65:407-425, 2008. |
---|
36 | ! ----------------------------------------------------------------------- |
---|
37 | |
---|
38 | SUBROUTINE calltherm_interface (ngrid,nlayer,nq, igcm_co2, & |
---|
39 | & zzlev,zzlay, & |
---|
40 | & ptimestep,pu,pv,pt,pq,pdu,pdv,pdt,pdq,q2, & |
---|
41 | & pplay,pplev,pphi,zpopsk, & |
---|
42 | & pdu_th,pdv_th,pdt_th,pdq_th,lmax,zmaxth,pbl_dtke, & |
---|
43 | & pdhdif,hfmax,wstar,sensibFlux) |
---|
44 | |
---|
45 | use comtherm_h |
---|
46 | use tracer_mod, only: nqmx,noms |
---|
47 | |
---|
48 | ! SHARED VARIABLES. This needs adaptations in another climate model. |
---|
49 | ! contains physical constant values such as |
---|
50 | ! "g" : gravitational acceleration (m.s-2) |
---|
51 | ! "r" : recuced gas constant (J.K-1.mol-1) |
---|
52 | ! "cpp" : specific heat of the atmosphere (J.kg-1.K-1) |
---|
53 | USE comcstfi_h |
---|
54 | |
---|
55 | implicit none |
---|
56 | |
---|
57 | !-------------------------------------------------------- |
---|
58 | ! Input Variables |
---|
59 | !-------------------------------------------------------- |
---|
60 | |
---|
61 | INTEGER, INTENT(IN) :: ngrid ! number of horizontal grid points |
---|
62 | INTEGER, INTENT(IN) :: nlayer ! number of vertical grid points |
---|
63 | INTEGER, INTENT(IN) :: nq ! number of tracer species |
---|
64 | REAL, INTENT(IN) :: ptimestep !timestep (s) |
---|
65 | REAL, INTENT(IN) :: pplev(ngrid,nlayer+1) !intermediate pressure levels (Pa) |
---|
66 | REAL, INTENT(IN) :: pplay(ngrid,nlayer) !Pressure at the middle of the layers (Pa) |
---|
67 | REAL, INTENT(IN) :: pphi(ngrid,nlayer) !Geopotential at the middle of the layers (m2s-2) |
---|
68 | REAL, INTENT(IN) :: pu(ngrid,nlayer),pv(ngrid,nlayer) !u,v components of the wind (ms-1) |
---|
69 | REAL, INTENT(IN) :: pt(ngrid,nlayer),pq(ngrid,nlayer,nq)!temperature (K) and tracer concentration (kg/kg) |
---|
70 | REAL, INTENT(IN) :: zzlay(ngrid,nlayer) ! altitude at the middle of the layers |
---|
71 | REAL, INTENT(IN) :: zzlev(ngrid,nlayer+1) ! altitude at layer boundaries |
---|
72 | INTEGER, INTENT(IN) :: igcm_co2 ! index of the CO2 tracer in mixing ratio array |
---|
73 | ! --> 0 if no tracer is CO2 (or no tracer at all) |
---|
74 | ! --> this prepares special treatment for polar night mixing |
---|
75 | ! (see thermcell_main_mars) |
---|
76 | REAL, INTENT(IN) :: pdu(ngrid,nlayer),pdv(ngrid,nlayer) ! wind velocity change from routines called |
---|
77 | ! before thermals du/dt (m/s/s) |
---|
78 | REAL, INTENT(IN) :: pdq(ngrid,nlayer,nq) ! tracer concentration change from routines called |
---|
79 | ! before thermals dq/dt (kg/kg/s) |
---|
80 | REAL, INTENT(IN) :: pdt(ngrid,nlayer) ! temperature change from routines called before thermals dT/dt (K/s) |
---|
81 | REAL, INTENT(IN) :: q2(ngrid,nlayer+1) ! turbulent kinetic energy |
---|
82 | REAL, INTENT(IN) :: zpopsk(ngrid,nlayer) ! ratio of pressure at middle of layer to surface pressure, |
---|
83 | ! to the power r/cp, i.e. zpopsk=(pplay(ig,l)/pplev(ig,1))**rcp |
---|
84 | REAL, INTENT(IN) :: pdhdif(ngrid,nlayer) ! potential temperature change from turbulent diffusion scheme dT/dt (K/s) |
---|
85 | REAL, INTENT(IN) :: sensibFlux(ngrid) ! sensible heat flux computed from surface layer scheme |
---|
86 | |
---|
87 | !-------------------------------------------------------- |
---|
88 | ! Output Variables |
---|
89 | !-------------------------------------------------------- |
---|
90 | |
---|
91 | REAL, INTENT(OUT) :: pdu_th(ngrid,nlayer) ! wind velocity change from thermals du/dt (m/s/s) |
---|
92 | REAL, INTENT(OUT) :: pdv_th(ngrid,nlayer) ! wind velocity change from thermals dv/dt (m/s/s) |
---|
93 | REAL, INTENT(OUT) :: pdt_th(ngrid,nlayer) ! temperature change from thermals dT/dt (K/s) |
---|
94 | REAL, INTENT(OUT) :: pdq_th(ngrid,nlayer,nq) ! tracer change from thermals dq/dt (kg/kg/s) |
---|
95 | INTEGER, INTENT(OUT) :: lmax(ngrid) ! layer number reacher by thermals in grid point |
---|
96 | REAL, INTENT(OUT) :: zmaxth(ngrid) ! equivalent to lmax, but in (m), interpolated |
---|
97 | REAL, INTENT(OUT) :: pbl_dtke(ngrid,nlayer+1) ! turbulent kinetic energy change from thermals dtke/dt |
---|
98 | REAL, INTENT(OUT) :: wstar(ngrid) ! free convection velocity (m/s) |
---|
99 | |
---|
100 | !-------------------------------------------------------- |
---|
101 | ! Thermals local variables |
---|
102 | !-------------------------------------------------------- |
---|
103 | REAL zu(ngrid,nlayer), zv(ngrid,nlayer) |
---|
104 | REAL zt(ngrid,nlayer) |
---|
105 | REAL d_t_ajs(ngrid,nlayer) |
---|
106 | REAL d_u_ajs(ngrid,nlayer), d_q_ajs(ngrid,nlayer,nq) |
---|
107 | REAL d_v_ajs(ngrid,nlayer) |
---|
108 | REAL fm_therm(ngrid,nlayer+1), entr_therm(ngrid,nlayer) |
---|
109 | REAL detr_therm(ngrid,nlayer),detrmod(ngrid,nlayer) |
---|
110 | REAL zw2(ngrid,nlayer+1) |
---|
111 | REAL fraca(ngrid,nlayer+1),zfraca(ngrid,nlayer+1) |
---|
112 | REAL q_therm(ngrid,nlayer), pq_therm(ngrid,nlayer,nq) |
---|
113 | REAL q2_therm(ngrid,nlayer), dq2_therm(ngrid,nlayer) |
---|
114 | REAL lmax_real(ngrid) |
---|
115 | REAL masse(ngrid,nlayer) |
---|
116 | |
---|
117 | INTEGER l,ig,iq,ii(1),k |
---|
118 | CHARACTER (LEN=20) modname |
---|
119 | |
---|
120 | !-------------------------------------------------------- |
---|
121 | ! Local variables for sub-timestep |
---|
122 | !-------------------------------------------------------- |
---|
123 | |
---|
124 | REAL d_t_the(ngrid,nlayer), d_q_the(ngrid,nlayer,nq) |
---|
125 | INTEGER isplit |
---|
126 | REAL fact |
---|
127 | REAL zfm_therm(ngrid,nlayer+1),zdt |
---|
128 | REAL zentr_therm(ngrid,nlayer),zdetr_therm(ngrid,nlayer) |
---|
129 | REAL zheatFlux(ngrid,nlayer) |
---|
130 | REAL zheatFlux_down(ngrid,nlayer) |
---|
131 | REAL zbuoyancyOut(ngrid,nlayer) |
---|
132 | REAL zbuoyancyEst(ngrid,nlayer) |
---|
133 | REAL zzw2(ngrid,nlayer+1) |
---|
134 | REAL zmax(ngrid) |
---|
135 | INTEGER ndt,limz |
---|
136 | |
---|
137 | !-------------------------------------------------------- |
---|
138 | ! Diagnostics |
---|
139 | !-------------------------------------------------------- |
---|
140 | |
---|
141 | REAL heatFlux(ngrid,nlayer) |
---|
142 | REAL heatFlux_down(ngrid,nlayer) |
---|
143 | REAL buoyancyOut(ngrid,nlayer) |
---|
144 | REAL buoyancyEst(ngrid,nlayer) |
---|
145 | REAL hfmax(ngrid),wmax(ngrid) |
---|
146 | REAL pbl_teta(ngrid),dteta(ngrid,nlayer) |
---|
147 | REAL rpdhd(ngrid,nlayer) |
---|
148 | REAL wtdif(ngrid,nlayer),rho(ngrid,nlayer) |
---|
149 | REAL wtth(ngrid,nlayer) |
---|
150 | |
---|
151 | ! ********************************************************************** |
---|
152 | ! Initializations |
---|
153 | ! ********************************************************************** |
---|
154 | |
---|
155 | lmax(:)=0 |
---|
156 | pdu_th(:,:)=0. |
---|
157 | pdv_th(:,:)=0. |
---|
158 | pdt_th(:,:)=0. |
---|
159 | entr_therm(:,:)=0. |
---|
160 | detr_therm(:,:)=0. |
---|
161 | q2_therm(:,:)=0. |
---|
162 | dq2_therm(:,:)=0. |
---|
163 | pbl_dtke(:,:)=0. |
---|
164 | fm_therm(:,:)=0. |
---|
165 | zw2(:,:)=0. |
---|
166 | fraca(:,:)=0. |
---|
167 | zfraca(:,:)=0. |
---|
168 | pdq_th(:,:,:)=0. |
---|
169 | d_t_ajs(:,:)=0. |
---|
170 | d_u_ajs(:,:)=0. |
---|
171 | d_v_ajs(:,:)=0. |
---|
172 | d_q_ajs(:,:,:)=0. |
---|
173 | heatFlux(:,:)=0. |
---|
174 | heatFlux_down(:,:)=0. |
---|
175 | buoyancyOut(:,:)=0. |
---|
176 | buoyancyEst(:,:)=0. |
---|
177 | zmaxth(:)=0. |
---|
178 | lmax_real(:)=0. |
---|
179 | |
---|
180 | |
---|
181 | ! ********************************************************************** |
---|
182 | ! Preparing inputs for the thermals: increment tendancies |
---|
183 | ! from other subroutines called before the thermals model |
---|
184 | ! ********************************************************************** |
---|
185 | |
---|
186 | zu(:,:)=pu(:,:)+pdu(:,:)*ptimestep ! u-component of wind velocity |
---|
187 | zv(:,:)=pv(:,:)+pdv(:,:)*ptimestep ! v-component of wind velocity |
---|
188 | zt(:,:)=pt(:,:)+pdt(:,:)*ptimestep ! temperature |
---|
189 | |
---|
190 | pq_therm(:,:,:)=0. ! tracer concentration |
---|
191 | |
---|
192 | if(qtransport_thermals) then |
---|
193 | pq_therm(:,:,:)=pq(:,:,:)+pdq(:,:,:)*ptimestep ! tracer concentration |
---|
194 | endif |
---|
195 | |
---|
196 | |
---|
197 | IF(dtke_thermals) THEN |
---|
198 | DO l=1,nlayer |
---|
199 | q2_therm(:,l)=0.5*(q2(:,l)+q2(:,l+1)) |
---|
200 | ENDDO |
---|
201 | ENDIF |
---|
202 | |
---|
203 | ! ********************************************************************** |
---|
204 | ! --> CALLTHERM |
---|
205 | ! SUB-TIMESTEP LOOP |
---|
206 | ! ********************************************************************** |
---|
207 | |
---|
208 | zdt=ptimestep/REAL(nsplit_thermals) !subtimestep |
---|
209 | |
---|
210 | DO isplit=1,nsplit_thermals !temporal loop on the subtimestep |
---|
211 | |
---|
212 | ! Initialization of intermediary variables |
---|
213 | |
---|
214 | zzw2(:,:)=0. |
---|
215 | zmax(:)=0. |
---|
216 | lmax(:)=0 |
---|
217 | |
---|
218 | if (nq .ne. 0 .and. igcm_co2 .ne. 0) then !initialize co2 tracer tendancy |
---|
219 | d_q_the(:,:,igcm_co2)=0. |
---|
220 | endif |
---|
221 | |
---|
222 | ! CALL to main thermal routine |
---|
223 | CALL thermcell_main_mars(ngrid,nlayer,nq & |
---|
224 | & ,igcm_co2 & |
---|
225 | & ,zdt & |
---|
226 | & ,pplay,pplev,pphi,zzlev,zzlay & |
---|
227 | & ,zu,zv,zt,pq_therm,q2_therm & |
---|
228 | & ,d_t_the,d_q_the & |
---|
229 | & ,zfm_therm,zentr_therm,zdetr_therm,lmax,zmax,limz & |
---|
230 | & ,zzw2,fraca,zpopsk & |
---|
231 | & ,zheatFlux,zheatFlux_down & |
---|
232 | & ,zbuoyancyOut,zbuoyancyEst) |
---|
233 | |
---|
234 | fact=1./REAL(nsplit_thermals) |
---|
235 | |
---|
236 | ! Update thermals tendancies |
---|
237 | d_t_the(:,:)=d_t_the(:,:)*ptimestep*fact !temperature |
---|
238 | if (igcm_co2 .ne. 0) then |
---|
239 | d_q_the(:,:,igcm_co2)=d_q_the(:,:,igcm_co2)*ptimestep*fact !co2 mass mixing ratio |
---|
240 | endif |
---|
241 | zmaxth(:)=zmaxth(:)+zmax(:)*fact !thermals height |
---|
242 | lmax_real(:)=lmax_real(:)+float(lmax(:))*fact !thermals height index |
---|
243 | fm_therm(:,:)=fm_therm(:,:) & !upward mass flux |
---|
244 | & +zfm_therm(:,:)*fact |
---|
245 | entr_therm(:,:)=entr_therm(:,:) & !entrainment mass flux |
---|
246 | & +zentr_therm(:,:)*fact |
---|
247 | detr_therm(:,:)=detr_therm(:,:) & !detrainment mass flux |
---|
248 | & +zdetr_therm(:,:)*fact |
---|
249 | zfraca(:,:)=zfraca(:,:) + fraca(:,:)*fact !updraft fractional coverage |
---|
250 | heatFlux(:,:)=heatFlux(:,:) & !upward heat flux |
---|
251 | & +zheatFlux(:,:)*fact |
---|
252 | heatFlux_down(:,:)=heatFlux_down(:,:) & !downward heat flux |
---|
253 | & +zheatFlux_down(:,:)*fact |
---|
254 | buoyancyOut(:,:)=buoyancyOut(:,:) & !plume final buoyancy |
---|
255 | & +zbuoyancyOut(:,:)*fact |
---|
256 | buoyancyEst(:,:)=buoyancyEst(:,:) & !plume estimated buoyancy used for vertical velocity computation |
---|
257 | & +zbuoyancyEst(:,:)*fact |
---|
258 | zw2(:,:)=zw2(:,:) + zzw2(:,:)*fact !vertical velocity |
---|
259 | |
---|
260 | ! Save tendancies |
---|
261 | d_t_ajs(:,:)=d_t_ajs(:,:)+d_t_the(:,:) !temperature tendancy (delta T) |
---|
262 | if (igcm_co2 .ne. 0) then |
---|
263 | d_q_ajs(:,:,igcm_co2)=d_q_ajs(:,:,igcm_co2)+d_q_the(:,:,igcm_co2) !tracer tendancy (delta q) |
---|
264 | endif |
---|
265 | |
---|
266 | ! Increment temperature and co2 concentration for next pass in subtimestep loop |
---|
267 | zt(:,:) = zt(:,:) + d_t_the(:,:) !temperature |
---|
268 | if (igcm_co2 .ne. 0) then |
---|
269 | pq_therm(:,:,igcm_co2) = & |
---|
270 | & pq_therm(:,:,igcm_co2) + d_q_the(:,:,igcm_co2) !co2 tracer |
---|
271 | endif |
---|
272 | |
---|
273 | |
---|
274 | ENDDO ! isplit |
---|
275 | !**************************************************************** |
---|
276 | |
---|
277 | ! Now that we have computed total entrainment and detrainment, we can |
---|
278 | ! advect u, v, and q in thermals. (potential temperature and co2 MMR |
---|
279 | ! have already been advected in thermcell_main because they are coupled |
---|
280 | ! to the determination of the thermals caracteristics). This is done |
---|
281 | ! separatly because u,v, and q are not used in thermcell_main for |
---|
282 | ! any thermals-related computation : they are purely passive. |
---|
283 | |
---|
284 | ! mass of cells |
---|
285 | do l=1,nlayer |
---|
286 | masse(:,l)=(pplev(:,l)-pplev(:,l+1))/g |
---|
287 | enddo |
---|
288 | |
---|
289 | ! recompute detrainment mass flux from entrainment and upward mass flux |
---|
290 | ! this ensure mass flux conservation |
---|
291 | detrmod(:,:)=0. |
---|
292 | do l=1,limz |
---|
293 | do ig=1,ngrid |
---|
294 | detrmod(ig,l)=fm_therm(ig,l)-fm_therm(ig,l+1) & |
---|
295 | & +entr_therm(ig,l) |
---|
296 | if (detrmod(ig,l).lt.0.) then |
---|
297 | entr_therm(ig,l)=entr_therm(ig,l)-detrmod(ig,l) |
---|
298 | detrmod(ig,l)=0. |
---|
299 | endif |
---|
300 | enddo |
---|
301 | enddo |
---|
302 | |
---|
303 | ! u component of wind velocity advection in thermals |
---|
304 | ! result is a derivative (d_u_ajs in m/s/s) |
---|
305 | ndt=10 |
---|
306 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
307 | & ,fm_therm,entr_therm,detrmod, & |
---|
308 | & masse,zu,d_u_ajs,ndt,limz) |
---|
309 | |
---|
310 | ! v component of wind velocity advection in thermals |
---|
311 | ! result is a derivative (d_v_ajs in m/s/s) |
---|
312 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
313 | & ,fm_therm,entr_therm,detrmod, & |
---|
314 | & masse,zv,d_v_ajs,ndt,limz) |
---|
315 | |
---|
316 | ! non co2 tracers advection in thermals |
---|
317 | ! result is a derivative (d_q_ajs in kg/kg/s) |
---|
318 | |
---|
319 | if (nq .ne. 0.) then |
---|
320 | DO iq=1,nq |
---|
321 | if (iq .ne. igcm_co2) then |
---|
322 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
323 | & ,fm_therm,entr_therm,detrmod, & |
---|
324 | & masse,pq_therm(:,:,iq),d_q_ajs(:,:,iq),ndt,limz) |
---|
325 | endif |
---|
326 | ENDDO |
---|
327 | endif |
---|
328 | |
---|
329 | ! tke advection in thermals |
---|
330 | ! result is a tendancy (d_u_ajs in J) |
---|
331 | if (dtke_thermals) then |
---|
332 | call thermcell_dqup(ngrid,nlayer,ptimestep & |
---|
333 | & ,fm_therm,entr_therm,detrmod, & |
---|
334 | & masse,q2_therm,dq2_therm,ndt,limz) |
---|
335 | endif |
---|
336 | |
---|
337 | ! compute wmax for diagnostics |
---|
338 | DO ig=1,ngrid |
---|
339 | wmax(ig)=MAXVAL(zw2(ig,:)) |
---|
340 | ENDDO |
---|
341 | |
---|
342 | ! ********************************************************************** |
---|
343 | ! ********************************************************************** |
---|
344 | ! ********************************************************************** |
---|
345 | ! CALLTHERM END |
---|
346 | ! ********************************************************************** |
---|
347 | ! ********************************************************************** |
---|
348 | ! ********************************************************************** |
---|
349 | |
---|
350 | |
---|
351 | ! ********************************************************************** |
---|
352 | ! Preparing outputs |
---|
353 | ! ********************************************************************** |
---|
354 | |
---|
355 | do l=1,limz |
---|
356 | pdu_th(:,l)=d_u_ajs(:,l) |
---|
357 | pdv_th(:,l)=d_v_ajs(:,l) |
---|
358 | enddo |
---|
359 | |
---|
360 | ! if tracers are transported in thermals, update output variables, else these are 0. |
---|
361 | if(qtransport_thermals) then |
---|
362 | do iq=1,nq |
---|
363 | if (iq .ne. igcm_co2) then |
---|
364 | do l=1,limz |
---|
365 | pdq_th(:,l,iq)=d_q_ajs(:,l,iq) !non-co2 tracers d_q_ajs are dq/dt (kg/kg/s) |
---|
366 | enddo |
---|
367 | else |
---|
368 | do l=1,limz |
---|
369 | pdq_th(:,l,iq)=d_q_ajs(:,l,iq)/ptimestep !co2 tracer d_q_ajs is dq (kg/kg) |
---|
370 | enddo |
---|
371 | endif |
---|
372 | enddo |
---|
373 | endif |
---|
374 | |
---|
375 | ! if tke is transported in thermals, update output variable, else this is 0. |
---|
376 | IF(dtke_thermals) THEN |
---|
377 | DO l=2,nlayer |
---|
378 | pbl_dtke(:,l)=0.5*(dq2_therm(:,l-1)+dq2_therm(:,l)) |
---|
379 | ENDDO |
---|
380 | |
---|
381 | pbl_dtke(:,1)=0.5*dq2_therm(:,1) |
---|
382 | pbl_dtke(:,nlayer+1)=0. |
---|
383 | ENDIF |
---|
384 | |
---|
385 | ! update output variable for temperature. d_t_ajs is delta T in (K), pdt_th is dT/dt in (K/s) |
---|
386 | do l=1,limz |
---|
387 | pdt_th(:,l)=d_t_ajs(:,l)/ptimestep |
---|
388 | enddo |
---|
389 | |
---|
390 | |
---|
391 | ! ********************************************************************** |
---|
392 | ! SURFACE LAYER INTERFACE |
---|
393 | ! Compute the free convection velocity w* scale for surface layer gustiness |
---|
394 | ! speed parameterization. The computed value of w* will be used at the next |
---|
395 | ! timestep to modify surface-atmosphere exchange fluxes, because the surface |
---|
396 | ! layer scheme and diffusion are called BEFORE the thermals. (outside of these |
---|
397 | ! routines) |
---|
398 | ! ********************************************************************** |
---|
399 | |
---|
400 | ! Potential temperature gradient |
---|
401 | |
---|
402 | dteta(:,nlayer)=0. |
---|
403 | DO l=1,nlayer-1 |
---|
404 | DO ig=1, ngrid |
---|
405 | dteta(ig,l) = ((zt(ig,l+1)-zt(ig,l))/zpopsk(ig,l)) & |
---|
406 | & /(zzlay(ig,l+1)-zzlay(ig,l)) |
---|
407 | ENDDO |
---|
408 | ENDDO |
---|
409 | |
---|
410 | ! Computation of the PBL mixed layer temperature |
---|
411 | |
---|
412 | DO ig=1, ngrid |
---|
413 | ii=MINLOC(abs(dteta(ig,1:lmax(ig)))) |
---|
414 | pbl_teta(ig) = zt(ig,ii(1))/zpopsk(ig,ii(1)) |
---|
415 | ENDDO |
---|
416 | |
---|
417 | ! In order to have an accurate w*, we must add the heat flux from the |
---|
418 | ! diffusion scheme to the computation of the maximum heat flux hfmax |
---|
419 | ! Here pdhdif is the potential temperature change from the diffusion |
---|
420 | ! scheme (Mellor and Yamada, see paper section 6, paragraph 57) |
---|
421 | |
---|
422 | ! compute rho as it is after the diffusion |
---|
423 | |
---|
424 | rho(:,:)=pplay(:,:) & |
---|
425 | & /(r*(pt(:,:)+pdhdif(:,:)*zpopsk(:,:)*ptimestep)) |
---|
426 | |
---|
427 | ! integrate -rho*pdhdif |
---|
428 | |
---|
429 | rpdhd(:,:)=0. |
---|
430 | |
---|
431 | DO ig=1,ngrid |
---|
432 | DO l=1,lmax(ig) |
---|
433 | rpdhd(ig,l)=0. |
---|
434 | DO k=1,l |
---|
435 | rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*pdhdif(ig,k)* & |
---|
436 | & (zzlev(ig,k+1)-zzlev(ig,k)) |
---|
437 | ENDDO |
---|
438 | rpdhd(ig,l)=rpdhd(ig,l)-sensibFlux(ig)/cpp |
---|
439 | ENDDO |
---|
440 | ENDDO |
---|
441 | |
---|
442 | ! compute w'theta' (vertical turbulent flux of temperature) from |
---|
443 | ! the diffusion scheme |
---|
444 | |
---|
445 | wtdif(:,:)=rpdhd(:,:)/rho(:,:) |
---|
446 | |
---|
447 | ! Now we compute the contribution of the thermals to w'theta': |
---|
448 | ! compute rho as it is after the thermals |
---|
449 | |
---|
450 | rho(:,:)=pplay(:,:)/(r*(zt(:,:))) |
---|
451 | |
---|
452 | ! integrate -rho*pdhdif |
---|
453 | |
---|
454 | DO ig=1,ngrid |
---|
455 | DO l=1,lmax(ig) |
---|
456 | rpdhd(ig,l)=0. |
---|
457 | DO k=1,l |
---|
458 | rpdhd(ig,l)=rpdhd(ig,l)-rho(ig,k)*(pdt_th(ig,k)/zpopsk(ig,k))* & |
---|
459 | & (zzlev(ig,k+1)-zzlev(ig,k)) |
---|
460 | ENDDO |
---|
461 | rpdhd(ig,l)=rpdhd(ig,l)+ & |
---|
462 | & rho(ig,1)*(heatFlux(ig,1)+heatFlux_down(ig,1)) |
---|
463 | ENDDO |
---|
464 | ENDDO |
---|
465 | rpdhd(:,nlayer)=0. |
---|
466 | |
---|
467 | ! compute w'teta' from thermals |
---|
468 | |
---|
469 | wtth(:,:)=rpdhd(:,:)/rho(:,:) |
---|
470 | |
---|
471 | ! Add vertical turbulent heat fluxes from the thermals and the diffusion scheme |
---|
472 | ! and compute the maximum |
---|
473 | |
---|
474 | DO ig=1,ngrid |
---|
475 | hfmax(ig)=MAXVAL(wtth(ig,:)+wtdif(ig,:)) |
---|
476 | ENDDO |
---|
477 | |
---|
478 | ! Finally we can compute the free convection velocity scale |
---|
479 | ! We follow Spiga et. al 2010 (QJRMS) |
---|
480 | ! ------------ |
---|
481 | |
---|
482 | DO ig=1, ngrid |
---|
483 | IF (zmax(ig) .gt. 0.) THEN |
---|
484 | wstar(ig)=(g*zmaxth(ig)*hfmax(ig)/pbl_teta(ig))**(1./3.) |
---|
485 | ELSE |
---|
486 | wstar(ig)=0. |
---|
487 | ENDIF |
---|
488 | ENDDO |
---|
489 | |
---|
490 | END |
---|