source: trunk/LMDZ.TITAN.old/libf/phytitan/n_ethane.F @ 3303

Last change on this file since 3303 was 1056, checked in by slebonnois, 11 years ago

SL: Titan runs ! see DOC/chantiers/commit_importants.log

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