1 | SUBROUTINE mass_redistribution(ngrid,nlayer,nq,ptimestep, & |
---|
2 | rnat,pcapcal,pplay,pplev,pt,ptsrf,pq,pqs, & |
---|
3 | pu,pv,pdt,pdtsrf,pdq,pdu,pdv,pdmassmr, & |
---|
4 | pdtmr,pdtsrfmr,pdpsrfmr,pdumr,pdvmr,pdqmr,pdqsmr) |
---|
5 | |
---|
6 | USE watercommon_h, Only: Tsat_water,RLVTT |
---|
7 | use surfdat_h |
---|
8 | use radcommon_h, only: glat |
---|
9 | USE tracer_h |
---|
10 | USE planete_mod, only: bp |
---|
11 | use comcstfi_mod, only: g |
---|
12 | USE callkeys_mod, ONLY: water |
---|
13 | |
---|
14 | IMPLICIT NONE |
---|
15 | !======================================================================= |
---|
16 | ! subject: |
---|
17 | ! -------- |
---|
18 | ! Mass and momentum fluxes through sigma levels as the surface pressure is modified are also taken into account |
---|
19 | ! |
---|
20 | ! author: Jeremy Leconte 2012 (from F.Forget 1998) |
---|
21 | ! ------ |
---|
22 | ! |
---|
23 | ! input: |
---|
24 | ! ----- |
---|
25 | ! ngrid nombre de points de verticales |
---|
26 | ! (toutes les boucles de la physique sont au |
---|
27 | ! moins vectorisees sur ngrid) |
---|
28 | ! nlayer nombre de couches |
---|
29 | ! pplay(ngrid,nlayer) Pressure levels |
---|
30 | ! pplev(ngrid,nlayer+1) Pressure levels |
---|
31 | ! nq Number of tracers |
---|
32 | ! |
---|
33 | ! pt(ngrid,nlayer) temperature (en K) |
---|
34 | ! pq(ngrid,nlayer,nq) tracer specific concentration (kg/kg of air) |
---|
35 | ! pu,pv (ngrid,nlayer) wind velocity (m/s) |
---|
36 | ! |
---|
37 | ! |
---|
38 | ! pdX physical tendency of X before mass redistribution |
---|
39 | ! |
---|
40 | ! pdmassmr air Mass added to the atmosphere in each layer (kg/m2/s) |
---|
41 | ! |
---|
42 | ! output: |
---|
43 | ! ------- |
---|
44 | ! |
---|
45 | ! pdXmr(ngrid) physical tendency of X after mass redistribution |
---|
46 | ! |
---|
47 | ! |
---|
48 | ! |
---|
49 | !======================================================================= |
---|
50 | ! |
---|
51 | ! 0. Declarations : |
---|
52 | ! ------------------ |
---|
53 | |
---|
54 | !----------------------------------------------------------------------- |
---|
55 | ! Arguments : |
---|
56 | ! --------- |
---|
57 | INTEGER ngrid, nlayer, nq |
---|
58 | REAL ptimestep |
---|
59 | REAL pcapcal(ngrid) |
---|
60 | INTEGER rnat(ngrid) |
---|
61 | REAL pplay(ngrid,nlayer),pplev(ngrid,nlayer+1) |
---|
62 | REAL pt(ngrid,nlayer),pdt(ngrid,nlayer) |
---|
63 | REAL ptsrf(ngrid),pdtsrf(ngrid) |
---|
64 | REAL pdtmr(ngrid,nlayer) |
---|
65 | REAL pu(ngrid,nlayer) , pv(ngrid,nlayer) |
---|
66 | REAL pdu(ngrid,nlayer) , pdv(ngrid,nlayer) |
---|
67 | REAL pdmassmr(ngrid,nlayer) |
---|
68 | REAL pdumr(ngrid,nlayer) , pdvmr(ngrid,nlayer) |
---|
69 | REAL pq(ngrid,nlayer,nq),pdq(ngrid,nlayer,nq) |
---|
70 | REAL pqs(ngrid,nq) |
---|
71 | REAL pdqmr(ngrid,nlayer,nq),pdqsmr(ngrid,nq) |
---|
72 | REAL pdpsrfmr(ngrid),pdtsrfmr(ngrid) |
---|
73 | ! |
---|
74 | ! Local variables : |
---|
75 | ! ----------------- |
---|
76 | |
---|
77 | ! Boiling/sublimation |
---|
78 | REAL Tsat(ngrid),zmassboil(ngrid) |
---|
79 | |
---|
80 | ! vertical reorganization of sigma levels |
---|
81 | REAL zzu(nlayer),zzv(nlayer) |
---|
82 | REAL zzq(nlayer,nq),zzt(nlayer) |
---|
83 | ! Dummy variables |
---|
84 | INTEGER n,l,ig,iq |
---|
85 | REAL zdtsig(ngrid,nlayer) |
---|
86 | REAL zmass(ngrid,nlayer),zzmass(nlayer),w(nlayer+1) |
---|
87 | REAL zdmass_sum(ngrid,nlayer+1) |
---|
88 | REAL zmflux(nlayer+1) |
---|
89 | REAL zq1(nlayer) |
---|
90 | REAL ztsrf(ngrid) |
---|
91 | REAL ztm(nlayer+1) |
---|
92 | REAL zum(nlayer+1) , zvm(nlayer+1) |
---|
93 | REAL zqm(nlayer+1,nq),zqm1(nlayer+1) |
---|
94 | REAL sigma(nlayer+1) |
---|
95 | |
---|
96 | ! local saved variables |
---|
97 | LOGICAL, SAVE :: firstcall=.true. |
---|
98 | !$OMP THREADPRIVATE(firstcall) |
---|
99 | |
---|
100 | !---------------------------------------------------------------------- |
---|
101 | |
---|
102 | ! Initialisation |
---|
103 | ! -------------- |
---|
104 | ! |
---|
105 | IF (firstcall) THEN |
---|
106 | firstcall=.false. |
---|
107 | ENDIF |
---|
108 | ! |
---|
109 | !====================================================================== |
---|
110 | ! Calcul of h2o condensation |
---|
111 | ! ============================================================ |
---|
112 | ! |
---|
113 | ! Used variable : |
---|
114 | ! pdmassmr : air Mass added to the atmosphere in each layer per unit time (kg/m2/s) |
---|
115 | ! zdmass_sum(ngrid,l) : total air mass added to the atm above layer l per unit time (kg/m2/s) |
---|
116 | ! |
---|
117 | ! |
---|
118 | ! Surface tracer Tendencies set to 0 |
---|
119 | ! ------------------------------------- |
---|
120 | pdqsmr(1:ngrid,1:nq)=0. |
---|
121 | |
---|
122 | ztsrf(1:ngrid) = ptsrf(1:ngrid) + pdtsrf(1:ngrid)*ptimestep |
---|
123 | |
---|
124 | |
---|
125 | DO ig=1,ngrid |
---|
126 | zdmass_sum(ig,nlayer+1)=0. |
---|
127 | DO l = nlayer, 1, -1 |
---|
128 | zmass(ig,l) = (pplev(ig,l)-pplev(ig,l+1))/glat(ig) |
---|
129 | zdmass_sum(ig,l)= zdmass_sum(ig,l+1)+pdmassmr(ig,l) |
---|
130 | END DO |
---|
131 | END DO |
---|
132 | |
---|
133 | |
---|
134 | if (water) then |
---|
135 | do ig=1,ngrid |
---|
136 | call Tsat_water(pplev(ig,1)+zdmass_sum(ig,1)*g*ptimestep,Tsat(ig)) |
---|
137 | enddo |
---|
138 | #ifndef MESOSCALE |
---|
139 | call writediagfi(ngrid,'Tsat','saturation temperature at surface','',2,Tsat) |
---|
140 | #endif |
---|
141 | |
---|
142 | do ig=1,ngrid |
---|
143 | if (ztsrf(ig).gt.Tsat(ig)) then |
---|
144 | zmassboil(ig)=(ptsrf(ig)-Tsat(ig))*pcapcal(ig)/RLVTT/ptimestep |
---|
145 | if ((zmassboil(ig)*ptimestep.gt.pqs(ig,igcm_h2o_vap)).and.(rnat(ig).eq.1)) then |
---|
146 | zmassboil(ig)=pqs(ig,igcm_h2o_vap)/ptimestep |
---|
147 | endif |
---|
148 | zmassboil(ig)=zmassboil(ig)*0.0 !momentary, should be 1. JL12 |
---|
149 | pdqsmr(ig,igcm_h2o_vap)=-zmassboil(ig) |
---|
150 | pdtsrfmr(ig)=-zmassboil(ig)*RLVTT/pcapcal(ig) |
---|
151 | ztsrf(ig)=ptsrf(ig)+pdtsrfmr(ig)*ptimestep |
---|
152 | else |
---|
153 | zmassboil(ig)=0. |
---|
154 | pdtsrfmr(ig)=0. |
---|
155 | endif |
---|
156 | enddo |
---|
157 | endif |
---|
158 | |
---|
159 | ! ************************* |
---|
160 | ! UPDATE SURFACE |
---|
161 | ! ************************* |
---|
162 | ! Changing pressure at the surface: |
---|
163 | ! """""""""""""""""""""""""""""""""""" |
---|
164 | |
---|
165 | pdpsrfmr(1:ngrid) = (zdmass_sum(1:ngrid,1)+zmassboil(1:ngrid))*g |
---|
166 | |
---|
167 | do ig = 1, ngrid |
---|
168 | IF(ABS(pdpsrfmr(ig)*ptimestep).GT.pplev(ig,1)) THEN |
---|
169 | PRINT*,'STOP in condens in mass_redistribution' |
---|
170 | PRINT*,'condensing more than total mass' |
---|
171 | PRINT*,'Grid point ',ig |
---|
172 | PRINT*,'Ps = ',pplev(ig,1) |
---|
173 | PRINT*,'d Ps = ',pdpsrfmr(ig)*ptimestep |
---|
174 | STOP |
---|
175 | ENDIF |
---|
176 | enddo ! of DO ig=1,ngrid |
---|
177 | |
---|
178 | |
---|
179 | ! *************************************************************** |
---|
180 | ! Correction to account for redistribution between sigma or hybrid |
---|
181 | ! layers when changing surface pressure |
---|
182 | ! zzx quantities have dimension (nlayer) to speed up calculation |
---|
183 | ! ************************************************************* |
---|
184 | |
---|
185 | DO ig=1,ngrid |
---|
186 | zzt(1:nlayer) = pt(ig,1:nlayer) + pdt(ig,1:nlayer) * ptimestep |
---|
187 | zzu(1:nlayer) = pu(ig,1:nlayer) + pdu(ig,1:nlayer) * ptimestep |
---|
188 | zzv(1:nlayer) = pv(ig,1:nlayer) + pdv(ig,1:nlayer) * ptimestep |
---|
189 | zzq(1:nlayer,1:nq)=pq(ig,1:nlayer,1:nq)+pdq(ig,1:nlayer,1:nq)*ptimestep ! must add the water that has fallen??? |
---|
190 | |
---|
191 | ! Mass fluxes of air through the sigma levels (kg.m-2.s-1) (>0 when up) |
---|
192 | ! """""""""""""""""""""""""""""""""""""""""""""""""""""""""""""" |
---|
193 | zmflux(1) = zmassboil(ig) |
---|
194 | sigma(1)=1 |
---|
195 | DO l=1,nlayer |
---|
196 | ! Ehouarn: shouldn't we rather compute sigma levels than use bp()? |
---|
197 | ! sigma(l+1)=pplev(ig,l+1)/pplev(ig,1) |
---|
198 | ! zmflux(l+1) = zmflux(l) + pdmassmr(ig,l) - & |
---|
199 | ! (sigma(l)-sigma(l+1))*(zdmass_sum(ig,1)+zmflux(1)) |
---|
200 | ! if (abs(zmflux(l+1)).lt.1E-13.OR.sigma(l+1).eq.0.) zmflux(l+1)=0. |
---|
201 | ! Ehouarn: but for now leave things as before |
---|
202 | zmflux(l+1) = zmflux(l) + pdmassmr(ig,l) - (bp(l)-bp(l+1))*(zdmass_sum(ig,1)+zmflux(1)) |
---|
203 | ! zmflux set to 0 if very low to avoid: top layer is disappearing in v1ld |
---|
204 | if (abs(zmflux(l+1)).lt.1E-13.OR.bp(l+1).eq.0.) zmflux(l+1)=0. |
---|
205 | END DO |
---|
206 | |
---|
207 | ! Mass of each layer |
---|
208 | ! ------------------ |
---|
209 | zzmass(1:nlayer)=zmass(ig,1:nlayer)*(1.+pdpsrfmr(ig)*ptimestep/pplev(ig,1)) |
---|
210 | |
---|
211 | |
---|
212 | ! Corresponding fluxes for T,U,V,Q |
---|
213 | ! """""""""""""""""""""""""""""""" |
---|
214 | |
---|
215 | ! averaging operator for TRANSPORT |
---|
216 | ! """""""""""""""""""""""""""""""" |
---|
217 | ! Value transfert at the surface interface when condensation |
---|
218 | ! sublimation: |
---|
219 | ztm(1) = ztsrf(ig) |
---|
220 | zum(1) = 0. |
---|
221 | zvm(1) = 0. |
---|
222 | zqm(1,1:nq)=0. ! most tracer do not condense ! |
---|
223 | if (water) zqm(1,igcm_h2o_vap)=1. ! flux is 100% h2o at surface |
---|
224 | |
---|
225 | ! Van Leer scheme: |
---|
226 | w(1:nlayer+1)=-zmflux(1:nlayer+1)*ptimestep |
---|
227 | call vl1d(nlayer,zzt,2.,zzmass,w,ztm) |
---|
228 | call vl1d(nlayer,zzu,2.,zzmass,w,zum) |
---|
229 | call vl1d(nlayer,zzv,2.,zzmass,w,zvm) |
---|
230 | do iq=1,nq |
---|
231 | zq1(1:nlayer)=zzq(1:nlayer,iq) |
---|
232 | zqm1(1)=zqm(1,iq) |
---|
233 | ! print*,iq |
---|
234 | ! print*,zq1 |
---|
235 | call vl1d(nlayer,zq1,2.,zzmass,w,zqm1) |
---|
236 | do l=2,nlayer |
---|
237 | zzq(l,iq)=zq1(l) |
---|
238 | zqm(l,iq)=zqm1(l) |
---|
239 | enddo |
---|
240 | enddo |
---|
241 | |
---|
242 | ! Surface condensation affects low winds |
---|
243 | if (zmflux(1).lt.0) then |
---|
244 | zum(1)= zzu(1) * (w(1)/zzmass(1)) |
---|
245 | zvm(1)= zzv(1) * (w(1)/zzmass(1)) |
---|
246 | if (w(1).gt.zzmass(1)) then ! ensure numerical stability |
---|
247 | zum(1)= (zzu(1)-zum(2))*zzmass(1)/w(1) + zum(2) |
---|
248 | zvm(1)= (zzv(1)-zvm(2))*zzmass(1)/w(1) + zvm(2) |
---|
249 | end if |
---|
250 | end if |
---|
251 | |
---|
252 | ztm(nlayer+1)= zzt(nlayer) ! should not be used, but... |
---|
253 | zum(nlayer+1)= zzu(nlayer) ! should not be used, but... |
---|
254 | zvm(nlayer+1)= zzv(nlayer) ! should not be used, but... |
---|
255 | zqm(nlayer+1,1:nq)= zzq(nlayer,1:nq) |
---|
256 | |
---|
257 | ! Tendencies on T, U, V, Q |
---|
258 | ! """""""""""""""""""""""" |
---|
259 | DO l=1,nlayer |
---|
260 | |
---|
261 | ! Tendencies on T |
---|
262 | pdtmr(ig,l) = (1/zzmass(l)) * & |
---|
263 | (zmflux(l)*(ztm(l) - zzt(l))-zmflux(l+1)*(ztm(l+1)-zzt(l))) |
---|
264 | !JL12 the last term in Newcondens has been set to zero because we are only dealing with redistribution here |
---|
265 | |
---|
266 | ! Tendencies on U |
---|
267 | pdumr(ig,l) = (1/zzmass(l)) *( zmflux(l)*(zum(l) - zzu(l)) - zmflux(l+1)*(zum(l+1) - zzu(l)) ) |
---|
268 | |
---|
269 | ! Tendencies on V |
---|
270 | pdvmr(ig,l) = (1/zzmass(l)) *( zmflux(l)*(zvm(l) - zzv(l)) - zmflux(l+1)*(zvm(l+1) - zzv(l)) ) |
---|
271 | |
---|
272 | END DO |
---|
273 | |
---|
274 | ! Tendencies on Q |
---|
275 | do iq=1,nq |
---|
276 | DO l=1,nlayer |
---|
277 | pdqmr(ig,l,iq)= (1/zzmass(l)) * & |
---|
278 | (zmflux(l)*(zqm(l,iq)-zzq(l,iq))- zmflux(l+1)*(zqm(l+1,iq)-zzq(l,iq)) - pdmassmr(ig,l)*zzq(l,iq)) |
---|
279 | END DO |
---|
280 | enddo |
---|
281 | |
---|
282 | END DO ! loop on ig |
---|
283 | |
---|
284 | CONTAINS |
---|
285 | |
---|
286 | ! ***************************************************************** |
---|
287 | SUBROUTINE vl1d(llm,q,pente_max,zzmass,w,qm) |
---|
288 | ! |
---|
289 | ! |
---|
290 | ! Operateur de moyenne inter-couche pour calcul de transport type |
---|
291 | ! Van-Leer " pseudo amont " dans la verticale |
---|
292 | ! q,w sont des arguments d'entree pour le s-pg .... |
---|
293 | ! masse : masse de la couche Dp/g |
---|
294 | ! w : masse d'atm ``transferee'' a chaque pas de temps (kg.m-2) |
---|
295 | ! pente_max = 2 conseillee |
---|
296 | ! |
---|
297 | ! |
---|
298 | ! -------------------------------------------------------------------- |
---|
299 | |
---|
300 | IMPLICIT NONE |
---|
301 | |
---|
302 | ! Arguments: |
---|
303 | ! ---------- |
---|
304 | integer,intent(in) :: llm |
---|
305 | real zzmass(llm),pente_max |
---|
306 | REAL q(llm),qm(llm+1) |
---|
307 | REAL w(llm+1) |
---|
308 | ! |
---|
309 | ! Local |
---|
310 | ! --------- |
---|
311 | ! |
---|
312 | INTEGER l |
---|
313 | ! |
---|
314 | real dzq(llm),dzqw(llm),adzqw(llm),dzqmax |
---|
315 | real sigw, Mtot, MQtot |
---|
316 | integer m |
---|
317 | ! integer ismax,ismin |
---|
318 | |
---|
319 | |
---|
320 | ! On oriente tout dans le sens de la pression |
---|
321 | ! W > 0 WHEN DOWN !!!!!!!!!!!!! |
---|
322 | |
---|
323 | do l=2,llm |
---|
324 | dzqw(l)=q(l-1)-q(l) |
---|
325 | adzqw(l)=abs(dzqw(l)) |
---|
326 | enddo |
---|
327 | |
---|
328 | do l=2,llm-1 |
---|
329 | if(dzqw(l)*dzqw(l+1).gt.0.) then |
---|
330 | dzq(l)=0.5*(dzqw(l)+dzqw(l+1)) |
---|
331 | else |
---|
332 | dzq(l)=0. |
---|
333 | endif |
---|
334 | dzqmax=pente_max*min(adzqw(l),adzqw(l+1)) |
---|
335 | dzq(l)=sign(min(abs(dzq(l)),dzqmax),dzq(l)) |
---|
336 | enddo |
---|
337 | |
---|
338 | dzq(1)=0. |
---|
339 | dzq(llm)=0. |
---|
340 | |
---|
341 | do l = 1,llm-1 |
---|
342 | |
---|
343 | ! Regular scheme (transfered mass < layer mass) |
---|
344 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
345 | if(w(l+1).gt.0. .and. w(l+1).le.zzmass(l+1)) then |
---|
346 | sigw=w(l+1)/zzmass(l+1) |
---|
347 | qm(l+1)=(q(l+1)+0.5*(1.-sigw)*dzq(l+1)) |
---|
348 | else if(w(l+1).le.0. .and. -w(l+1).le.zzmass(l)) then |
---|
349 | sigw=w(l+1)/zzmass(l) |
---|
350 | qm(l+1)=(q(l)-0.5*(1.+sigw)*dzq(l)) |
---|
351 | |
---|
352 | ! Extended scheme (transfered mass > layer mass) |
---|
353 | ! ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
354 | else if(w(l+1).gt.0.) then |
---|
355 | m=l+1 |
---|
356 | Mtot = zzmass(m) |
---|
357 | MQtot = zzmass(m)*q(m) |
---|
358 | do while ((m.lt.llm).and.(w(l+1).gt.(Mtot+zzmass(m+1)))) |
---|
359 | m=m+1 |
---|
360 | Mtot = Mtot + zzmass(m) |
---|
361 | MQtot = MQtot + zzmass(m)*q(m) |
---|
362 | end do |
---|
363 | if (m.lt.llm) then |
---|
364 | sigw=(w(l+1)-Mtot)/zzmass(m+1) |
---|
365 | qm(l+1)= (1/w(l+1))*(MQtot + (w(l+1)-Mtot)*(q(m+1)+0.5*(1.-sigw)*dzq(m+1)) ) |
---|
366 | else |
---|
367 | ! w(l+1) = Mtot |
---|
368 | ! qm(l+1) = Mqtot / Mtot |
---|
369 | write(*,*) 'top layer is disappearing !',l,Mtot,w(l+1),qm(l+1) |
---|
370 | print*,zzmass |
---|
371 | print*,w |
---|
372 | print*,q |
---|
373 | print*,qm |
---|
374 | stop |
---|
375 | end if |
---|
376 | else ! if(w(l+1).lt.0) |
---|
377 | m = l-1 |
---|
378 | Mtot = zzmass(m+1) |
---|
379 | MQtot = zzmass(m+1)*q(m+1) |
---|
380 | if (m.gt.0) then ! because some compilers will have problems |
---|
381 | ! evaluating zzmass(0) |
---|
382 | do while ((m.gt.0).and.(-w(l+1).gt.(Mtot+zzmass(m)))) |
---|
383 | m=m-1 |
---|
384 | Mtot = Mtot + zzmass(m+1) |
---|
385 | MQtot = MQtot + zzmass(m+1)*q(m+1) |
---|
386 | if (m.eq.0) exit |
---|
387 | end do |
---|
388 | endif |
---|
389 | if (m.gt.0) then |
---|
390 | sigw=(w(l+1)+Mtot)/zzmass(m) |
---|
391 | qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)*(q(m)-0.5*(1.+sigw)*dzq(m)) ) |
---|
392 | else |
---|
393 | qm(l+1)= (-1/w(l+1))*(MQtot + (-w(l+1)-Mtot)*qm(1)) |
---|
394 | end if |
---|
395 | end if |
---|
396 | enddo |
---|
397 | |
---|
398 | ! boundary conditions (not used in newcondens !!) |
---|
399 | ! qm(llm+1)=0. |
---|
400 | ! if(w(1).gt.0.) then |
---|
401 | ! qm(1)=q(1) |
---|
402 | ! else |
---|
403 | ! qm(1)=0. |
---|
404 | ! end if |
---|
405 | |
---|
406 | END SUBROUTINE vl1d |
---|
407 | |
---|
408 | END SUBROUTINE mass_redistribution |
---|