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

Last change on this file since 3026 was 3008, checked in by emillour, 18 months ago

Mars PCM:
Some code cleanup around microphysics. Turn microphys.h into module
microphys_h.F90, and while at it also turn nuclea.F, growthrate.F90 and
massflowrateco2.F90 into modules.
EM

  • Property svn:executable set to *
File size: 3.8 KB
RevLine 
[2562]1module nucleaco2_mod
2  implicit none
3  contains
[3008]4
[2562]5subroutine nucleaco2(pco2,temp,sat,n_ccn,nucrate,vo2co2, teta)
6
[3008]7use comcstfi_h, only: pi
8use microphys_h, only: nbinco2_cld, rad_cldco2, desorpco2, m0co2, kbz
9use microphys_h, only: nusco2, sigco2, surfdifco2
10
[2562]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
32include "callkeys.h"
33
34!     Inputs
35DOUBLE PRECISION, intent(in) :: pco2,sat,vo2co2
36DOUBLE PRECISION, intent(in) :: n_ccn(nbinco2_cld)
37REAL, intent(in) :: temp !temperature
38REAL, intent(in) :: teta
39!     Output
40DOUBLE PRECISION, intent(out) :: nucrate(nbinco2_cld)
41
42!     Local variables
43DOUBLE PRECISION nco2
44DOUBLE PRECISION rstar    ! Radius of the critical germ (m)
45DOUBLE PRECISION gstar    ! # of molecules forming a critical embryo
46DOUBLE PRECISION fistar   ! Activation energy required to form a critical embryo (J)
47!DOUBLE PRECISION fshapeco2   ! function defined at the end of the file
48DOUBLE PRECISION deltaf
49double precision mtetalocal ! local teta in double precision
50double precision fshapeco2simple,zefshapeco2
51integer i
52!*************************************************
53
54mtetalocal = dble(teta)
55
56nco2   = pco2 / kbz / temp
57rstar  = 2. * sigco2 * vo2co2 / (kbz*temp*log(sat))
58gstar  = 4. * pi * (rstar * rstar * rstar) / (3.*vo2co2)
59
60fshapeco2simple = (2.+mtetalocal)*(1.-mtetalocal)*(1.-mtetalocal) / 4.
61
62!c       Loop over size bins
63do i=1,nbinco2_cld
64
65  if ( n_ccn(i) .lt. 1e-10 ) then
66    !c no dust, no need to compute nucleation!
67    nucrate(i)=0.
68  else
69    if (rad_cldco2(i).gt.3000.*rstar) then
70      zefshapeco2 = fshapeco2simple
71    else
72      zefshapeco2 = fshapeco2(mtetalocal, rad_cldco2(i)/rstar)
73    endif
74
75    fistar = (4./3.*pi) * sigco2 * (rstar * rstar) * zefshapeco2
76    deltaf = (2.*desorpco2-surfdifco2-fistar) / (kbz*temp)
77    deltaf = min( max(deltaf, -100.d0), 100.d0)
78
79    if (deltaf.eq.-100.) then
80        nucrate(i) = 0.
81    else
82        nucrate(i) = dble(sqrt ( fistar / (3.*pi*kbz*temp*(gstar*gstar)) ) * kbz * temp * rstar * rstar * 4. * pi &
83            * ( nco2*rad_cldco2(i) ) * ( nco2*rad_cldco2(i) ) / ( zefshapeco2 * nusco2 * m0co2 )  * exp (deltaf))
84    endif
85  endif ! if n_ccn(i) .lt. 1e-10
86
87enddo
88
89end subroutine nucleaco2
90
91
92!*********************************************************
93double precision function fshapeco2(cost,rap)
94implicit none
95!*        function computing the f(m,x) factor           *
96!* related to energy required to form a critical embryo  *
97!*********************************************************
98
[3008]99double precision, intent(in) :: cost,rap
[2562]100double precision yeah
101
102!! PHI
103yeah = sqrt( 1. - 2.*cost*rap + rap*rap )
104!! FSHAPECO2 = TERM A
105fshapeco2 = (1.-cost*rap) / yeah
106fshapeco2 = fshapeco2 * fshapeco2 * fshapeco2
107fshapeco2 = 1. + fshapeco2
108!! ... + TERM B
109yeah = (rap-cost)/yeah
110fshapeco2 = fshapeco2 + rap*rap*rap*(2.-3.*yeah+yeah*yeah*yeah)
111!! ... + TERM C
112fshapeco2 = fshapeco2 + 3. * cost * rap * rap * (yeah-1.)
113!! FACTOR 1/2
114fshapeco2 = 0.5*fshapeco2
115
116end function fshapeco2
117
118end module nucleaco2_mod
Note: See TracBrowser for help on using the repository browser.