1 | subroutine condense_cloud(ngrid,nlayer,nq,ptimestep, & |
---|
2 | pcapcal,pplay,pplev,ptsrf,pt, & |
---|
3 | pphi,pdt,pdu,pdv,pdtsrf,pu,pv,pq,pdq, & |
---|
4 | piceco2,psolaralb,pemisurf, & |
---|
5 | pdtc,pdtsrfc,pdpsrf,pduc,pdvc, & |
---|
6 | pdqc,reffrad,cpp3D) |
---|
7 | |
---|
8 | use radinc_h, only : naerkind |
---|
9 | use gases_h |
---|
10 | |
---|
11 | implicit none |
---|
12 | |
---|
13 | !================================================================== |
---|
14 | ! Purpose |
---|
15 | ! ------- |
---|
16 | ! Condense and/or sublime CO2 ice on the ground and in the |
---|
17 | ! atmosphere, and sediment the ice. |
---|
18 | ! |
---|
19 | ! Inputs |
---|
20 | ! ------ |
---|
21 | ! ngrid Number of vertical columns |
---|
22 | ! nlayer Number of layers |
---|
23 | ! pplay(ngrid,nlayer) Pressure layers |
---|
24 | ! pplev(ngrid,nlayer+1) Pressure levels |
---|
25 | ! pt(ngrid,nlayer) Temperature (in K) |
---|
26 | ! ptsrf(ngrid) Surface temperature |
---|
27 | ! |
---|
28 | ! pdt(ngrid,nlayermx) Time derivative before condensation/sublimation of pt |
---|
29 | ! pdtsrf(ngrid) Time derivative before condensation/sublimation of ptsrf |
---|
30 | ! pqsurf(ngrid,nq) Sedimentation flux at the surface (kg.m-2.s-1) |
---|
31 | ! |
---|
32 | ! Outputs |
---|
33 | ! ------- |
---|
34 | ! pdpsrf(ngrid) \ Contribution of condensation/sublimation |
---|
35 | ! pdtc(ngrid,nlayermx) / to the time derivatives of Ps, pt, and ptsrf |
---|
36 | ! pdtsrfc(ngrid) / |
---|
37 | ! |
---|
38 | ! Both |
---|
39 | ! ---- |
---|
40 | ! piceco2(ngrid) CO2 ice at the surface (kg/m2) |
---|
41 | ! psolaralb(ngrid) Albedo at the surface |
---|
42 | ! pemisurf(ngrid) Emissivity of the surface |
---|
43 | ! |
---|
44 | ! Authors |
---|
45 | ! ------- |
---|
46 | ! Francois Forget (1996) |
---|
47 | ! Converted to Fortran 90 and slightly modified by R. Wordsworth (2009) |
---|
48 | ! |
---|
49 | !================================================================== |
---|
50 | |
---|
51 | #include "dimensions.h" |
---|
52 | #include "dimphys.h" |
---|
53 | #include "comcstfi.h" |
---|
54 | #include "surfdat.h" |
---|
55 | #include "comgeomfi.h" |
---|
56 | #include "comvert.h" |
---|
57 | #include "callkeys.h" |
---|
58 | #include "tracer.h" |
---|
59 | |
---|
60 | !----------------------------------------------------------------------- |
---|
61 | ! Arguments |
---|
62 | |
---|
63 | INTEGER ngrid, nlayer, nq |
---|
64 | |
---|
65 | REAL ptimestep |
---|
66 | REAL pplay(ngrid,nlayer),pplev(ngrid,nlayer+1) |
---|
67 | REAL pcapcal(ngrid) |
---|
68 | REAL pt(ngrid,nlayer) |
---|
69 | REAL ptsrf(ngrid) |
---|
70 | REAL pphi(ngrid,nlayer) |
---|
71 | REAL pdt(ngrid,nlayer),pdtsrf(ngrid),pdtc(ngrid,nlayer) |
---|
72 | REAL pdtsrfc(ngrid),pdpsrf(ngrid) |
---|
73 | ! REAL piceco2(ngrid),psolaralb(ngrid,2),pemisurf(ngrid) |
---|
74 | REAL piceco2(ngrid),psolaralb(ngrid),pemisurf(ngrid) |
---|
75 | |
---|
76 | REAL pu(ngrid,nlayer) , pv(ngrid,nlayer) |
---|
77 | REAL pdu(ngrid,nlayer) , pdv(ngrid,nlayer) |
---|
78 | REAL pduc(ngrid,nlayer) , pdvc(ngrid,nlayer) |
---|
79 | REAL pq(ngridmx,nlayer,nq),pdq(ngrid,nlayer,nq) |
---|
80 | REAL pdqc(ngrid,nlayer,nq), pdqsc(ngrid,nq) |
---|
81 | |
---|
82 | REAL reffrad(ngrid,nlayer,naerkind) |
---|
83 | |
---|
84 | ! Allow variations in cp with location |
---|
85 | REAL cpp3D(ngridmx,nlayermx) ! specific heat capacity at const. pressure |
---|
86 | |
---|
87 | |
---|
88 | |
---|
89 | !----------------------------------------------------------------------- |
---|
90 | ! Local variables |
---|
91 | |
---|
92 | INTEGER l,ig,icap,ilay,it,iq |
---|
93 | |
---|
94 | REAL*8 zt(ngridmx,nlayermx) |
---|
95 | REAL zq(ngridmx,nlayermx,nqmx) |
---|
96 | REAL zcpi |
---|
97 | REAL ztcond (ngridmx,nlayermx) |
---|
98 | REAL ztcondsol(ngridmx) |
---|
99 | REAL zdiceco2(ngridmx) |
---|
100 | REAL zcondicea(ngridmx,nlayermx), zcondices(ngridmx) |
---|
101 | REAL zfallice(ngridmx), Mfallice(ngridmx) |
---|
102 | REAL zmflux(nlayermx+1) |
---|
103 | REAL zu(nlayermx),zv(nlayermx) |
---|
104 | REAL ztsrf(ngridmx) |
---|
105 | REAL ztc(nlayermx), ztm(nlayermx+1) |
---|
106 | REAL zum(nlayermx+1) , zvm(nlayermx+1) |
---|
107 | LOGICAL condsub(ngridmx) |
---|
108 | REAL subptimestep |
---|
109 | Integer Ntime |
---|
110 | real masse (ngridmx,nlayermx), w(ngridmx,nlayermx,nqmx) |
---|
111 | real wq(ngridmx,nlayermx+1) |
---|
112 | real vstokes,reff |
---|
113 | |
---|
114 | ! Special diagnostic variables |
---|
115 | real tconda1(ngridmx,nlayermx) |
---|
116 | real tconda2(ngridmx,nlayermx) |
---|
117 | real zdtsig (ngridmx,nlayermx) |
---|
118 | real zdt (ngridmx,nlayermx) |
---|
119 | |
---|
120 | !----------------------------------------------------------------------- |
---|
121 | ! Saved local variables |
---|
122 | |
---|
123 | REAL emisref(ngridmx) |
---|
124 | REAL latcond |
---|
125 | REAL ccond |
---|
126 | REAL cpice |
---|
127 | SAVE emisref, cpice |
---|
128 | SAVE latcond,ccond |
---|
129 | |
---|
130 | LOGICAL firstcall |
---|
131 | SAVE firstcall |
---|
132 | REAL SSUM |
---|
133 | EXTERNAL SSUM |
---|
134 | |
---|
135 | DATA latcond /5.9e5/ |
---|
136 | DATA cpice /1000./ |
---|
137 | DATA firstcall/.true./ |
---|
138 | |
---|
139 | REAL CBRT |
---|
140 | EXTERNAL CBRT |
---|
141 | |
---|
142 | INTEGER,SAVE :: i_co2ice=0 ! co2 ice |
---|
143 | CHARACTER(LEN=20) :: tracername ! to temporarily store text |
---|
144 | |
---|
145 | integer igas |
---|
146 | integer,save :: igasco2=0 |
---|
147 | character(len=3) :: gasname |
---|
148 | |
---|
149 | real reffradmin, reffradmax |
---|
150 | |
---|
151 | real ppco2 |
---|
152 | |
---|
153 | !----------------------------------------------------------------------- |
---|
154 | ! Initializations |
---|
155 | |
---|
156 | call zerophys(ngrid*nlayer*nq,pdqc) |
---|
157 | call zerophys(ngrid*nlayer*nq,pdtc) |
---|
158 | call zerophys(ngridmx*nlayermx*nqmx,zq) |
---|
159 | call zerophys(ngridmx*nlayermx,zt) |
---|
160 | |
---|
161 | !reffradmin=1.e-7 |
---|
162 | !reffradmax=5.e-4 |
---|
163 | !reffradmin=0.5e-7 |
---|
164 | !reffradmax=0.1e-3 ! FF data |
---|
165 | reffradmin=0.1e-5 |
---|
166 | reffradmax=0.1e-3 ! JB data |
---|
167 | ! improve this in the future... |
---|
168 | |
---|
169 | ! Initialisation |
---|
170 | IF (firstcall) THEN |
---|
171 | |
---|
172 | ! find CO2 ice tracer |
---|
173 | do iq=1,nqmx |
---|
174 | tracername=noms(iq) |
---|
175 | if (tracername.eq."co2_ice") then |
---|
176 | i_co2ice=iq |
---|
177 | endif |
---|
178 | enddo |
---|
179 | |
---|
180 | ! find CO2 gas |
---|
181 | do igas=1,ngasmx |
---|
182 | gasname=gnom(igas) |
---|
183 | ! gasname=noms(igas) ! was a bug |
---|
184 | if (gasname.eq."CO2") then |
---|
185 | igasco2=igas |
---|
186 | endif |
---|
187 | enddo |
---|
188 | |
---|
189 | write(*,*) "condense_cloud: i_co2ice=",i_co2ice |
---|
190 | |
---|
191 | if((i_co2ice.lt.1))then |
---|
192 | print*,'In condens_cloud but no CO2 ice tracer, exiting.' |
---|
193 | print*,'Still need generalisation to arbitrary species!' |
---|
194 | stop |
---|
195 | endif |
---|
196 | |
---|
197 | ccond=cpp/(g*latcond) |
---|
198 | print*,'In condens_cloud: ccond=',ccond,' latcond=',latcond |
---|
199 | |
---|
200 | ! Prepare special treatment if gas is not pure CO2 |
---|
201 | !if (addn2) then |
---|
202 | ! m_co2 = 44.01E-3 ! CO2 molecular mass (kg/mol) |
---|
203 | ! m_noco2 = 28.02E-3 ! N2 molecular mass (kg/mol) |
---|
204 | ! Compute A and B coefficient use to compute |
---|
205 | ! mean molecular mass Mair defined by |
---|
206 | ! 1/Mair = q(ico2)/m_co2 + (1-q(ico2))/m_noco2 |
---|
207 | ! 1/Mair = A*q(ico2) + B |
---|
208 | ! A = (1/m_co2 - 1/m_noco2) |
---|
209 | ! B = 1/m_noco2 |
---|
210 | !endif |
---|
211 | |
---|
212 | ! Minimum CO2 mixing ratio below which mixing occurs with layer above: |
---|
213 | !qco2min =0.75 |
---|
214 | |
---|
215 | firstcall=.false. |
---|
216 | ENDIF |
---|
217 | zcpi=1./cpp |
---|
218 | |
---|
219 | !----------------------------------------------------------------------- |
---|
220 | ! Calculation of CO2 condensation / sublimation |
---|
221 | ! |
---|
222 | ! Variables used: |
---|
223 | ! piceco2(ngrid) amount of co2 ice on the ground (kg/m2) |
---|
224 | ! zcondicea(ngrid,l) condensation rate in layer l (kg/m2/s) |
---|
225 | ! zcondices(ngrid) condensation rate on the ground (kg/m2/s) |
---|
226 | ! zfallice(ngrid) flux of ice falling on surface (kg/m2/s) |
---|
227 | ! pdtc(ngrid,nlayermx) dT/dt due to phase changes (K/s) |
---|
228 | |
---|
229 | |
---|
230 | ! Tendencies initially set to 0 (except pdtc) |
---|
231 | DO l=1,nlayer |
---|
232 | DO ig=1,ngrid |
---|
233 | zcondicea(ig,l) = 0. |
---|
234 | pduc(ig,l) = 0 |
---|
235 | pdvc(ig,l) = 0 |
---|
236 | pdqc(ig,l,i_co2ice) = 0 |
---|
237 | END DO |
---|
238 | END DO |
---|
239 | |
---|
240 | DO ig=1,ngrid |
---|
241 | Mfallice(ig) = 0. |
---|
242 | zfallice(ig) = 0. |
---|
243 | zcondices(ig) = 0. |
---|
244 | pdtsrfc(ig) = 0. |
---|
245 | pdpsrf(ig) = 0. |
---|
246 | condsub(ig) = .false. |
---|
247 | zdiceco2(ig) = 0. |
---|
248 | ENDDO |
---|
249 | |
---|
250 | !----------------------------------------------------------------------- |
---|
251 | ! Atmospheric condensation |
---|
252 | |
---|
253 | |
---|
254 | ! Compute CO2 Volume mixing ratio |
---|
255 | ! ------------------------------- |
---|
256 | ! if (addn2) then |
---|
257 | ! DO l=1,nlayer |
---|
258 | ! DO ig=1,ngrid |
---|
259 | ! qco2=pq(ig,l,ico2)+pdq(ig,l,ico2)*ptimestep |
---|
260 | !! Mean air molecular mass = 1/(q(ico2)/m_co2 + (1-q(ico2))/m_noco2) |
---|
261 | ! mmean=1/(A*qco2 +B) |
---|
262 | ! vmr_co2(ig,l) = qco2*mmean/m_co2 |
---|
263 | ! ENDDO |
---|
264 | ! ENDDO |
---|
265 | ! else |
---|
266 | ! DO l=1,nlayer |
---|
267 | ! DO ig=1,ngrid |
---|
268 | ! vmr_co2(ig,l)=0.5 |
---|
269 | ! ENDDO |
---|
270 | ! ENDDO |
---|
271 | ! end if |
---|
272 | |
---|
273 | ! Forecast the atmospheric frost temperature 'ztcond' |
---|
274 | DO l=1,nlayer |
---|
275 | DO ig=1,ngrid |
---|
276 | ppco2=gfrac(igasco2)*pplay(ig,l) |
---|
277 | call get_tcond_co2(ppco2,ztcond(ig,l)) |
---|
278 | ENDDO |
---|
279 | ENDDO |
---|
280 | |
---|
281 | ! Initialize zq and zt at the beginning of the sub-timestep loop |
---|
282 | DO l=1,nlayer |
---|
283 | DO ig=1,ngrid |
---|
284 | zt(ig,l)=pt(ig,l) |
---|
285 | zq(ig,l,i_co2ice)=pq(ig,l,i_co2ice) |
---|
286 | IF( zq(ig,l,i_co2ice).lt.-1.e-6 ) THEN |
---|
287 | print*,'Uh-oh, zq = ',zq(ig,l,i_co2ice),'at ig,l=',ig,l |
---|
288 | if(l.eq.1)then |
---|
289 | print*,'Perhaps the atmosphere is collapsing on surface...?' |
---|
290 | endif |
---|
291 | END IF |
---|
292 | ENDDO |
---|
293 | ENDDO |
---|
294 | |
---|
295 | ! Calculate the mass of each atmospheric layer (kg.m-2) |
---|
296 | do ilay=1,nlayer |
---|
297 | do ig=1, ngrid |
---|
298 | masse(ig,ilay)=(pplev(ig,ilay) - pplev(ig,ilay+1)) /g |
---|
299 | end do |
---|
300 | end do |
---|
301 | |
---|
302 | ! ----------------------------------------------- |
---|
303 | ! START CONDENSATION/SEDIMENTATION SUB-TIME LOOP |
---|
304 | ! ----------------------------------------------- |
---|
305 | Ntime = 20 ! number of sub-timestep |
---|
306 | subptimestep = ptimestep/float(Ntime) |
---|
307 | |
---|
308 | DO it=1,Ntime |
---|
309 | |
---|
310 | ! Add the tendencies from other physical processes at each subtimstep |
---|
311 | DO l=1,nlayer |
---|
312 | DO ig=1,ngrid |
---|
313 | zt(ig,l) = zt(ig,l) + pdt(ig,l) * subptimestep |
---|
314 | zq(ig,l,i_co2ice) = zq(ig,l,i_co2ice) + pdq(ig,l,i_co2ice) * subptimestep |
---|
315 | END DO |
---|
316 | END DO |
---|
317 | |
---|
318 | |
---|
319 | ! Gravitational sedimentation |
---|
320 | |
---|
321 | ! sedimentation computed from radius computed from q |
---|
322 | ! assuming that the ice is splitted in Nmix particle |
---|
323 | do ilay=1,nlayer |
---|
324 | do ig=1, ngrid |
---|
325 | |
---|
326 | reff = CBRT( 3*zq(ig,ilay,i_co2ice)/( 4*Nmix_co2*pi*rho_co2 )) |
---|
327 | |
---|
328 | ! there should be a more elegant way of doing this... |
---|
329 | if(reff.lt.1.e-16) reff=1.e-16 |
---|
330 | |
---|
331 | ! update reffrad for radiative transfer |
---|
332 | if(reff.lt.reffradmin)then |
---|
333 | reffrad(ig,ilay,1)=reffradmin |
---|
334 | !print*,'reff below optical limit' |
---|
335 | elseif(reff.gt.reffradmax)then |
---|
336 | reffrad(ig,ilay,1)=reffradmax |
---|
337 | !print*,'reff above optical limit' |
---|
338 | else |
---|
339 | reffrad(ig,ilay,1)=reff |
---|
340 | endif |
---|
341 | |
---|
342 | call stokes & |
---|
343 | (pplev(ig,ilay),pt(ig,ilay), & |
---|
344 | reff,vstokes,rho_co2) |
---|
345 | |
---|
346 | !w(ig,ilay,i_co2ice) = 0.0 |
---|
347 | w(ig,ilay,i_co2ice) = vstokes * subptimestep * & |
---|
348 | pplev(ig,ilay)/(r*pt(ig,ilay)) |
---|
349 | |
---|
350 | end do |
---|
351 | end do |
---|
352 | |
---|
353 | ! Computing q after sedimentation |
---|
354 | |
---|
355 | call vlz_fi(ngrid,zq(1,1,i_co2ice),2.,masse,w(1,1,i_co2ice),wq) |
---|
356 | |
---|
357 | |
---|
358 | ! Progressively accumulating the flux to the ground |
---|
359 | ! Mfallice is the total amount of ice fallen to the ground |
---|
360 | do ig=1,ngrid |
---|
361 | Mfallice(ig) = Mfallice(ig) + wq(ig,i_co2ice) |
---|
362 | end do |
---|
363 | |
---|
364 | |
---|
365 | ! Condensation / sublimation in the atmosphere |
---|
366 | ! -------------------------------------------- |
---|
367 | ! (calculation of zcondicea, zfallice and pdtc) |
---|
368 | ! (MODIFICATIONS FOR EARLY MARS: falling heat neglected, condensation |
---|
369 | ! of CO2 into tracer i_co2ice) |
---|
370 | |
---|
371 | DO l=nlayer , 1, -1 |
---|
372 | DO ig=1,ngrid |
---|
373 | pdtc(ig,l)=0. |
---|
374 | |
---|
375 | ! tweak ccond if the gas is non-ideal |
---|
376 | if(nonideal)then |
---|
377 | ccond=cpp3D(ig,l)/(g*latcond) |
---|
378 | endif |
---|
379 | |
---|
380 | |
---|
381 | IF ((zt(ig,l).LT.ztcond(ig,l)).or.(zq(ig,l,i_co2ice).gt.0)) THEN |
---|
382 | condsub(ig)=.true. |
---|
383 | pdtc(ig,l) = (ztcond(ig,l) - zt(ig,l))/subptimestep |
---|
384 | pdqc(ig,l,i_co2ice) = pdtc(ig,l)*ccond*g |
---|
385 | |
---|
386 | ! Case when the ice from above sublimes entirely |
---|
387 | IF ((zq(ig,l,i_co2ice).lt.-pdqc(ig,l,i_co2ice)*subptimestep) & |
---|
388 | .AND. (zq(ig,l,i_co2ice).gt.0)) THEN |
---|
389 | |
---|
390 | pdqc(ig,l,i_co2ice) = -zq(ig,l,i_co2ice)/subptimestep |
---|
391 | pdtc(ig,l) =-zq(ig,l,i_co2ice)/(ccond*g*subptimestep) |
---|
392 | |
---|
393 | END IF |
---|
394 | |
---|
395 | ! Temperature and q after condensation |
---|
396 | zt(ig,l) = zt(ig,l) + pdtc(ig,l) * subptimestep |
---|
397 | zq(ig,l,i_co2ice) = zq(ig,l,i_co2ice) + pdqc(ig,l,i_co2ice) * subptimestep |
---|
398 | END IF |
---|
399 | |
---|
400 | ENDDO |
---|
401 | ENDDO |
---|
402 | ENDDO ! end of subtimestep loop |
---|
403 | |
---|
404 | ! Computing global tendencies after the subtimestep |
---|
405 | DO l=1,nlayer |
---|
406 | DO ig=1,ngrid |
---|
407 | pdtc(ig,l) = & |
---|
408 | (zt(ig,l) - (pt(ig,l) + pdt(ig,l)*ptimestep))/ptimestep |
---|
409 | pdqc(ig,l,i_co2ice) = & |
---|
410 | (zq(ig,l,i_co2ice)-(pq(ig,l,i_co2ice)+pdq(ig,l,i_co2ice)*ptimestep))/ptimestep |
---|
411 | END DO |
---|
412 | END DO |
---|
413 | DO ig=1,ngrid |
---|
414 | zfallice(ig) = Mfallice(ig)/ptimestep |
---|
415 | END DO |
---|
416 | |
---|
417 | |
---|
418 | !----------------------------------------------------------------------- |
---|
419 | ! Condensation/sublimation on the ground |
---|
420 | ! (calculation of zcondices and pdtsrfc) |
---|
421 | |
---|
422 | ! forecast of ground temperature ztsrf and frost temperature ztcondsol |
---|
423 | DO ig=1,ngrid |
---|
424 | ppco2=gfrac(igasco2)*pplay(ig,1) |
---|
425 | call get_tcond_co2(ppco2,ztcondsol(ig)) |
---|
426 | |
---|
427 | ztsrf(ig) = ptsrf(ig) |
---|
428 | |
---|
429 | if((ztsrf(ig).le.ztcondsol(ig)+2.0).and.(ngrid.eq.1))then |
---|
430 | print*,'CO2 is condensing on the surface in 1D. This atmosphere is doomed.' |
---|
431 | print*,'T_surf = ',ztsrf,'K' |
---|
432 | print*,'T_cond = ',ztcondsol,'K' |
---|
433 | open(116,file='surf_vals.out') |
---|
434 | write(116,*) 0.0, pplev(1,1), 0.0, 0.0 |
---|
435 | close(116) |
---|
436 | call abort |
---|
437 | endif |
---|
438 | |
---|
439 | ztsrf(ig) = ptsrf(ig) + pdtsrf(ig)*ptimestep |
---|
440 | |
---|
441 | ENDDO |
---|
442 | |
---|
443 | DO ig=1,ngrid |
---|
444 | IF(ig.GT.ngrid/2+1) THEN |
---|
445 | icap=2 |
---|
446 | ELSE |
---|
447 | icap=1 |
---|
448 | ENDIF |
---|
449 | |
---|
450 | ! Loop over where we have condensation / sublimation |
---|
451 | IF ((ztsrf(ig) .LT. ztcondsol(ig)) .OR. & ! ground condensation |
---|
452 | (zfallice(ig).NE.0.) .OR. & ! falling snow |
---|
453 | ((ztsrf(ig) .GT. ztcondsol(ig)) .AND. & ! ground sublimation |
---|
454 | ((piceco2(ig)+zfallice(ig)*ptimestep) .NE. 0.))) THEN |
---|
455 | condsub(ig) = .true. |
---|
456 | |
---|
457 | ! Condensation or partial sublimation of CO2 ice |
---|
458 | zcondices(ig)=pcapcal(ig)*(ztcondsol(ig)-ztsrf(ig)) & |
---|
459 | /(latcond*ptimestep) |
---|
460 | pdtsrfc(ig) = (ztcondsol(ig) - ztsrf(ig))/ptimestep |
---|
461 | |
---|
462 | ! If the entire CO_2 ice layer sublimes |
---|
463 | ! (including what has just condensed in the atmosphere) |
---|
464 | IF((piceco2(ig)/ptimestep+zfallice(ig)).LE. & |
---|
465 | -zcondices(ig))THEN |
---|
466 | zcondices(ig) = -piceco2(ig)/ptimestep - zfallice(ig) |
---|
467 | pdtsrfc(ig)=(latcond/pcapcal(ig))* & |
---|
468 | (zcondices(ig)) |
---|
469 | END IF |
---|
470 | |
---|
471 | ! Changing CO2 ice amount and pressure |
---|
472 | |
---|
473 | zdiceco2(ig) = zcondices(ig) + zfallice(ig) |
---|
474 | piceco2(ig) = piceco2(ig) + zdiceco2(ig)*ptimestep |
---|
475 | pdpsrf(ig) = -zdiceco2(ig)*g |
---|
476 | |
---|
477 | IF(ABS(pdpsrf(ig)*ptimestep).GT.pplev(ig,1)) THEN |
---|
478 | PRINT*,'STOP in condens' |
---|
479 | PRINT*,'condensing more than total mass' |
---|
480 | PRINT*,'Grid point ',ig |
---|
481 | PRINT*,'Ps = ',pplev(ig,1) |
---|
482 | PRINT*,'d Ps = ',pdpsrf(ig) |
---|
483 | STOP |
---|
484 | ENDIF |
---|
485 | END IF |
---|
486 | ENDDO |
---|
487 | |
---|
488 | ! Surface albedo and emissivity of the ground below the snow (emisref) |
---|
489 | ! -------------------------------------------------------------------- |
---|
490 | do ig =1,ngrid |
---|
491 | IF(ig.GT.ngrid/2+1) THEN |
---|
492 | icap=2 |
---|
493 | ELSE |
---|
494 | icap=1 |
---|
495 | ENDIF |
---|
496 | |
---|
497 | if(.not.piceco2(ig).ge.0.) THEN |
---|
498 | if(piceco2(ig).le.-1.e-8) print*, & |
---|
499 | 'WARNING in condense_co2cloud: piceco2(',ig,')=', piceco2(ig) |
---|
500 | piceco2(ig)=0. |
---|
501 | endif |
---|
502 | if (piceco2(ig).gt.0) then |
---|
503 | psolaralb(ig) = albedice(icap) |
---|
504 | emisref(ig) = emisice(icap) |
---|
505 | else |
---|
506 | psolaralb(ig) = albedodat(ig) |
---|
507 | emisref(ig) = emissiv |
---|
508 | pemisurf(ig) = emissiv |
---|
509 | end if |
---|
510 | end do |
---|
511 | |
---|
512 | return |
---|
513 | end subroutine condense_cloud |
---|
514 | |
---|
515 | !------------------------------------------------------------------------- |
---|
516 | subroutine get_tcond_co2(p,tcond) |
---|
517 | ! Calculates the condensation temperature for CO2 |
---|
518 | |
---|
519 | implicit none |
---|
520 | |
---|
521 | #include "callkeys.h" |
---|
522 | |
---|
523 | real p, peff, tcond |
---|
524 | real, parameter :: ptriple=518000.0 |
---|
525 | |
---|
526 | peff=p!/co2supsat |
---|
527 | |
---|
528 | if(peff.lt.ptriple)then |
---|
529 | tcond = (-3167.8)/(log(.01*peff)-23.23) ! Fanale's formula |
---|
530 | else |
---|
531 | tcond = 684.2-92.3*log(peff)+4.32*log(peff)**2 |
---|
532 | ! liquid-vapour transition (based on CRC handbook 2003 data) |
---|
533 | endif |
---|
534 | return |
---|
535 | |
---|
536 | end subroutine get_tcond_co2 |
---|