source: trunk/LMDZ.MARS/libf/phymars/nucleaco2.F90 @ 3739

Last change on this file since 3739 was 3739, checked in by emillour, 7 weeks ago

Mars PCM:
Code tidying: put routines in modules, remove useless "return" statements, etc.
EM

File size: 3.8 KB
Line 
1module nucleaco2_mod
2  implicit none
3  contains
4
5subroutine nucleaco2(pco2,temp,sat,n_ccn,nucrate,vo2co2, teta)
6
7use comcstfi_h, only: pi
8use microphys_h, only: nbinco2_cld, rad_cldco2, desorpco2, m0co2, kbz
9use microphys_h, only: nusco2, sigco2, surfdifco2
10
11implicit none
12!*                                                     *
13!*   This subroutine computes the nucleation rate      *
14!*   as given in Pruppacher & Klett (1978) in the      *
15!*   case of water ice forming on a solid substrate.   *
16!*     Definition refined by Keese (jgr,1989)          *
17!*   Authors: F. Montmessin                            *
18!*     Adapted for the LMD/GCM by J.-B. Madeleine      *
19!*     (October 2011)                                  *
20!*     Optimisation by A. Spiga (February 2012)        *
21!*   CO2 nucleation routine dev. by Constantino        *
22!*     Listowski and Joachim Audouard (2016-2017),     *
23!*     adapted from the water ice nucleation
24!* It computes two different nucleation rates : one
25!* on the dust CCN distribution and the other one on
26!* the water ice particles distribution
27!*******************************************************
28 ! nucrate = output
29 ! nucrate_h2o en sortie aussi :
30!nucleation sur dust et h2o separement ici
31
32!     Inputs
33DOUBLE PRECISION, intent(in) :: pco2,sat,vo2co2
34DOUBLE PRECISION, intent(in) :: n_ccn(nbinco2_cld)
35REAL, intent(in) :: temp !temperature
36REAL, intent(in) :: teta
37!     Output
38DOUBLE PRECISION, intent(out) :: nucrate(nbinco2_cld)
39
40!     Local variables
41DOUBLE PRECISION nco2
42DOUBLE PRECISION rstar    ! Radius of the critical germ (m)
43DOUBLE PRECISION gstar    ! # of molecules forming a critical embryo
44DOUBLE PRECISION fistar   ! Activation energy required to form a critical embryo (J)
45!DOUBLE PRECISION fshapeco2   ! function defined at the end of the file
46DOUBLE PRECISION deltaf
47double precision mtetalocal ! local teta in double precision
48double precision fshapeco2simple,zefshapeco2
49integer i
50!*************************************************
51
52mtetalocal = dble(teta)
53
54nco2   = pco2 / kbz / temp
55rstar  = 2. * sigco2 * vo2co2 / (kbz*temp*log(sat))
56gstar  = 4. * pi * (rstar * rstar * rstar) / (3.*vo2co2)
57
58fshapeco2simple = (2.+mtetalocal)*(1.-mtetalocal)*(1.-mtetalocal) / 4.
59
60!c       Loop over size bins
61do i=1,nbinco2_cld
62
63  if ( n_ccn(i) .lt. 1e-10 ) then
64    !c no dust, no need to compute nucleation!
65    nucrate(i)=0.
66  else
67    if (rad_cldco2(i).gt.3000.*rstar) then
68      zefshapeco2 = fshapeco2simple
69    else
70      zefshapeco2 = fshapeco2(mtetalocal, rad_cldco2(i)/rstar)
71    endif
72
73    fistar = (4./3.*pi) * sigco2 * (rstar * rstar) * zefshapeco2
74    deltaf = (2.*desorpco2-surfdifco2-fistar) / (kbz*temp)
75    deltaf = min( max(deltaf, -100.d0), 100.d0)
76
77    if (deltaf.eq.-100.) then
78        nucrate(i) = 0.
79    else
80        nucrate(i) = dble(sqrt ( fistar / (3.*pi*kbz*temp*(gstar*gstar)) ) * kbz * temp * rstar * rstar * 4. * pi &
81            * ( nco2*rad_cldco2(i) ) * ( nco2*rad_cldco2(i) ) / ( zefshapeco2 * nusco2 * m0co2 )  * exp (deltaf))
82    endif
83  endif ! if n_ccn(i) .lt. 1e-10
84
85enddo
86
87end subroutine nucleaco2
88
89
90!*********************************************************
91double precision function fshapeco2(cost,rap)
92implicit none
93!*        function computing the f(m,x) factor           *
94!* related to energy required to form a critical embryo  *
95!*********************************************************
96
97double precision, intent(in) :: cost,rap
98double precision yeah
99
100!! PHI
101yeah = sqrt( 1. - 2.*cost*rap + rap*rap )
102!! FSHAPECO2 = TERM A
103fshapeco2 = (1.-cost*rap) / yeah
104fshapeco2 = fshapeco2 * fshapeco2 * fshapeco2
105fshapeco2 = 1. + fshapeco2
106!! ... + TERM B
107yeah = (rap-cost)/yeah
108fshapeco2 = fshapeco2 + rap*rap*rap*(2.-3.*yeah+yeah*yeah*yeah)
109!! ... + TERM C
110fshapeco2 = fshapeco2 + 3. * cost * rap * rap * (yeah-1.)
111!! FACTOR 1/2
112fshapeco2 = 0.5*fshapeco2
113
114end function fshapeco2
115
116end module nucleaco2_mod
Note: See TracBrowser for help on using the repository browser.