source: trunk/LMDZ.TITAN/libf/phytitan/n_acethylene.F @ 537

Last change on this file since 537 was 175, checked in by slebonnois, 13 years ago

S.LEBONNOIS:

  • Revision majeure de la physique Titan => ajout des nuages version 10 bins (Jeremie Burgalat) Cette version reste a tester mais avec clouds=0, on reste sur l'ancienne.
  • Quelques ajouts dans la doc.
File size: 16.7 KB
Line 
1      subroutine n_acethylene(ngrid,nq,nbin,
2     *                     dt,pl,tl,aerad,
3     *                     q,qprime)
4
5      implicit none
6#include "dimensions.h"
7#include "microtab.h"
8#include "varmuphy.h"
9
10c  Arguments
11c  ---------
12
13      integer ngrid,nq,nbin
14
15      REAL dt              ! physical time step (s)
16      REAL pl(ngrid,nz)    ! pressure at each level (mbar)
17      REAL tl(ngrid,nz)    ! temperature at each level (K)
18      REAL aerad(nbin)     ! Radius array
19
20c    Tracers :
21      REAL q(ngrid,nz,nq)         ! tracer (kg/kg)
22      REAL qprime(ngrid,nz,nbin)  ! tracer (kg/kg)
23
24
25c     Local variables
26c     ---------------
27
28      integer ntyp
29      parameter (ntyp=3)
30
31      real n_aer(nz,nbin,ntyp)
32      real c2h2vap(nz)
33
34      integer itrac
35      integer ig,i,j,k,l,n   ! Loop integers
36      integer ilay,iq
37
38c  Treatment
39c  ---------
40
41      DO ig = 1 , NGRID
42
43c     Set up the aerosol array
44        do j = 1, ntyp
45          do k = 1, nbin
46            itrac = (j-1) * nbin + k
47            do l = 1, nz
48              n_aer(l,k,j) = max(q(ig,l,itrac),0.)
49            enddo
50          enddo
51        enddo
52
53c     Set up the methane vapor array
54        do l = 1, nz
55          c2h2vap(l) = q(ig,l,nq)
56        enddo
57
58        call nucleacond3(ngrid,nbin,dt,ig,pl,tl,aerad,
59     &                             n_aer,qprime,c2h2vap)
60
61c       Update q arrays
62        do j = 1, ntyp
63          do k = 1, nbin
64            itrac = (j-1) * nbin + k
65            do l = 1, nz
66              q(ig,l,itrac) = n_aer(l,k,j)
67            enddo
68          enddo
69        enddo
70
71c     Update methane vapor array
72        do l = 1, nz
73          q(ig,l,nq) = c2h2vap(l)
74        enddo
75
76      ENDDO
77
78      return
79      END
80
81****************************************************************
82      subroutine nucleacond3(ngrid,nbin,dt,ig,
83     *                      pl,tl,aerad,n_aer,qprime,c2h2vap)
84*                                                              *
85*     This routine updates species concentrations due          *
86*     to both nucleation and condensation-induced variations.  *
87*     Gain and loss rates associated to each one of these      *
88*     processes are computed separately in other routines.     *
89*                                                              *
90****************************************************************
91
92      implicit none
93#include "dimensions.h"
94#include "microtab.h"
95#include "varmuphy.h"
96
97
98      integer ng,nalt
99      parameter(ng=1,nalt=llm)
100
101
102      real lv
103
104      common/lheat/lv
105
106
107
108c  Arguments
109c  ---------
110
111      integer ngrid,nbin
112      integer ig
113      integer ntyp
114      parameter (ntyp=3)
115
116      real dt                    ! Global time step
117      real pl(ngrid,nz),tl(ngrid,nz)
118      real aerad(nbin)
119      real c2h2vap(nz)            ! Methane vapor mass mixing ratio (kg/m3)
120      real c2h2vap_old
121      real n_aer(nz,nbin,ntyp)  ! number concentrations of particle/each size bin
122      real qprime(ngrid,nz,nbin)  ! tracer (kg/kg)
123      REAL total1(nz),total11(nz),total2(nz),total22(nz)
124      REAL dmsm,mtot
125
126
127
128
129
130c  Local
131c  -----
132
133      integer i,j,k,l,n,iindice,iselec
134
135      real dQc           ! Amount of condensed methane (kg/m3) during timestep
136      real*8 sat_ratio     ! Methane saturation ratio over liquid
137      real*8 sat_ratmix  ! Methane saturation ratio over liquid
138      real*8 pc2h2         ! Methane partial pressure (Pa)
139      real qsat          ! Methane mass mixing ratio at saturation (kg/kg of air)
140      real qsatmix          ! Methane mass mixing ratio at saturation (kg/kg of air)
141      real*8 rate(nbin)    ! Heterogeneous Nucleation rate (s-1)
142      real*8 elim         
143
144      real nsav(nbin,ntyp)
145      real dn(nbin,ntyp)
146      real rad(nbin)     ! Radius of droplets in each size bin
147      real*8 gr(nbin)      ! Growth rate in each bin
148      real radius        ! Radius of droplets after growth
149      real Qs            ! Mass of condensate required to reach saturation
150      real newsat
151      real vol(nbin)
152
153      real press
154      real sig3,temp,seq(nbin)
155      real Ctot,up,dwn,newvap,gltot
156
157      real temp0,temp1,temp2,last_temp
158      real qsat1,sat_ratio1,tempf(0:10),sat_ratiof(0:10)
159      real rho_a,cap
160      real tempref
161      real xtime,xtime_prime
162
163
164c     Variables for latent heat release
165      real lw,cpp
166      data lw / 581.e+3/
167      data cpp/1050./ ! pour etre cohérent avec le reste...
168      save lw,cpp
169
170
171c  Treatment
172c  ---------
173      do i = 1, nbin
174        vol(i) = 4./3. * pi * aerad(i)**3.
175      enddo
176
177      do l = 1, nz
178      total1(l)=0. !solide
179      do k = 1, nbin
180      total1(l)=total1(l)+n_aer(l,k,2)*rhoi_c2h2
181      enddo
182      total2(l)=c2h2vap(l)
183      enddo
184
185
186c     Start loop over heights
187      DO 100 l = 1, nz
188
189        iindice=0                ! mettre l'indice à 0
190
191        temp   = tl(ig,l)
192        press  = pl(ig,l)
193        tempref=temp
194
195c       Save the values of the particle arrays before condensation
196        do j = 1, ntyp
197          do i = 1, nbin
198             nsav(i,j) = n_aer(l,i,j)
199          enddo
200        enddo
201
202
203 99     continue
204
205
206        call c2h2sat(temp,press,qsat)
207        qsatmix=qsat
208 
209c  quantité  pmixc2h2(l) déjà calculé dans cnuages.F et passé dans un common
210 
211c       Get the partial presure of methane vapor and its saturation ratio
212        pc2h2      = c2h2vap(l) * (Mn2/Mc2h2) * press
213        sat_ratio  = c2h2vap(l) / qsat
214        sat_ratmix = c2h2vap(l) / qsatmix
215
216c       Get the rates of nucleation
217        call nuclea3(nbin,aerad,pc2h2,temp,sat_ratio,rate)
218
219c       Get the growth rates of condensation/sublimation
220        up   = c2h2vap(l)
221        dwn  = 1.
222        Ctot = c2h2vap(l)
223        DO i = 1, nbin
224
225        if (n_aer(l,i,3).eq.0) then
226         rad(i) = aerad(i)
227        else
228         rad(i) = ((n_aer(l,i,2)/n_aer(l,i,3) +
229     &   qprime(ig,l,i)/n_aer(l,i,3)
230     &   +vol(i))*0.75/pi)**(1./3.)
231        endif
232
233
234*       Equilibrium saturation ratio (due to curvature effect)
235        seq(i) = exp( 2.*sig3(temp)*Mc2h2 /(rhoi_c2h2*rgp*temp*rad(i)))
236
237        call growthrate3(dt,temp,press,pc2h2,
238     &                   sat_ratmix,seq(i),rad(i),gr(i))
239
240        up = up + dt * gr(i) * 4. * pi * rhoi_c2h2 * rad(i) * seq(i)
241     *                 * nsav(i,3)
242        dwn= dwn+ dt * gr(i) * 4. * pi * rhoi_c2h2 * rad(i) / qsat
243     *                 * nsav(i,3)
244        Ctot= Ctot + rhoi_c2h2 * nsav(i,2)
245
246        ENDDO
247
248        newvap = min(up/dwn,Ctot)
249        newvap = max(newvap,0.)
250
251        gltot = 0.
252        DO i = 1, nbin
253          gr(i)  = gr(i) * ( newvap/qsat - seq(i) )
254          if(nsav(i,2).le.0. .and. gr(i).le.0.) then
255              n_aer(l,i,2) = 0.
256          else
257          n_aer(l,i,2) = nsav(i,2) + dt * gr(i) * 4. * pi * rad(i)
258     *                               * n_aer(l,i,3)
259          if (n_aer(l,i,2).le.0.) then
260            n_aer(l,i,1) = n_aer(l,i,1) + n_aer(l,i,3)
261            n_aer(l,i,2) = 0.
262            n_aer(l,i,3) = 0.
263          endif
264          gltot=n_aer(l,i,2)*rhoi_c2h2+gltot
265          endif
266
267        ENDDO
268   
269c       Determine the mass of exchanged methane
270
271        dQc = 0.
272        DO i = 1, nbin
273          dQc = dQc - rhoi_c2h2 * ( n_aer(l,i,2) - nsav(i,2) )
274        ENDDO
275
276c       Update the methane vapor mixing ratio implied by
277c       the cond/eva processes.
278
279
280c       Arrays resetted to their initial value before condensation
281        do j = 1, ntyp
282          do i = 1, nbin
283            dn(i,j)      = n_aer(l,i,j) - nsav(i,j)
284            n_aer(l,i,j) = nsav(i,j)
285          enddo
286        enddo
287
288c       Update the c arrays.
289c       nucleation & cond/eva tendencies added together.
290
291        do i=1,nbin
292          elim         = dt * rate(i)
293          n_aer(l,i,1) = n_aer(l,i,1) / (1.+elim)
294          n_aer(l,i,3) = n_aer(l,i,3) + elim * n_aer(l,i,1) + dn(i,3)
295          n_aer(l,i,1) = n_aer(l,i,1) + dn(i,1)
296          n_aer(l,i,2) = n_aer(l,i,2) + dn(i,2)
297          if(n_aer(l,i,2).lt.0.) n_aer(l,i,2)=0.
298        enddo
299
300        dQc = 0.
301        DO i = 1, nbin
302          dQc = dQc - rhoi_c2h2 * ( n_aer(l,i,2) - nsav(i,2) )
303        ENDDO
304
305
306        c2h2vap(l)  = c2h2vap(l) + dQc
307
308100   CONTINUE
309
310      do l = 1, nz
311      total11(l)=0.
312      do k = 1, nbin
313      total11(l)=total11(l)+n_aer(l,k,2)*rhoi_c2h2
314      enddo
315      total22(l)=c2h2vap(l)
316      enddo
317
318      return
319      end
320
321
322*******************************************************
323*                                                     *
324      subroutine nuclea3(nbin,aerad,pc2h2,temp,sat,nucrate)
325*   This subroutine computes the nucleation rate      *
326*   as given in Pruppacher & Klett (1978) in the      *
327*   case of water ice forming on a solid substrate.   *
328*     Definition refined by Keese (jgr,1989)          *
329*                                                     *
330*******************************************************
331
332      implicit none
333#include "dimensions.h"
334#include "microtab.h"
335#include "varmuphy.h"
336
337      integer nbin
338      real aerad(nbin)
339
340      real*8 nucrate(nbin)
341      real*8 pc2h2
342      real   temp
343      real*8 sat
344
345      integer l,i
346      real*8 nc2h2
347      real sig3            ! Water-ice/air surface tension  (N.m)
348      real*8 rstar        ! Radius of the critical germ (m)
349      real*8 gstar        ! # of molecules forming a critical embryo
350      real*8 x            ! Ratio rstar/radius of the nucleating dust particle
351      real fistar         ! Activation energy required to form a critical embryo (J)
352      real*8 zeldov       ! Zeldovitch factor (no dim)
353      real*8 fshape3       ! function defined at the end of the file
354      real*8 deltaf
355
356      real nus
357      data nus/1.e+13/       ! Jump frequency of a molecule (s-1)
358      real m0
359      data m0/4.31894e-26/     ! Weight of a methane molecule (kg)
360      real vo1
361      data vo1/4.22764e-5/    ! Volume molaire (masse molaire/masse volumique = MolWt/LDEN)
362      real desorp
363      data desorp/0.288e-19/ ! Activation energy for desorption of water on a dust-like substrate (J/molecule)
364      real surfdif
365      data surfdif/0.288e-20/! Estimated activation energy for surface diffusion of water molecules (J/molecule)
366
367      IF (sat .GT. 1.) then    ! minimum condition to activate nucleation
368
369        nc2h2    = pc2h2 / kbz / temp
370        rstar  = 2. * sig3(temp) * vo1 / (rgp*temp*log(sat))
371        gstar  = 4. * nav * pi * (rstar**3) / (3.*vo1)
372c       Loop over size bins
373        do i=1,nbin
374          x      = aerad(i) / rstar
375          x      = aerad(imono) / rstar  ! r(5)=monomere
376          fistar = (4./3.*pi) * sig3(temp) * (rstar**2.)
377     &     *fshape3(mtetac2h2,x)
378          deltaf = min( max((2.*desorp-surfdif-fistar)/(kbz*temp)
379     &           , -100.), 100.)
380          if (deltaf.eq.-100.) then
381            nucrate(i) = 0.
382          else
383            zeldov = sqrt ( fistar / (3.*pi*kbz*temp*(gstar**2.)) )
384            nucrate(i)  = zeldov * kbz* temp * rstar**2.
385     &                  * 4. * pi * ( nc2h2*aerad(i) )**2.
386     &                  / ( fshape3(mtetac2h2,x) * nus * m0 )
387     &                  * dexp(deltaf)
388
389
390            if(i.gt.imono)  nucrate(i)= zeldov * kbz* temp * rstar**2.
391     &          * 4. * pi * vrat_e**(i-imono)*(nc2h2*aerad(imono) )**2.
392     &                  / (fshape3(mtetac2h2,x) * nus * m0 )
393     &                  * dexp(deltaf)
394
395          endif
396        enddo
397      ELSE
398        do i=1,nbin
399          nucrate(i) = 0.
400        enddo
401
402      ENDIF
403
404      return
405      end
406
407******************************************************************
408        subroutine c2h2sat(t,p,qsat)
409*                                                                 *
410* cette fonction calcule la pression de vapeur saturante de l'    *
411* ethane a une altitude donnee z par Reid et al., p657            *
412*                                                                 *
413* Compatible avec Barth et al., dans l'intervalle 30-90K          *
414*                                                                 *
415*                                                                 *
416******************************************************************
417
418        real rgp
419        data rgp/8.3143/
420
421* declaration des variables internes
422* ----------------------------------
423
424        real qsat,t,p
425         
426
427        a=-6.90128
428        b=1.26873
429        c=-2.09113
430        d=-2.75601
431        pc=61.4*1.013e5
432        tc=308.3
433
434        x=(1.-t/tc)
435        if(x.gt. 0.) qsat=(1-x)**(-1)*(a*x+b*x**1.5+c*x**3.+d*x**6.)
436        if(x.le. 0.) qsat=a*x/abs(1.-x)     ! approx pour  t > tc
437        qsat=pc*exp(qsat)
438
439        qsat=qsat* 26.0 / (28.0*p)  ! kg/kg
440
441        return
442        end
443
444c=======================================================================
445      subroutine growthrate3(timestep,temp,press,pc2h2,sat,seq,r,Cste)
446c
447c     Determination of the droplet growth rate
448c
449c=======================================================================
450
451      IMPLICIT NONE
452#include "dimensions.h"
453#include "microtab.h"
454#include "varmuphy.h"
455
456c-----------------------------------------------------------------------
457C   DECLARATIONS:
458c   -------------
459
460      common/lheat/Lv
461
462c
463c   arguments:
464c   ----------
465
466      REAL timestep
467      REAL temp    ! temperature in the middle of the layer (K)
468      REAL press    ! pressure in the middle of the layer (K)
469      REAL*8 pc2h2 ! Methane vapor partial pressure (Pa)
470      REAL*8 sat  ! saturation ratio
471      REAL r    ! crystal radius before condensation (m)
472      REAL seq  ! Equilibrium saturation ratio
473
474c   local:
475c   ------
476
477      REAL psat
478      REAL moln2,molc2h2
479      REAL To,tc2h2,wc2h2       ! Reid et al., (eq 7-9.4 + Appendix compound  [168])
480      REAL fte
481
482c     Effective gas molecular radius (m)
483      data moln2/1.75e-10/   ! N2
484c     Effective gas molecular radius (m)
485      data molc2h2/2.015e-10/   ! C2H2
486c     Temperature critique  + omega
487      data tc2h2/308.3/
488      data wc2h2/19.0e-2/
489
490      REAL k,Lv                 
491      REAL knudsen           ! Knudsen number (gas mean free path/particle radius)
492      REAL a,Dv,lambda,Rk,Rd ! Intermediate computations for growth rate
493      REAL*8 Cste
494
495c-----------------------------------------------------------------------
496c      Ice particle growth rate by diffusion/impegement of molecules
497c                r.dr/dt = (S-Seq) / (Seq*Rk+Rd)
498c        with r the crystal radius, Rk and Rd the resistances due to
499c        latent heat release and to vapor diffusion respectively
500c-----------------------------------------------------------------------
501
502      psat = pc2h2 / sat
503
504c     - Thermal conductibility of N2
505     
506      k = ( 2.857e-2 * temp - 0.5428  ) * 4.184e-3
507     
508     
509c     - Latent heat of c2h2 (J.kg-1)
510      Lv =581.e3                      ! eq (7-9.4) Reid et al.
511      fte=(1.-temp/tc2h2)
512      if (fte.le.1.e-3)  fte=1.e-3
513      Lv=8.314*tc2h2*(7.08*fte**0.354+10.95*wc2h2*fte**0.456)/26.e-3
514
515     
516
517c     - Constant to compute gas mean free path
518c     l= (T/P)*a, with a = (  0.707*8.31/(4*pi*molrad**2 * avogadro))
519
520      a = 0.707*rgp/(4 * pi* (moln2*1.e10)**2  * (nav*1.e-20))
521
522c     - Compute Dv, methane vapor diffusion coefficient
523c       accounting for both kinetic and continuum regime of diffusion,
524c       the nature of which depending on the Knudsen number.
525
526      Dv = 1./3. * sqrt( 8*rgp*temp/(pi*Mc2h2) )* (kbz*1.e20) * temp/
527     & (pi*press*(moln2*1.e10+molc2h2*1.e10)**2 * sqrt(1.+Mc2h2/Mn2))
528
529      knudsen = temp / press * a / r
530      lambda  = (1.333+0.71/knudsen) / (1.+1./knudsen)
531      Dv      = Dv / (1. + lambda * knudsen)
532
533c     - Compute Rk
534      Rk = Lv**2 * rhoi_c2h2 * Mc2h2 / (k*rgp*temp**2.)
535*     print*,'Cste Rk :',Lv,k,rgp,t
536
537c     - Compute Rd
538      Rd = rgp * temp *rhoi_c2h2 / (Dv*psat*Mc2h2)
539*     print*,'Cste Rd :',Dv,psat,Mc2h2
540
541c     - Compute:      rdr/dt = Cste * (S-Seq)
542      Cste = 1. / (seq*Rk+Rd)
543*     print*,'Cste Cste :',seq,Rk,Rd
544
545
546      RETURN
547      END
548
549
550*********************************************************
551      real function sig3(t)
552*   this function computes the surface tension (N.m)   *
553*   between acethylene and air as a function of temp.    *
554*********************************************************
555
556      real t
557      pc=61.4*1.01325e5
558      tc=308.3
559      tb=188.4
560      tr=t/tc
561      tbr=tb/tc
562      if(t.gt.308.0) then
563         tr=308./tc
564      endif
565
566
567
568      sig3=0.1196*(1.+(tbr*alog(pc/1.01325))/(1.-tbr))-0.279
569      sig3=pc**(2./3.)*tc**(1./3.)*sig3*(1.-tr)**(11./9.)
570      sig3=sig3*1.e-8
571
572      return
573      end
574
575*********************************************************
576      real*8 function fshape3(cost,rap)
577*        function computing the f(m,x) factor           *
578* related to energy required to form a critical embryo  *
579*********************************************************
580
581      implicit none
582
583      real cost
584      real*8 rap
585      real*8 phi
586      real*8 a,b,c
587
588
589      phi = sqrt( 1. - 2.*cost*rap + rap**2 )
590      a = 1. + ( (1.-cost*rap)/phi )**3
591      b = (rap**3) * (2.-3.*(rap-cost)/phi+((rap-cost)/phi)**3)
592      c = 3. * cost * (rap**2) * ((rap-cost)/phi-1.)
593
594      fshape3 = 0.5*(a+b+c)
595
596      if (rap.gt.3000.) fshape3 = ((2.+cost)*(1.-cost)**2)/4.
597
598      return
599      end
600
Note: See TracBrowser for help on using the repository browser.