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