1 | MODULE watercloud_mod |
---|
2 | |
---|
3 | IMPLICIT NONE |
---|
4 | |
---|
5 | REAL,SAVE,ALLOCATABLE :: zdqcloud(:,:,:) ! tendencies on pq due to condensation of H2O(kg/kg.s-1) |
---|
6 | REAL,SAVE,ALLOCATABLE :: zdqscloud(:,:) ! tendencies on qsurf (calculated only by calchim but declared here) |
---|
7 | |
---|
8 | CONTAINS |
---|
9 | |
---|
10 | SUBROUTINE watercloud(ngrid,nlay,ptimestep, |
---|
11 | & pplev,pplay,pdpsrf,pzlay,pt,pdt, |
---|
12 | & pq,pdq,pdqcloud,pdtcloud, |
---|
13 | & nq,tau,tauscaling,rdust,rice,nuice, |
---|
14 | & rsedcloud,rhocloud,totcloudfrac) |
---|
15 | USE ioipsl_getin_p_mod, ONLY : getin_p |
---|
16 | USE updaterad, ONLY: updaterdust, updaterice_micro, |
---|
17 | & updaterice_typ |
---|
18 | USE improvedclouds_mod, ONLY: improvedclouds |
---|
19 | USE watersat_mod, ONLY: watersat |
---|
20 | use tracer_mod, only: nqmx, igcm_h2o_vap, igcm_h2o_ice, |
---|
21 | & igcm_hdo_vap, igcm_hdo_ice, |
---|
22 | & igcm_dust_mass, igcm_dust_number, |
---|
23 | & igcm_ccn_mass, igcm_ccn_number, |
---|
24 | & rho_dust, nuice_sed, nuice_ref |
---|
25 | use dimradmars_mod, only: naerkind |
---|
26 | IMPLICIT NONE |
---|
27 | |
---|
28 | |
---|
29 | c======================================================================= |
---|
30 | c Water-ice cloud formation |
---|
31 | c |
---|
32 | c Includes two different schemes: |
---|
33 | c - A simplified scheme (see simpleclouds.F) |
---|
34 | c - An improved microphysical scheme (see improvedclouds.F) |
---|
35 | c |
---|
36 | c There is a time loop specific to cloud formation |
---|
37 | c due to timescales smaller than the GCM integration timestep. |
---|
38 | c |
---|
39 | c Authors: Franck Montmessin, Francois Forget, Ehouarn Millour, |
---|
40 | c J.-B. Madeleine, Thomas Navarro |
---|
41 | c |
---|
42 | c 2004 - 2012 |
---|
43 | c======================================================================= |
---|
44 | |
---|
45 | c----------------------------------------------------------------------- |
---|
46 | c declarations: |
---|
47 | c ------------- |
---|
48 | |
---|
49 | include "callkeys.h" |
---|
50 | |
---|
51 | c Inputs/outputs: |
---|
52 | c ------ |
---|
53 | |
---|
54 | INTEGER, INTENT(IN) :: ngrid,nlay |
---|
55 | INTEGER, INTENT(IN) :: nq ! nombre de traceurs |
---|
56 | REAL, INTENT(IN) :: ptimestep ! pas de temps physique (s) |
---|
57 | REAL, INTENT(IN) :: pplev(ngrid,nlay+1) ! pression aux inter-couches (Pa) |
---|
58 | REAL, INTENT(IN) :: pplay(ngrid,nlay) ! pression au milieu des couches (Pa) |
---|
59 | REAL, INTENT(IN) :: pdpsrf(ngrid) ! tendence surf pressure |
---|
60 | REAL, INTENT(IN) :: pzlay(ngrid,nlay) ! altitude at the middle of the layers |
---|
61 | REAL, INTENT(IN) :: pt(ngrid,nlay) ! temperature at the middle of the layers (K) |
---|
62 | REAL, INTENT(IN) :: pdt(ngrid,nlay) ! tendence temperature des autres param. |
---|
63 | |
---|
64 | REAL, INTENT(IN) :: pq(ngrid,nlay,nq) ! traceur (kg/kg) |
---|
65 | rEAL, INTENT(IN) :: pdq(ngrid,nlay,nq) ! tendence avant condensation (kg/kg.s-1) |
---|
66 | |
---|
67 | REAL, INTENT(IN) :: tau(ngrid,naerkind) ! Column dust optical depth at each point |
---|
68 | REAL, INTENT(IN) :: tauscaling(ngrid) ! Convertion factor for dust amount |
---|
69 | REAL, INTENT(INOUT) :: rdust(ngrid,nlay) ! Dust geometric mean radius (m) |
---|
70 | |
---|
71 | REAL, INTENT(OUT) :: pdqcloud(ngrid,nlay,nq) ! tendence de la condensation H2O(kg/kg.s-1) |
---|
72 | REAL, INTENT(OUT) :: pdtcloud(ngrid,nlay) ! tendence temperature due |
---|
73 | ! a la chaleur latente |
---|
74 | REAL, INTENT(INOUT) :: rice(ngrid,nlay) ! Ice mass mean radius (m) |
---|
75 | ! (r_c in montmessin_2004) |
---|
76 | REAL, INTENT(OUT) :: nuice(ngrid,nlay) ! Estimated effective variance |
---|
77 | ! of the size distribution |
---|
78 | REAL, INTENT(OUT) :: rsedcloud(ngrid,nlay) ! Cloud sedimentation radius |
---|
79 | REAL, INTENT(OUT) :: rhocloud(ngrid,nlay) ! Cloud density (kg.m-3) |
---|
80 | |
---|
81 | REAL, INTENT(INOUT):: totcloudfrac(ngrid) ! Cloud fraction (A. Pottier 2013) |
---|
82 | |
---|
83 | c Locals: |
---|
84 | c ------ |
---|
85 | |
---|
86 | ! for ice radius computation |
---|
87 | REAL Mo,No |
---|
88 | REAl ccntyp |
---|
89 | |
---|
90 | ! for time loop |
---|
91 | INTEGER microstep ! time subsampling step variable |
---|
92 | INTEGER,SAVE :: imicro ! time subsampling for coupled water microphysics & sedimentation |
---|
93 | REAL,SAVE :: microtimestep ! integration timestep for coupled water microphysics & sedimentation |
---|
94 | REAL,SAVE :: microtimestep_prev=-999 |
---|
95 | |
---|
96 | ! tendency given by clouds (inside the micro loop) |
---|
97 | REAL subpdqcloud(ngrid,nlay,nq) ! cf. pdqcloud |
---|
98 | REAL subpdtcloud(ngrid,nlay) ! cf. pdtcloud |
---|
99 | |
---|
100 | ! global tendency (clouds+physics) |
---|
101 | REAL sum_subpdq(ngrid,nlay,nq) ! cf. pdqcloud |
---|
102 | REAL sum_subpdt(ngrid,nlay) ! cf. pdtcloud |
---|
103 | |
---|
104 | ! no supersaturation when option supersat is false |
---|
105 | REAL zt(ngrid,nlay) ! local value of temperature |
---|
106 | REAL zqsat(ngrid,nlay) ! saturation |
---|
107 | |
---|
108 | INTEGER iq,ig,l |
---|
109 | LOGICAL,SAVE :: firstcall=.true. |
---|
110 | |
---|
111 | ! Representation of sub-grid water ice clouds A. Pottier 2013 |
---|
112 | REAL :: ztclf(ngrid, nlay) |
---|
113 | REAL :: zqclf(ngrid, nlay,nq) |
---|
114 | REAL :: zdelt |
---|
115 | REAL :: norm |
---|
116 | REAL :: ponder |
---|
117 | REAL :: tcond(ngrid,nlay) |
---|
118 | REAL :: zqvap(ngrid,nlay) |
---|
119 | REAL :: zqice(ngrid,nlay) |
---|
120 | REAL :: spant ! delta T for the temperature distribution |
---|
121 | ! REAL :: zqsat(ngrid,nlay) ! saturation |
---|
122 | REAL :: pteff(ngrid, nlay)! effective temperature in the cloud,neb |
---|
123 | REAL :: pqeff(ngrid, nlay, nq)! effective tracers quantities in the cloud |
---|
124 | REAL :: cloudfrac(ngrid,nlay) ! cloud fraction |
---|
125 | REAL :: mincloud ! min cloud frac |
---|
126 | LOGICAL, save :: flagcloud=.true. |
---|
127 | |
---|
128 | c ** un petit test de coherence |
---|
129 | c -------------------------- |
---|
130 | |
---|
131 | IF (firstcall) THEN |
---|
132 | |
---|
133 | if (nq.gt.nqmx) then |
---|
134 | write(*,*) 'stop in watercloud (nq.gt.nqmx)!' |
---|
135 | write(*,*) 'nq=',nq,' nqmx=',nqmx |
---|
136 | call abort_physic("watercloud","nq.gt.nqmx",1) |
---|
137 | endif |
---|
138 | |
---|
139 | write(*,*) "watercloud: igcm_h2o_vap=",igcm_h2o_vap |
---|
140 | write(*,*) " igcm_h2o_ice=",igcm_h2o_ice |
---|
141 | |
---|
142 | write(*,*) "time subsampling for microphysic ?" |
---|
143 | #ifdef MESOSCALE |
---|
144 | imicro = 2 |
---|
145 | #else |
---|
146 | imicro = 30 |
---|
147 | #endif |
---|
148 | call getin_p("imicro",imicro) |
---|
149 | write(*,*)"watercloud: imicro = ",imicro |
---|
150 | |
---|
151 | firstcall=.false. |
---|
152 | ENDIF ! of IF (firstcall) |
---|
153 | |
---|
154 | !! AS: moved out of firstcall to allow nesting+evoluting timestep |
---|
155 | !! TBD: consider possible diff imicro with domains? |
---|
156 | microtimestep = ptimestep/real(imicro) |
---|
157 | if (microtimestep/=microtimestep_prev) then |
---|
158 | ! only tell the world if microtimestep has changed |
---|
159 | write(*,*)"watercloud: Physical timestep is ",ptimestep |
---|
160 | write(*,*)"watercloud: Microphysics timestep is ",microtimestep |
---|
161 | microtimestep_prev=microtimestep |
---|
162 | endif |
---|
163 | |
---|
164 | c-----Initialization |
---|
165 | sum_subpdq(1:ngrid,1:nlay,1:nq) = 0 |
---|
166 | sum_subpdt(1:ngrid,1:nlay) = 0 |
---|
167 | |
---|
168 | ! default value if no ice |
---|
169 | rhocloud(1:ngrid,1:nlay) = rho_dust |
---|
170 | |
---|
171 | c------------------------------------------------------------------- |
---|
172 | c 0. Representation of sub-grid water ice clouds |
---|
173 | c------------------ |
---|
174 | c-----Initialization |
---|
175 | pteff(1:ngrid,1:nlay) = 0 |
---|
176 | pqeff(1:ngrid,1:nlay,1:nq) = 0 |
---|
177 | DO l=1,nlay |
---|
178 | DO ig=1,ngrid |
---|
179 | pteff(ig,l)=pt(ig,l) |
---|
180 | END DO |
---|
181 | END DO |
---|
182 | DO l=1,nlay |
---|
183 | DO ig=1,ngrid |
---|
184 | DO iq=1,nq |
---|
185 | pqeff(ig,l,iq)=pq(ig,l,iq) |
---|
186 | ENDDO |
---|
187 | ENDDO |
---|
188 | ENDDO |
---|
189 | c-----Tendencies |
---|
190 | DO l=1,nlay |
---|
191 | DO ig=1,ngrid |
---|
192 | ztclf(ig,l)=pt(ig,l)+ pdt(ig,l)*ptimestep |
---|
193 | ENDDO |
---|
194 | ENDDO |
---|
195 | DO l=1,nlay |
---|
196 | DO ig=1,ngrid |
---|
197 | DO iq=1,nq |
---|
198 | zqclf(ig,l,iq)=pq(ig,l,iq)+pdq(ig,l,iq)*ptimestep |
---|
199 | ENDDO |
---|
200 | ENDDO |
---|
201 | ENDDO |
---|
202 | c-----Effective temperature calculation |
---|
203 | IF (CLFvarying) THEN |
---|
204 | spant=3.0 ! delta T for the temprature distribution |
---|
205 | mincloud=0.1 ! min cloudfrac when there is ice |
---|
206 | IF (flagcloud) THEN |
---|
207 | write(*,*) "Delta T", spant |
---|
208 | write(*,*) "mincloud", mincloud |
---|
209 | flagcloud=.false. |
---|
210 | END IF |
---|
211 | !CALL watersat(ngrid*nlay,ztclf,pplay,zqsat) !MV17: we dont need zqsat in the CLFvarying scheme |
---|
212 | zqvap=zqclf(:,:,igcm_h2o_vap) |
---|
213 | zqice=zqclf(:,:,igcm_h2o_ice) |
---|
214 | CALL tcondwater(ngrid*nlay,pplay,zqvap+zqice,tcond) |
---|
215 | DO l=1,nlay |
---|
216 | DO ig=1,ngrid |
---|
217 | zdelt=spant !MAX(spant*ztclf(ig,l),1.e-12), now totally in K. Fixed width |
---|
218 | IF (tcond(ig,l) .ge. (ztclf(ig,l)+zdelt)) THEN |
---|
219 | pteff(ig,l)=ztclf(ig,l) |
---|
220 | cloudfrac(ig,l)=1. |
---|
221 | ELSE IF (tcond(ig,l) .le. (ztclf(ig,l)-zdelt)) THEN |
---|
222 | pteff(ig,l)=ztclf(ig,l)-zdelt |
---|
223 | cloudfrac(ig,l)=mincloud |
---|
224 | ELSE |
---|
225 | cloudfrac(ig,l)=(tcond(ig,l)-ztclf(ig,l)+zdelt)/ |
---|
226 | & (2.0*zdelt) |
---|
227 | pteff(ig,l)=(tcond(ig,l)+ztclf(ig,l)-zdelt)/2. |
---|
228 | END IF |
---|
229 | pteff(ig,l)=pteff(ig,l)-pdt(ig,l)*ptimestep |
---|
230 | IF (cloudfrac(ig,l).le.mincloud) THEN !MV17: replaced .le.0 by .le.mincloud |
---|
231 | cloudfrac(ig,l)=mincloud |
---|
232 | ELSE IF (cloudfrac(ig,l).gt.1) THEN |
---|
233 | cloudfrac(ig,l)=1. |
---|
234 | END IF |
---|
235 | ENDDO |
---|
236 | ENDDO |
---|
237 | c-----Calculation of the total cloud coverage of the column |
---|
238 | DO ig=1,ngrid |
---|
239 | totcloudfrac(ig) = 0. |
---|
240 | norm=0. |
---|
241 | DO l=1,nlay |
---|
242 | ponder=zqice(ig,l)*(pplev(ig,l) - pplev(ig,l+1)) |
---|
243 | totcloudfrac(ig) = totcloudfrac(ig) |
---|
244 | & + cloudfrac(ig,l)*ponder |
---|
245 | norm=norm+ponder |
---|
246 | ENDDO |
---|
247 | totcloudfrac(ig)=MAX(totcloudfrac(ig)/norm,1.e-12) ! min value if NaNs |
---|
248 | ENDDO |
---|
249 | c-----Effective tracers quantities in the cloud fraction |
---|
250 | IF (microphys) THEN |
---|
251 | pqeff(:,:,igcm_ccn_mass)=pq(:,:,igcm_ccn_mass)/ |
---|
252 | & cloudfrac(:,:) |
---|
253 | pqeff(:,:,igcm_ccn_number)=pq(:,:,igcm_ccn_number)/ |
---|
254 | & cloudfrac(:,:) |
---|
255 | END IF ! end if (microphys) |
---|
256 | pqeff(:,:,igcm_h2o_ice)=pq(:,:,igcm_h2o_ice)/ |
---|
257 | & cloudfrac(:,:) |
---|
258 | !! CLFvarying outputs |
---|
259 | CALL WRITEDIAGFI(ngrid,'pqeffice','pqeffice', |
---|
260 | & 'kg/kg',3,pqeff(:,:,igcm_h2o_ice)) |
---|
261 | CALL WRITEDIAGFI(ngrid,'pteff','pteff', |
---|
262 | & 'K',3,pteff(:,:)) |
---|
263 | CALL WRITEDIAGFI(ngrid,'tcond','tcond', |
---|
264 | & 'K',3,tcond(:,:)) |
---|
265 | CALL WRITEDIAGFI(ngrid,'cloudfrac','cloudfrac', |
---|
266 | & 'K',3,cloudfrac(:,:)) |
---|
267 | END IF ! end if (CLFvarying) |
---|
268 | c------------------------------------------------------------------ |
---|
269 | c Time subsampling for microphysics |
---|
270 | c------------------------------------------------------------------ |
---|
271 | rhocloud(1:ngrid,1:nlay) = rho_dust |
---|
272 | DO microstep=1,imicro |
---|
273 | |
---|
274 | c------------------------------------------------------------------- |
---|
275 | c 1. Tendencies: |
---|
276 | c------------------ |
---|
277 | |
---|
278 | |
---|
279 | c------ Temperature tendency subpdt |
---|
280 | ! Each microtimestep we give the cloud scheme a stepped entry subpdt instead of pdt |
---|
281 | ! If imicro=1 subpdt is the same as pdt |
---|
282 | DO l=1,nlay |
---|
283 | DO ig=1,ngrid |
---|
284 | sum_subpdt(ig,l) = sum_subpdt(ig,l) |
---|
285 | & + pdt(ig,l) ! At each micro timestep we add pdt in order to have a stepped entry |
---|
286 | ENDDO |
---|
287 | ENDDO |
---|
288 | c------ Tracers tendencies subpdq are additionned |
---|
289 | c------ At each micro timestep we add pdq in order to have a stepped entry |
---|
290 | IF (microphys) THEN |
---|
291 | DO l=1,nlay |
---|
292 | DO ig=1,ngrid |
---|
293 | sum_subpdq(ig,l,igcm_dust_mass) = |
---|
294 | & sum_subpdq(ig,l,igcm_dust_mass) |
---|
295 | & + pdq(ig,l,igcm_dust_mass) |
---|
296 | sum_subpdq(ig,l,igcm_dust_number) = |
---|
297 | & sum_subpdq(ig,l,igcm_dust_number) |
---|
298 | & + pdq(ig,l,igcm_dust_number) |
---|
299 | sum_subpdq(ig,l,igcm_ccn_mass) = |
---|
300 | & sum_subpdq(ig,l,igcm_ccn_mass) |
---|
301 | & + pdq(ig,l,igcm_ccn_mass) |
---|
302 | sum_subpdq(ig,l,igcm_ccn_number) = |
---|
303 | & sum_subpdq(ig,l,igcm_ccn_number) |
---|
304 | & + pdq(ig,l,igcm_ccn_number) |
---|
305 | ENDDO |
---|
306 | ENDDO |
---|
307 | ENDIF |
---|
308 | DO l=1,nlay |
---|
309 | DO ig=1,ngrid |
---|
310 | sum_subpdq(ig,l,igcm_h2o_ice) = |
---|
311 | & sum_subpdq(ig,l,igcm_h2o_ice) |
---|
312 | & + pdq(ig,l,igcm_h2o_ice) |
---|
313 | sum_subpdq(ig,l,igcm_h2o_vap) = |
---|
314 | & sum_subpdq(ig,l,igcm_h2o_vap) |
---|
315 | & + pdq(ig,l,igcm_h2o_vap) |
---|
316 | IF (hdo) THEN |
---|
317 | sum_subpdq(ig,l,igcm_hdo_ice) = |
---|
318 | & sum_subpdq(ig,l,igcm_hdo_ice) |
---|
319 | & + pdq(ig,l,igcm_hdo_ice) |
---|
320 | sum_subpdq(ig,l,igcm_hdo_vap) = |
---|
321 | & sum_subpdq(ig,l,igcm_hdo_vap) |
---|
322 | & + pdq(ig,l,igcm_hdo_vap) |
---|
323 | ENDIF !hdo |
---|
324 | ENDDO |
---|
325 | ENDDO |
---|
326 | |
---|
327 | c------------------------------------------------------------------- |
---|
328 | c 2. Main call to the different cloud schemes: |
---|
329 | c------------------------------------------------ |
---|
330 | IF (microphys) THEN |
---|
331 | CALL improvedclouds(ngrid,nlay,microtimestep, |
---|
332 | & pplay,pteff,sum_subpdt, |
---|
333 | & pqeff,sum_subpdq,subpdqcloud,subpdtcloud, |
---|
334 | & nq,tauscaling) |
---|
335 | |
---|
336 | ELSE |
---|
337 | CALL simpleclouds(ngrid,nlay,microtimestep, |
---|
338 | & pplay,pzlay,pteff,sum_subpdt, |
---|
339 | & pqeff,sum_subpdq,subpdqcloud,subpdtcloud, |
---|
340 | & nq,tau,rice) |
---|
341 | ENDIF |
---|
342 | |
---|
343 | c------------------------------------------------------------------- |
---|
344 | c 3. Updating tendencies after cloud scheme: |
---|
345 | c----------------------------------------------- |
---|
346 | |
---|
347 | IF (microphys) THEN |
---|
348 | DO l=1,nlay |
---|
349 | DO ig=1,ngrid |
---|
350 | sum_subpdq(ig,l,igcm_dust_mass) = |
---|
351 | & sum_subpdq(ig,l,igcm_dust_mass) |
---|
352 | & + subpdqcloud(ig,l,igcm_dust_mass) |
---|
353 | sum_subpdq(ig,l,igcm_dust_number) = |
---|
354 | & sum_subpdq(ig,l,igcm_dust_number) |
---|
355 | & + subpdqcloud(ig,l,igcm_dust_number) |
---|
356 | sum_subpdq(ig,l,igcm_ccn_mass) = |
---|
357 | & sum_subpdq(ig,l,igcm_ccn_mass) |
---|
358 | & + subpdqcloud(ig,l,igcm_ccn_mass) |
---|
359 | sum_subpdq(ig,l,igcm_ccn_number) = |
---|
360 | & sum_subpdq(ig,l,igcm_ccn_number) |
---|
361 | & + subpdqcloud(ig,l,igcm_ccn_number) |
---|
362 | ENDDO |
---|
363 | ENDDO |
---|
364 | ENDIF |
---|
365 | DO l=1,nlay |
---|
366 | DO ig=1,ngrid |
---|
367 | sum_subpdq(ig,l,igcm_h2o_ice) = |
---|
368 | & sum_subpdq(ig,l,igcm_h2o_ice) |
---|
369 | & + subpdqcloud(ig,l,igcm_h2o_ice) |
---|
370 | sum_subpdq(ig,l,igcm_h2o_vap) = |
---|
371 | & sum_subpdq(ig,l,igcm_h2o_vap) |
---|
372 | & + subpdqcloud(ig,l,igcm_h2o_vap) |
---|
373 | |
---|
374 | IF (hdo) THEN |
---|
375 | sum_subpdq(ig,l,igcm_hdo_ice) = |
---|
376 | & sum_subpdq(ig,l,igcm_hdo_ice) |
---|
377 | & + subpdqcloud(ig,l,igcm_hdo_ice) |
---|
378 | sum_subpdq(ig,l,igcm_hdo_vap) = |
---|
379 | & sum_subpdq(ig,l,igcm_hdo_vap) |
---|
380 | & + subpdqcloud(ig,l,igcm_hdo_vap) |
---|
381 | ENDIF ! hdo |
---|
382 | |
---|
383 | ENDDO |
---|
384 | ENDDO |
---|
385 | |
---|
386 | IF (activice) THEN |
---|
387 | DO l=1,nlay |
---|
388 | DO ig=1,ngrid |
---|
389 | sum_subpdt(ig,l) = |
---|
390 | & sum_subpdt(ig,l) + subpdtcloud(ig,l) |
---|
391 | ENDDO |
---|
392 | ENDDO |
---|
393 | ENDIF |
---|
394 | |
---|
395 | ! !! Example of how to use writediagmicrofi useful to |
---|
396 | ! !! get outputs at each microphysical sub-timestep (better to be used in 1D) |
---|
397 | ! CALL WRITEDIAGMICROFI(ngrid,imicro,microstep, |
---|
398 | ! & microtimestep,'subpdtcloud', |
---|
399 | ! & 'subpdtcloud','K/s',1,subpdtcloud(:,:)) |
---|
400 | |
---|
401 | ENDDO ! of DO microstep=1,imicro |
---|
402 | |
---|
403 | c------------------------------------------------------------------- |
---|
404 | c 6. Compute final tendencies after time loop: |
---|
405 | c------------------------------------------------ |
---|
406 | |
---|
407 | c------ Temperature tendency pdtcloud |
---|
408 | DO l=1,nlay |
---|
409 | DO ig=1,ngrid |
---|
410 | pdtcloud(ig,l) = |
---|
411 | & sum_subpdt(ig,l)/real(imicro)-pdt(ig,l) |
---|
412 | ENDDO |
---|
413 | ENDDO |
---|
414 | |
---|
415 | c------ Tracers tendencies pdqcloud |
---|
416 | DO l=1,nlay |
---|
417 | DO ig=1,ngrid |
---|
418 | pdqcloud(ig,l,igcm_h2o_ice) = |
---|
419 | & sum_subpdq(ig,l,igcm_h2o_ice)/real(imicro) |
---|
420 | & - pdq(ig,l,igcm_h2o_ice) |
---|
421 | pdqcloud(ig,l,igcm_h2o_vap) = |
---|
422 | & sum_subpdq(ig,l,igcm_h2o_vap)/real(imicro) |
---|
423 | & - pdq(ig,l,igcm_h2o_vap) |
---|
424 | IF (hdo) THEN |
---|
425 | pdqcloud(ig,l,igcm_hdo_ice) = |
---|
426 | & sum_subpdq(ig,l,igcm_hdo_ice)/real(imicro) |
---|
427 | & - pdq(ig,l,igcm_hdo_ice) |
---|
428 | pdqcloud(ig,l,igcm_hdo_vap) = |
---|
429 | & sum_subpdq(ig,l,igcm_hdo_vap)/real(imicro) |
---|
430 | & - pdq(ig,l,igcm_hdo_vap) |
---|
431 | ENDIF !hdo |
---|
432 | ENDDO |
---|
433 | ENDDO |
---|
434 | |
---|
435 | IF(microphys) THEN |
---|
436 | DO l=1,nlay |
---|
437 | DO ig=1,ngrid |
---|
438 | pdqcloud(ig,l,igcm_ccn_mass) = |
---|
439 | & sum_subpdq(ig,l,igcm_ccn_mass)/real(imicro) |
---|
440 | & - pdq(ig,l,igcm_ccn_mass) |
---|
441 | pdqcloud(ig,l,igcm_ccn_number) = |
---|
442 | & sum_subpdq(ig,l,igcm_ccn_number)/real(imicro) |
---|
443 | & - pdq(ig,l,igcm_ccn_number) |
---|
444 | ENDDO |
---|
445 | ENDDO |
---|
446 | ENDIF |
---|
447 | |
---|
448 | IF(scavenging) THEN |
---|
449 | DO l=1,nlay |
---|
450 | DO ig=1,ngrid |
---|
451 | pdqcloud(ig,l,igcm_dust_mass) = |
---|
452 | & sum_subpdq(ig,l,igcm_dust_mass)/real(imicro) |
---|
453 | & - pdq(ig,l,igcm_dust_mass) |
---|
454 | pdqcloud(ig,l,igcm_dust_number) = |
---|
455 | & sum_subpdq(ig,l,igcm_dust_number)/real(imicro) |
---|
456 | & - pdq(ig,l,igcm_dust_number) |
---|
457 | ENDDO |
---|
458 | ENDDO |
---|
459 | ENDIF |
---|
460 | |
---|
461 | c------- Due to stepped entry, other processes tendencies can add up to negative values |
---|
462 | c------- Therefore, enforce positive values and conserve mass |
---|
463 | IF(microphys) THEN |
---|
464 | DO l=1,nlay |
---|
465 | DO ig=1,ngrid |
---|
466 | IF ((pq(ig,l,igcm_ccn_number) + |
---|
467 | & ptimestep* (pdq(ig,l,igcm_ccn_number) + |
---|
468 | & pdqcloud(ig,l,igcm_ccn_number)) .le. 1.) |
---|
469 | & .or. (pq(ig,l,igcm_ccn_mass) + |
---|
470 | & ptimestep* (pdq(ig,l,igcm_ccn_mass) + |
---|
471 | & pdqcloud(ig,l,igcm_ccn_mass)) .le. 1.e-20)) THEN |
---|
472 | pdqcloud(ig,l,igcm_ccn_number) = |
---|
473 | & - pq(ig,l,igcm_ccn_number)/ptimestep |
---|
474 | & - pdq(ig,l,igcm_ccn_number) + 1. |
---|
475 | pdqcloud(ig,l,igcm_dust_number) = |
---|
476 | & -pdqcloud(ig,l,igcm_ccn_number) |
---|
477 | pdqcloud(ig,l,igcm_ccn_mass) = |
---|
478 | & - pq(ig,l,igcm_ccn_mass)/ptimestep |
---|
479 | & - pdq(ig,l,igcm_ccn_mass) + 1.e-20 |
---|
480 | pdqcloud(ig,l,igcm_dust_mass) = |
---|
481 | & -pdqcloud(ig,l,igcm_ccn_mass) |
---|
482 | ENDIF |
---|
483 | ENDDO |
---|
484 | ENDDO |
---|
485 | ENDIF |
---|
486 | |
---|
487 | IF(scavenging) THEN |
---|
488 | DO l=1,nlay |
---|
489 | DO ig=1,ngrid |
---|
490 | IF ((pq(ig,l,igcm_dust_number) + |
---|
491 | & ptimestep* (pdq(ig,l,igcm_dust_number) + |
---|
492 | & pdqcloud(ig,l,igcm_dust_number)) .le. 1.) |
---|
493 | & .or. (pq(ig,l,igcm_dust_mass) + |
---|
494 | & ptimestep* (pdq(ig,l,igcm_dust_mass) + |
---|
495 | & pdqcloud(ig,l,igcm_dust_mass)) .le. 1.e-20)) THEN |
---|
496 | pdqcloud(ig,l,igcm_dust_number) = |
---|
497 | & - pq(ig,l,igcm_dust_number)/ptimestep |
---|
498 | & - pdq(ig,l,igcm_dust_number) + 1. |
---|
499 | pdqcloud(ig,l,igcm_ccn_number) = |
---|
500 | & -pdqcloud(ig,l,igcm_dust_number) |
---|
501 | pdqcloud(ig,l,igcm_dust_mass) = |
---|
502 | & - pq(ig,l,igcm_dust_mass)/ptimestep |
---|
503 | & - pdq(ig,l,igcm_dust_mass) + 1.e-20 |
---|
504 | pdqcloud(ig,l,igcm_ccn_mass) = |
---|
505 | & -pdqcloud(ig,l,igcm_dust_mass) |
---|
506 | ENDIF |
---|
507 | ENDDO |
---|
508 | ENDDO |
---|
509 | ENDIF |
---|
510 | |
---|
511 | DO l=1,nlay |
---|
512 | DO ig=1,ngrid |
---|
513 | |
---|
514 | IF (pq(ig,l,igcm_h2o_ice) + ptimestep* |
---|
515 | & (pdq(ig,l,igcm_h2o_ice) + pdqcloud(ig,l,igcm_h2o_ice)) |
---|
516 | & .le. 1.e-8) THEN |
---|
517 | pdqcloud(ig,l,igcm_h2o_ice) = |
---|
518 | & - pq(ig,l,igcm_h2o_ice)/ptimestep - pdq(ig,l,igcm_h2o_ice) |
---|
519 | pdqcloud(ig,l,igcm_h2o_vap) = -pdqcloud(ig,l,igcm_h2o_ice) |
---|
520 | if (hdo) then |
---|
521 | pdqcloud(ig,l,igcm_hdo_ice) = |
---|
522 | & - pq(ig,l,igcm_hdo_ice)/ptimestep - pdq(ig,l,igcm_hdo_ice) |
---|
523 | pdqcloud(ig,l,igcm_hdo_vap) = -pdqcloud(ig,l,igcm_hdo_ice) |
---|
524 | endif |
---|
525 | ENDIF |
---|
526 | IF (pq(ig,l,igcm_h2o_vap) + ptimestep* |
---|
527 | & (pdq(ig,l,igcm_h2o_vap) + pdqcloud(ig,l,igcm_h2o_vap)) |
---|
528 | & .le. 1.e-8) THEN |
---|
529 | pdqcloud(ig,l,igcm_h2o_vap) = |
---|
530 | & - pq(ig,l,igcm_h2o_vap)/ptimestep - pdq(ig,l,igcm_h2o_vap) |
---|
531 | pdqcloud(ig,l,igcm_h2o_ice) = -pdqcloud(ig,l,igcm_h2o_vap) |
---|
532 | if (hdo) then |
---|
533 | pdqcloud(ig,l,igcm_hdo_vap) = |
---|
534 | & - pq(ig,l,igcm_hdo_vap)/ptimestep - pdq(ig,l,igcm_hdo_vap) |
---|
535 | pdqcloud(ig,l,igcm_hdo_ice) = -pdqcloud(ig,l,igcm_hdo_vap) |
---|
536 | endif |
---|
537 | ENDIF |
---|
538 | |
---|
539 | ENDDO |
---|
540 | ENDDO |
---|
541 | |
---|
542 | c------Update the ice and dust particle size "rice" for output or photochemistry |
---|
543 | c------Only rsedcloud is used for the water cycle |
---|
544 | |
---|
545 | IF(scavenging) THEN |
---|
546 | DO l=1, nlay |
---|
547 | DO ig=1,ngrid |
---|
548 | |
---|
549 | call updaterdust( |
---|
550 | & pq(ig,l,igcm_dust_mass) + ! dust mass |
---|
551 | & (pdq(ig,l,igcm_dust_mass) + ! dust mass |
---|
552 | & pdqcloud(ig,l,igcm_dust_mass))*ptimestep, ! dust mass |
---|
553 | & pq(ig,l,igcm_dust_number) + ! dust number |
---|
554 | & (pdq(ig,l,igcm_dust_number) + ! dust number |
---|
555 | & pdqcloud(ig,l,igcm_dust_number))*ptimestep, ! dust number |
---|
556 | & rdust(ig,l)) |
---|
557 | |
---|
558 | ENDDO |
---|
559 | ENDDO |
---|
560 | ENDIF |
---|
561 | |
---|
562 | IF(microphys) THEN |
---|
563 | |
---|
564 | ! In case one does not want to allow supersatured water when using microphysics. |
---|
565 | ! Not done by default. |
---|
566 | IF(.not.supersat) THEN |
---|
567 | zt = pt + (pdt+pdtcloud)*ptimestep |
---|
568 | call watersat(ngrid*nlay,zt,pplay,zqsat) |
---|
569 | DO l=1, nlay |
---|
570 | DO ig=1,ngrid |
---|
571 | IF (pq(ig,l,igcm_h2o_vap) |
---|
572 | & + (pdq(ig,l,igcm_h2o_vap) + pdqcloud(ig,l,igcm_h2o_vap)) |
---|
573 | & * ptimestep .ge. zqsat(ig,l)) THEN |
---|
574 | pdqcloud(ig,l,igcm_h2o_vap) = |
---|
575 | & (zqsat(ig,l) - pq(ig,l,igcm_h2o_vap))/ptimestep |
---|
576 | & - pdq(ig,l,igcm_h2o_vap) |
---|
577 | pdqcloud(ig,l,igcm_h2o_ice) = |
---|
578 | & -pdqcloud(ig,l,igcm_h2o_vap) |
---|
579 | ! no need to correct ccn_number, updaterad can handle this properly. |
---|
580 | ENDIF |
---|
581 | ENDDO |
---|
582 | ENDDO |
---|
583 | ENDIF |
---|
584 | |
---|
585 | DO l=1, nlay |
---|
586 | DO ig=1,ngrid |
---|
587 | |
---|
588 | call updaterice_micro( |
---|
589 | & pq(ig,l,igcm_h2o_ice) + ! ice mass |
---|
590 | & (pdq(ig,l,igcm_h2o_ice) + ! ice mass |
---|
591 | & pdqcloud(ig,l,igcm_h2o_ice))*ptimestep, ! ice mass |
---|
592 | & pq(ig,l,igcm_ccn_mass) + ! ccn mass |
---|
593 | & (pdq(ig,l,igcm_ccn_mass) + ! ccn mass |
---|
594 | & pdqcloud(ig,l,igcm_ccn_mass))*ptimestep, ! ccn mass |
---|
595 | & pq(ig,l,igcm_ccn_number) + ! ccn number |
---|
596 | & (pdq(ig,l,igcm_ccn_number) + ! ccn number |
---|
597 | & pdqcloud(ig,l,igcm_ccn_number))*ptimestep, ! ccn number |
---|
598 | & tauscaling(ig),rice(ig,l),rhocloud(ig,l)) |
---|
599 | |
---|
600 | ENDDO |
---|
601 | ENDDO |
---|
602 | |
---|
603 | ELSE ! no microphys |
---|
604 | |
---|
605 | DO l=1,nlay |
---|
606 | DO ig=1,ngrid |
---|
607 | |
---|
608 | call updaterice_typ( |
---|
609 | & pq(ig,l,igcm_h2o_ice) + ! ice mass |
---|
610 | & (pdq(ig,l,igcm_h2o_ice) + ! ice mass |
---|
611 | & pdqcloud(ig,l,igcm_h2o_ice))*ptimestep, ! ice mass |
---|
612 | & tau(ig,1),pzlay(ig,l),rice(ig,l)) |
---|
613 | |
---|
614 | ENDDO |
---|
615 | ENDDO |
---|
616 | |
---|
617 | ENDIF ! of IF(microphys) |
---|
618 | |
---|
619 | |
---|
620 | |
---|
621 | c A correction if a lot of subliming CO2 fills the 1st layer FF04/2005 |
---|
622 | c ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ |
---|
623 | c Then that should not affect the ice particle radius |
---|
624 | do ig=1,ngrid |
---|
625 | if(pdpsrf(ig)*ptimestep.gt.0.9*(pplev(ig,1)-pplev(ig,2)))then |
---|
626 | if(pdpsrf(ig)*ptimestep.gt.0.9*(pplev(ig,1)-pplev(ig,3))) |
---|
627 | & rice(ig,2)=rice(ig,3) |
---|
628 | rice(ig,1)=rice(ig,2) |
---|
629 | end if |
---|
630 | end do |
---|
631 | |
---|
632 | |
---|
633 | DO l=1,nlay |
---|
634 | DO ig=1,ngrid |
---|
635 | rsedcloud(ig,l)=max(rice(ig,l)* |
---|
636 | & (1.+nuice_sed)*(1.+nuice_sed)*(1.+nuice_sed), |
---|
637 | & rdust(ig,l)) |
---|
638 | ! rsedcloud(ig,l)=min(rsedcloud(ig,l),1.e-4) |
---|
639 | ENDDO |
---|
640 | ENDDO |
---|
641 | |
---|
642 | ! used for rad. transfer calculations |
---|
643 | ! nuice is constant because a lognormal distribution is prescribed |
---|
644 | nuice(1:ngrid,1:nlay)=nuice_ref |
---|
645 | |
---|
646 | c------Update tendencies for sub-grid water ice clouds |
---|
647 | IF (CLFvarying) THEN |
---|
648 | DO ig=1,ngrid |
---|
649 | DO l=1,nlay |
---|
650 | pdqcloud(ig,l,igcm_dust_mass)=pdqcloud(ig,l,igcm_dust_mass) |
---|
651 | & *cloudfrac(ig,l) |
---|
652 | pdqcloud(ig,l,igcm_ccn_mass)=pdqcloud(ig,l,igcm_ccn_mass) |
---|
653 | & *cloudfrac(ig,l) |
---|
654 | pdqcloud(ig,l,igcm_dust_number)=pdqcloud(ig,l, |
---|
655 | & igcm_dust_number) *cloudfrac(ig,l) |
---|
656 | pdqcloud(ig,l,igcm_ccn_number)=pdqcloud(ig,l, |
---|
657 | & igcm_ccn_number) *cloudfrac(ig,l) |
---|
658 | pdqcloud(ig,l,igcm_h2o_vap)=pdqcloud(ig,l, |
---|
659 | & igcm_h2o_vap) *cloudfrac(ig,l) |
---|
660 | pdqcloud(ig,l,igcm_h2o_ice)=pdqcloud(ig,l, |
---|
661 | & igcm_h2o_ice) *cloudfrac(ig,l) |
---|
662 | ENDDO |
---|
663 | ENDDO |
---|
664 | pdtcloud(:,:)=pdtcloud(:,:)*cloudfrac(:,:) |
---|
665 | ENDIF |
---|
666 | #ifndef MESOSCALE |
---|
667 | c======================================================================= |
---|
668 | call WRITEDIAGFI(ngrid,"pdqice2","pdqcloudice apres microphysique" |
---|
669 | & ,"kg/kg.s-1",3,pdqcloud(1:ngrid,1:nlay,igcm_h2o_ice)) |
---|
670 | call WRITEDIAGFI(ngrid,"pdqvap2","pdqcloudvap apres microphysique" |
---|
671 | & ,"kg/kg.s-1",3,pdqcloud(1:ngrid,1:nlay, |
---|
672 | & igcm_h2o_vap)) |
---|
673 | if (hdo) then |
---|
674 | call WRITEDIAGFI(ngrid,"pdqiceD","pdqiceD apres microphysique" |
---|
675 | & ,"kg/kg.s-1",3,pdqcloud(1:ngrid,1:nlay,igcm_hdo_ice)) |
---|
676 | call WRITEDIAGFI(ngrid,"pdqvapD","pdqvapD apres microphysique" |
---|
677 | & ,"kg/kg.s-1",3,pdqcloud(1:ngrid,1:nlay, |
---|
678 | & igcm_hdo_vap)) |
---|
679 | endif |
---|
680 | call WRITEDIAGFI(ngrid,"pdqccn2","pdqcloudccn apres microphysique" |
---|
681 | & ,"kg/kg.s-1",3,pdqcloud(1:ngrid,1:nlay, |
---|
682 | & igcm_ccn_mass)) |
---|
683 | call WRITEDIAGFI(ngrid,"pdqccnN2","pdqcloudccnN apres "// |
---|
684 | & "microphysique","nb/kg.s-1",3,pdqcloud(1:ngrid,1:nlay, |
---|
685 | & igcm_ccn_number)) |
---|
686 | call WRITEDIAGFI(ngrid,"pdqdust2", "pdqclouddust apres "// |
---|
687 | & "microphysique","kg/kg.s-1",3,pdqcloud(1:ngrid,1:nlay, |
---|
688 | & igcm_dust_mass)) |
---|
689 | call WRITEDIAGFI(ngrid,"pdqdustN2", "pdqclouddustN apres "// |
---|
690 | & "microphysique","nb/kg.s-1",3,pdqcloud(1:ngrid,1:nlay, |
---|
691 | & igcm_dust_number)) |
---|
692 | c======================================================================= |
---|
693 | #endif |
---|
694 | |
---|
695 | END SUBROUTINE watercloud |
---|
696 | |
---|
697 | subroutine ini_watercloud_mod(ngrid,nlayer,nq) |
---|
698 | implicit none |
---|
699 | |
---|
700 | integer,intent(in) :: ngrid ! number of atmospheric columns |
---|
701 | integer,intent(in) :: nlayer ! number of atmospheric layers |
---|
702 | integer,intent(in) :: nq ! number of tracers |
---|
703 | |
---|
704 | allocate(zdqcloud(ngrid,nlayer,nq)) |
---|
705 | zdqcloud(:,:,:)=0 |
---|
706 | allocate(zdqscloud(ngrid,nq)) |
---|
707 | zdqscloud(:,:)=0 |
---|
708 | |
---|
709 | end subroutine ini_watercloud_mod |
---|
710 | |
---|
711 | |
---|
712 | subroutine end_watercloud_mod |
---|
713 | implicit none |
---|
714 | |
---|
715 | if (allocated(zdqcloud)) deallocate(zdqcloud) |
---|
716 | if (allocated(zdqscloud)) deallocate(zdqscloud) |
---|
717 | |
---|
718 | end subroutine end_watercloud_mod |
---|
719 | |
---|
720 | END MODULE watercloud_mod |
---|