source: trunk/LMDZ.MARS/libf/phymars/massflowrateCO2.F @ 1818

Last change on this file since 1818 was 1816, checked in by jaudouard, 7 years ago

Commit for CO2 clouds microphysics.

  • Property svn:executable set to *
File size: 17.0 KB
Line 
1c=======================================================================
2      subroutine massflowrateCO2(P,T,Sat,Radius,Matm,Ic)
3c
4c     Determination of the mass transfer rate for CO2 condensation &
5c     sublimation
6c
7c     newton-raphson method
8c     CLASSICAL  (no SF etc.)
9c
10c   inputs: Pressure (P), Temperature (T), saturation ratio (Sat),
11c           particle radius (Radius), molecular mass of the atm (Matm)
12c   output:  MASS FLUX Ic
13c
14c   Authors: C. Listowski (2014) then J. Audouard (2016-2017)
15c=======================================================================
16      USE comcstfi_h
17
18      implicit none
19
20      include "microphys.h"
21
22
23
24c   arguments: INPUT
25c   ----------
26      REAL T,Matm
27      REAL*8 SAT
28      real P
29      DOUBLE PRECISION Radius
30c   arguments: OUTPUT
31c   ----------
32      DOUBLE PRECISION   Ic
33c   Local Variables
34c   ----------
35      DOUBLE PRECISION   Tcm
36      DOUBLE PRECISION   T_inf, T_sup, T_dT
37      DOUBLE PRECISION   C0,C1,C2
38      DOUBLE PRECISION   kmix,Lsub,cond
39      DOUBLE PRECISION   rtsafe
40      DOUBLE PRECISION   left, fval, dfval
41c    function for newton-raphson iterative method
42c    --------------------------
43      EXTERNAL classical
44         
45      Tcm = 110!dble(T)    ! initialize pourquoi 0 et pas t(i)
46      T_inf = 0d0
47      T_sup = 800d0
48      T_dT = 0.001  ! precision - mettre petit et limiter nb iteration?
49     
50666   CONTINUE
51
52      call  coefffunc(P,T,Sat,Radius,Matm,kmix,Lsub,C0,C1,C2)
53
54      if (isnan(C0) .eqv. .true.)  C0=0d0   
55c     FIND SURFACE TEMPERATURE (Tc) : iteration sur t
56      cond = 4.*pi*Radius*kmix
57      Tcm = rtsafe(classical,T_inf,T_sup,T_dT,Radius,C0,C1,C2)
58      if (Tcm.LE.0d0 .or. isnan(Tcm)) then ! unsignificant cases where S<<<Seq and Ncores <<1e-10
59         Tcm = T
60      endif
61c     THEN COMPUTE MASS FLUX Ic from FINAL Tsurface (Tc)
62      Ic = (Tcm-T)
63      Ic = cond*Ic/Lsub
64c regarder de combien varie la solution Ic entre Tcm et Tcm+T_dT
65 
66      RETURN
67
68      END
69
70
71c****************************************************************
72      FUNCTION rtsafe(funcd,x1,x2,xacc,Radius,C0,C1,C2)
73*
74*     Newton Raphsen routine (Numerical Recipe)
75*
76c****************************************************************
77
78      implicit none
79
80      INTEGER MAXIT
81       DOUBLE PRECISION x1,x2,xacc
82      DOUBLE PRECISION rtsafe
83      DOUBLE PRECISION Radius
84      DOUBLE PRECISION C0,C1,C2
85
86
87      EXTERNAL funcd
88
89      PARAMETER (MAXIT=10000)   !Maximum allowed number of iterations. Using a combination of Newton-Raphson and bisection,
90                                !find the root of a function bracketed between x1 and x2. The root, returned as the function value rtsafe,
91                                !will be refined until its accuracy is known within !!±xacc. funcd is a user-supplied subroutine which
92                                !returns both the function value and the first derivative of the function.
93
94      INTEGER j
95      DOUBLE PRECISION df,dx,dxold
96      DOUBLE PRECISION f,fh,fl,temp,xh,xl
97
98      call funcd(x1,fl,df,C0,C1,C2)
99      call funcd(x2,fh,df,C0,C1,C2)
100
101
102      if ((fl.gt.0..and.fh.gt.0.).or.(fl.lt.0..and.fh.lt.0.) ) then
103         x1=0d0
104         x2=1200d0         
105         call funcd(x1,fl,df,C0,C1,C2)
106         call funcd(x2,fh,df,C0,C1,C2)
107c         write(*,*) 'root must be bracketed in rtsafe'
108      endif
109
110
111      if (fl.eq.0.) then
112         rtsafe=x1
113         return
114      else if (fh.eq.0.) then
115         rtsafe=x2
116         return
117      else if (fl.lt.0.) then   !Orient the search so that f(xl) < 0.                                                                   
118         xl=x1
119         xh=x2
120      else
121        xh=x1
122        xl=x2
123      endif
124      rtsafe = .5*(x1+x2)       !Initialize the guess for root,
125      dxold  = abs(x2-x1)       !the stepsize before last,
126      dx     = dxold            ! and the last step.
127      call funcd(rtsafe,f,df,C0,C1,C2)
128      DO 11 j=1,MAXIT           !Loop over allowed iterations.
129
130         if (((rtsafe-xh)*df-f)*((rtsafe-xl)*df-f).gt.0.   ! Bisect if Newton out of range
131     *    .or. abs(2.*f).gt.abs(dxold*df) ) then               ! or not decreasing fst enough
132
133             dxold=dx
134             dx=0.5*(xh-xl)
135             rtsafe=xl+dx
136             if (xl.eq.rtsafe) return                   !Change in root is negligible. Newton step acceptable. Take it.
137         else
138            dxold=dx
139            dx=f/df
140            temp=rtsafe
141            rtsafe=rtsafe-dx
142            if(temp.eq.rtsafe) return
143         endif
144
145        if(abs(dx).lt.xacc) return !Convergence criterion. The one new function evaluation per iteration. Maintain the bracket on the root.
146         call funcd(rtsafe,f,df,C0,C1,C2)
147        if(f.lt.0.) then
148         xl=rtsafe
149        else
150         xh=rtsafe
151        endif
15211    ENDDO
153
154      rtsafe=0d0
155      write(*,*) 'rtsafe exceeding maximum iterations,Tcm=',rtsafe
156      return
157
158      END
159
160
161c********************************************************************************
162      subroutine classical(x,f,df,C0,C1,C2)
163c     Function to give as input to RTSAFE (NEWTON-RAPHOEN)
164c********************************************************************************
165
166      implicit none
167
168      DOUBLE PRECISION   x
169      DOUBLE PRECISION  C0,C1,C2
170      DOUBLE PRECISION f
171      DOUBLE PRECISION  df
172
173      f   = x + C0*exp(C1*x) - C2   ! start f
174      df  = 1. + C0*C1*exp(C1*x)    ! start df
175 
176      return
177      END 
178
179c********************************************************************************
180
181      subroutine coefffunc(P,T,S,rc,Matm,kmix,Lsub,C0,C1,C2)
182
183c********************************************************************************
184c defini la fonction eq 6 papier 2014  (Listowski et al., 2014)
185      use tracer_mod, only: rho_ice_co2
186      USE comcstfi_h
187
188      implicit none
189      include "microphys.h"
190
191c   arguments: INPUT
192c   ----------------
193      REAL P
194      real T
195      REAL*8 S
196      double precision rc
197      REAL Matm !g.mol-1 ( = mmean(ig,l) )
198
199c   local:
200c   ------
201      DOUBLE PRECISION Cpatm,Cpn2,Cpco2
202      DOUBLE PRECISION psat, xinf, pco2
203      DOUBLE PRECISION Dv     
204      DOUBLE PRECISION l0,l1,l2,l3,l4           
205      DOUBLE PRECISION knudsen, a, lambda      ! F and S correction
206      DOUBLE PRECISION Ak                       ! kelvin factor   
207      DOUBLE PRECISION vthatm,lpmt,rhoatm, vthco2 ! for Kn,th
208
209c   arguments: OUTPUT
210c   ----------
211      DOUBLE PRECISION  C0,C1,C2
212      DOUBLE PRECISION  kmix,Lsub
213
214c     DEFINE heat cap. J.kg-1.K-1 and To
215
216      data Cpco2/0.7e3/
217      data Cpn2/1e3/
218
219      kmix = 0d0
220      Lsub = 0d0
221
222      C0 = 0d0
223      C1 = 0d0
224      C2 = 0d0
225
226c     Equilibirum pressure over a flat surface
227      psat = 1.382 * 1.00e12 * exp(-3182.48/dble(T))  ! (Pa)
228c     Compute transport coefficient
229      pco2 = psat * dble(S)
230c     Latent heat of sublimation if CO2  co2 (J.kg-1)
231c     version Azreg_Ainou (J/kg) :
232      l0=595594.     
233      l1=903.111     
234      l2=-11.5959   
235      l3=0.0528288
236      l4=-0.000103183
237      Lsub = l0 + l1 * dble(T) + l2 * dble(T)**2 + l3 *
238     &     dble(T)**3 + l4 * dble(T)**4 ! J/kg
239c     atmospheric density
240      rhoatm = dble(P*Matm)/(rgp*dble(T))   ! g.m-3
241      rhoatm = rhoatm * 1.00e-3 !kg.m-3
242      call  KthMixNEW(kmix,T,pco2/dble(P),rhoatm) ! compute thermal cond of mixture co2/N2
243      call  Diffcoeff(P, T, Dv) 
244
245      Dv = Dv * 1.00e-4         !!! cm2.s-1  to m2.s-1
246
247c     ----- FS correction for Diff
248      vthco2  = sqrt(8d0*kbz*dble(T)/(dble(pi) * mco2/nav)) ! units OK: m.s-1
249      knudsen = 3*Dv / (vthco2 * rc)
250      lambda  = (1.333+0.71/knudsen) / (1.+1./knudsen) ! pas adaptée, Dahneke 1983? en fait si (Monschick&Black)
251      Dv      = Dv / (1. + lambda * knudsen)
252c     ----- FS correction for Kth
253      vthatm = sqrt(8d0*kbz*dble(T)/(pi * 1.00e-3*dble(Matm)/nav)) ! Matm/nav = mass of "air molecule" in G , *1e-3 --> kg
254      Cpatm = Cpco2 * pco2/dble(P) + Cpn2 * (1d0 - pco2/dble(P)) !J.kg-1.K-1
255      lpmt = 3 * kmix / (rhoatm * vthatm * (Cpatm - 0.5*rgp/
256     &     (dble(Matm)*1.00e-3))) ! mean free path related to heat transfer
257      knudsen = lpmt / rc
258      lambda  = (1.333+0.71/knudsen) / (1.+1./knudsen) ! pas adaptée, Dahneke 1983? en fait si (Monschick&Black)
259      kmix    = kmix /  (1. + lambda * knudsen)
260c     --------------------- ASSIGN coeff values for FUNCTION
261      xinf = dble(S) * psat / dble(P)
262      Ak = exp(2d0*sigco2*mco2/(rgp* dble(rho_ice_co2*T* rc) ))
263      C0 = mco2*Dv*psat*Lsub/(rgp*dble(T)*kmix)*Ak*exp(-Lsub*mco2/
264     &     (rgp*dble(T)))
265      C1 = Lsub*mco2/(rgp*dble(T)**2)
266      C2 = dble(T) + dble(P)*mco2*Dv*Lsub*xinf/(kmix*rgp*dble(T))
267      RETURN
268      END
269
270
271c======================================================================
272      subroutine Diffcoeff(P, T, Diff)
273c     Compute diffusion coefficient CO2/N2
274c     cited in Ilona's lecture - from Reid et al. 1987
275c======================================================================
276       IMPLICIT NONE
277
278       include "microphys.h"
279
280c      arguments
281c     -----------
282     
283      REAL P
284      REAL Pbar                     !!! has to be in bar for the formula
285      REAL T
286     
287c     output
288c     -----------
289
290      DOUBLE PRECISION Diff
291   
292c      local
293c     -----------
294
295      DOUBLE PRECISION dva, dvb, Mab  ! Mab has to be in g.mol-1
296       
297        Pbar = P * 1d-5
298   
299        Mab = 2. / ( 1./mn2 + 1./mco2 ) * 1000.
300
301        dva = 26.9        ! diffusion volume of CO2,  Reid et al. 1987 (cited in Ilona's lecture)
302        dvb = 18.5        ! diffusion volume of N2
303   
304        Diff  = 0.00143 * dble(T)**(1.75) / (dble(Pbar) * sqrt(Mab)
305     &       * (dble(dva)**(1./3.) + dble(dvb)**(1./3.))**2.) !!! in cm2.s-1
306 
307       RETURN
308
309       END
310
311
312c======================================================================
313
314         subroutine KthMixNEW(Kthmix,T,x,rho)
315
316c        Compute thermal conductivity of CO2/N2 mixture
317c         (***WITHOUT*** USE OF VISCOSITY)
318
319c          (Mason & Saxena, 1958 - Wassiljeva 1904)
320c======================================================================
321
322       implicit none
323       
324       include "microphys.h"
325c      arguments
326c     -----------
327         
328         REAL T
329         DOUBLE PRECISION x
330         DOUBLE PRECISION rho !kg.m-3
331
332c     outputs
333c     -----------
334
335         DOUBLE PRECISION Kthmix
336
337c     local
338c    ------------
339
340         DOUBLE PRECISION x1,x2
341
342         DOUBLE PRECISION  Tc1, Tc2, Pc1, Pc2
343
344         DOUBLE PRECISION  A12, A11, A22, A21
345
346         DOUBLE PRECISION  Gamma1, Gamma2, M1, M2
347         DOUBLE PRECISION  lambda_trans1, lambda_trans2,epsilon
348
349         DOUBLE PRECISION  kco2, kn2
350
351      x1 = x
352      x2 = 1d0 - x
353
354      M1 = mco2
355      M2 = mn2
356
357      Tc1 =  304.1282 !(Scalabrin et al. 2006)
358      Tc2 =  126.192  ! (Lemmon & Jacobsen 2003)
359
360      Pc1 =  73.773   ! (bars)
361      Pc2 =  33.958   ! (bars)
362   
363      Gamma1 = 210.*(Tc1*M1**(3.)/Pc1**(4.))**(1./6.)
364      Gamma2 = 210.*(Tc2*M2**(3.)/Pc2**(4.))**(1./6.)
365
366c Translational conductivities
367
368      lambda_trans1 = ( exp(0.0464 * T/Tc1) - exp(-0.2412 * T/Tc1) )
369     &                                                          /Gamma1
370
371      lambda_trans2 = ( exp(0.0464 * T/Tc2) - exp(-0.2412 * T/Tc2) )
372     &                                                          /Gamma2
373     
374c     Coefficient of Mason and Saxena
375      epsilon = 1.
376
377      A11 = 1.
378         
379      A22 = 1.
380
381      A12 = epsilon * (1. + sqrt(lambda_trans1/lambda_trans2)*
382     &                    (M1/M2)**(1./4.))**(2.) / sqrt(8*(1.+ M1/M2))
383
384      A21 = epsilon * (1. + sqrt(lambda_trans2/lambda_trans1)*
385     &                    (M2/M1)**(1./4.))**(2.) / sqrt(8*(1.+ M2/M1))
386
387c     INDIVIDUAL COND.
388
389         call KthCO2Scalab(kco2,T,rho)
390         call KthN2LemJac(kn2,T,rho)
391
392c     MIXTURE COND.
393        Kthmix = kco2*x1 /(x1*A11 + x2*A12) + kn2*x2 /(x1*A21 + x2*A22)
394        Kthmix = Kthmix*1e-3   ! from mW.m-1.K-1 to  W.m-1.K-1
395        RETURN
396
397        END
398
399c======================================================================
400         subroutine KthN2LemJac(kthn2,T,rho)
401c        Compute thermal cond of N2 (Lemmon and Jacobsen, 2003)
402cWITH viscosity
403c======================================================================
404
405       implicit none
406
407        include "microphys.h"
408c        include "microphysCO2.h"
409
410
411c      arguments
412c     -----------
413         
414         REAL T
415         DOUBLE PRECISION rho !kg.m-3
416
417c     outputs
418c     -----------
419
420         DOUBLE PRECISION kthn2
421
422c     local
423c    ------------
424
425        DOUBLE PRECISION g1,g2,g3,g4,g5,g6,g7,g8,g9,g10
426        DOUBLE PRECISION h1,h2,h3,h4,h5,h6,h7,h8,h9,h10
427        DOUBLE PRECISION n1,n2,n3,n4,n5,n6,n7,n8,n9,n10
428        DOUBLE PRECISION d4,d5,d6,d7,d8,d9
429        DOUBLE PRECISION l4,l5,l6,l7,l8,l9
430        DOUBLE PRECISION t2,t3,t4,t5,t6,t7,t8,t9
431        DOUBLE PRECISION gamma4,gamma5,gamma6,gamma7,gamma8,gamma9
432
433        DOUBLE PRECISION Tc,rhoc
434
435        DOUBLE PRECISION tau, delta
436
437        DOUBLE PRECISION visco
438
439        DOUBLE PRECISION k1, k2
440
441         N1 = 1.511d0
442         N2 = 2.117d0
443         N3 = -3.332d0
444
445         N4 = 8.862
446         N5 = 31.11
447         N6 = -73.13
448         N7 = 20.03
449         N8 = -0.7096
450         N9 = 0.2672
451
452         t2 = -1.0d0
453         t3 = -0.7d0
454         t4 = 0.0d0
455         t5 = 0.03
456         t6 = 0.2
457         t7 = 0.8
458         t8 = 0.6
459         t9 = 1.9
460   
461         d4 =  1.
462         d5 =  2.
463         d6 =  3.
464         d7 =  4.
465         d8 =  8.
466         d9 = 10.
467   
468         l4 = 0.
469         gamma4 = 0.
470       
471         l5 = 0.
472         gamma5 = 0.
473       
474         l6 = 1.
475         gamma6 = 1.
476
477         l7 = 2.
478         gamma7 = 1.
479
480         l8 = 2.
481         gamma8 = 1.
482
483         l9 = 2.
484         gamma9 = 1.
485
486c----------------------------------------------------------------------           
487         call viscoN2(T,visco)  !! v given in microPa.s
488     
489         Tc   = 126.192d0
490         rhoc = 11.1839  * 1000 * mn2   !!!from mol.dm-3 to kg.m-3
491
492         tau  = Tc / T
493         delta = rho/rhoc
494
495         k1 =  N1 * visco + N2 * tau**t2 + N3 * tau**t3  !!! mW m-1 K-1     
496c--------- residual thermal conductivity
497
498         k2 = N4 * tau**t4 * delta**d4 * exp(-gamma4*delta**l4)         
499     &  +     N5 * tau**t5 * delta**d5 * exp(-gamma5*delta**l5)
500     &  +     N6 * tau**t6 * delta**d6 * exp(-gamma6*delta**l6)
501     &  +     N7 * tau**t7 * delta**d7 * exp(-gamma7*delta**l7)
502     &  +     N8 * tau**t8 * delta**d8 * exp(-gamma8*delta**l8)
503     &  +     N9 * tau**t9 * delta**d9 * exp(-gamma9*delta**l9)
504
505         kthn2 = k1 + k2
506
507         RETURN
508
509         END
510
511
512c======================================================================
513
514         subroutine viscoN2(T,visco)
515
516c        Compute viscosity of N2 (Lemmon and Jacobsen, 2003)
517
518c======================================================================
519
520         implicit none
521
522       include "microphys.h"
523c       include "microphysCO2.h"
524c      arguments
525c     -----------
526         
527      REAL T
528
529c     outputs
530c     -----------
531
532      DOUBLE PRECISION visco
533
534
535c     local
536c    ------------
537
538      DOUBLE PRECISION a0,a1,a2,a3,a4
539      DOUBLE PRECISION Tstar,factor,sigma,M2
540      DOUBLE PRECISION RGCS
541     
542
543c---------------------------------------------------------------------- 
544
545 
546      factor = 98.94   ! (K)   
547 
548      sigma  = 0.3656  ! (nm)
549 
550      a0 =  0.431
551      a1 = -0.4623
552      a2 =  0.08406
553      a3 =  0.005341
554      a4 = -0.00331
555 
556      M2 = mn2 * 1.00e3   !!! to g.mol-1
557 
558      Tstar = T*1./factor
559
560      RGCS = exp( a0 + a1 * log(Tstar) + a2 * (log(Tstar))**2. +
561     &                a3 * (log(Tstar))**3. + a4 * (log(Tstar))**4. )
562 
563 
564      visco = 0.0266958 * sqrt(M2*T) / ( sigma**2. * RGCS )  !!! microPa.s
565
566
567      RETURN
568
569      END
570
571
572c======================================================================
573
574         subroutine KthCO2Scalab(kthco2,T,rho)
575
576c        Compute thermal cond of CO2 (Scalabrin et al. 2006)
577
578c======================================================================
579
580         implicit none
581
582
583
584c      arguments
585c     -----------
586
587      REAL T
588      DOUBLE PRECISION rho
589
590c     outputs
591c     -----------
592
593      DOUBLE PRECISION kthco2
594
595c     LOCAL
596c     -----------
597
598      DOUBLE PRECISION Tc,Pc,rhoc, Lambdac
599     
600      DOUBLE PRECISION Tr, rhor, k1, k2
601
602      DOUBLE PRECISION g1,g2,g3,g4,g5,g6,g7,g8,g9,g10
603      DOUBLE PRECISION h1,h2,h3,h4,h5,h6,h7,h8,h9,h10
604      DOUBLE PRECISION n1,n2,n3,n4,n5,n6,n7,n8,n9,n10
605
606      Tc   = 304.1282   !(K)
607      Pc   = 7.3773e6   !(MPa)
608      rhoc = 467.6      !(kg.m-3)
609      Lambdac = 4.81384 !(mW.m-1K-1)
610 
611      g1 = 0.
612      g2 = 0.
613      g3 = 1.5
614      g4 = 0.0
615      g5 = 1.0
616      g6 = 1.5
617      g7 = 1.5
618      g8 = 1.5
619      g9 = 3.5
620      g10 = 5.5
621
622      h1 = 1.
623      h2 = 5.
624      h3 = 1.
625      h4 = 1.
626      h5 = 2.
627      h6 = 0.
628      h7 = 5.0
629      h8 = 9.0
630      h9 = 0.
631      h10 = 0.
632
633      n1 = 7.69857587
634      n2 = 0.159885811
635      n3 = 1.56918621
636      n4 = -6.73400790
637      n5 = 16.3890156
638      n6 = 3.69415242
639      n7 = 22.3205514
640      n8 = 66.1420950
641      n9 = -0.171779133
642      n10 = 0.00433043347
643
644      Tr   = T/Tc
645      rhor = rho/rhoc
646
647      k1 = n1*Tr**(g1)*rhor**(h1) + n2*Tr**(g2)*rhor**(h2)
648     &     + n3*Tr**(g3)*rhor**(h3)
649
650      k2 = n4*Tr**(g4)*rhor**(h4) + n5*Tr**(g5)*rhor**(h5) 
651     &    + n6*Tr**(g6)*rhor**(h6) + n7*Tr**(g7)*rhor**(h7)
652     &    + n8*Tr**(g8)*rhor**(h8) + n9*Tr**(g9)*rhor**(h9)
653     &    + n10*Tr**(g10)*rhor**(h10)
654   
655      k2  = exp(-5.*rhor**(2.)) * k2
656               
657      kthco2 = (k1 + k2) *  Lambdac   ! mW
658
659      RETURN
660
661      END
Note: See TracBrowser for help on using the repository browser.