source: trunk/LMDZ.MARS/libf/phymars/nucleaCO2.F @ 1770

Last change on this file since 1770 was 1720, checked in by jaudouard, 8 years ago

Update on CO2 ice clouds scheme

  • Property svn:executable set to *
File size: 6.2 KB
Line 
1*******************************************************
2*                                                     *
3      subroutine nucleaCO2(pco2,temp,sat,n_ccn,nucrate,
4     &           n_ccn_h2oice,rad_h2oice,nucrate_h2oice,
5     &           vo2co2)
6      USE comcstfi_h
7
8      implicit none
9*                                                     *
10*   This subroutine computes the nucleation rate      *
11*   as given in Pruppacher & Klett (1978) in the      *
12*   case of water ice forming on a solid substrate.   *
13*     Definition refined by Keese (jgr,1989)          *
14*   Authors: F. Montmessin                            *
15*     Adapted for the LMD/GCM by J.-B. Madeleine      *
16*     (October 2011)                                  *
17*     Optimisation by A. Spiga (February 2012)        * 
18*******************************************************
19 ! nucrate = output
20      ! nucrate_h2o en sortie aussi :
21!nucleation sur dust et h2o separement ici
22!#include "tracer.h"
23#include "microphys.h"
24c#include "microphysCO2.h"
25
26c     Inputs
27      DOUBLE PRECISION pco2,sat,vo2co2
28      DOUBLE PRECISION n_ccn(nbinco2_cld), n_ccn_h2oice(nbinco2_cld)
29      REAL temp
30
31c     Output
32   !   DOUBLE PRECISION nucrate(nbinco2_cld)
33      DOUBLE PRECISION nucrate(nbinco2_cld)
34      DOUBLE PRECISION nucrate_h2oice(nbinco2_cld) ! h2o as substrate
35
36      double precision rad_h2oice(nbinco2_cld) ! h2o ice grid (as substrate)
37
38c     Local variables
39      DOUBLE PRECISION nco2
40c      DOUBLE PRECISION sigco2      ! Water-ice/air surface tension  (N.m)
41c      external sigco2
42      DOUBLE PRECISION rstar    ! Radius of the critical germ (m)
43      DOUBLE PRECISION gstar    ! # of molecules forming a critical embryo
44      DOUBLE PRECISION fistar   ! Activation energy required to form a critical embryo (J)
45!      DOUBLE PRECISION zeldov   ! Zeldovitch factor (no dim)
46      DOUBLE PRECISION fshapeco2   ! function defined at the end of the file
47      DOUBLE PRECISION deltaf
48
49c     Ratio rstar/radius of the nucleating dust particle
50c     double precision xratio
51     
52      double precision mtetalocal,mtetalocalh ! local mteta in double precision
53
54      double precision fshapeco2simple,zefshapeco2
55
56
57      integer i
58     
59      LOGICAL firstcall
60      DATA firstcall/.true./
61      SAVE firstcall
62
63c     *************************************************
64
65      mtetalocal = dble(mtetaco2)  !! use mtetalocal for better performance
66      mtetalocalh=dble(mteta)
67
68
69      if (sat .gt. 1.) then    ! minimum condition to activate nucleation
70
71        nco2   = pco2 / kbz / temp
72        rstar  = 2. * sigco2 * vo2co2 / (kbz*temp*dlog(sat))
73        gstar  = 4. * pi * (rstar * rstar * rstar) / (3.*vo2co2)
74       
75       fshapeco2simple = (2.+mtetalocal)*(1.-mtetalocal)*(1.-mtetalocal)
76     &                   / 4.
77
78c       Loop over size bins
79        do 200 i=1,nbinco2_cld
80c            write(*,*) "IN NUCLEA, i, RAD_CLDCO2(i) = ",i, rad_cldco2(i),
81c     &          n_ccn(i)
82
83          if ( n_ccn(i) .lt. 1e-10 ) then
84c           no dust, no need to compute nucleation!
85            nucrate(i)=0.
86            goto 210
87          endif
88
89          if (rad_cldco2(i).gt.3000.*rstar) then
90            zefshapeco2 = fshapeco2simple
91          else
92            zefshapeco2 = fshapeco2(mtetalocal,rad_cldco2(i)/rstar)
93          endif
94
95          fistar = (4./3.*pi) * sigco2 * (rstar * rstar) *
96     &             zefshapeco2
97          deltaf = (2.*desorpco2-surfdifco2-fistar)/
98     &             (kbz*temp)
99          deltaf = min( max(deltaf, -100.d0), 100.d0)
100
101          if (deltaf.eq.-100.) then
102            nucrate(i) = 0.
103          else
104            nucrate(i)= dble(sqrt ( fistar /
105     &               (3.*pi*kbz*temp*(gstar*gstar)) )
106     &                  * kbz * temp * rstar
107     &                  * rstar * 4. * pi
108     &                  * ( nco2*rad_cldco2(i) )
109     &                  * ( nco2*rad_cldco2(i) )
110     &                  / ( zefshapeco2 * nusco2 * m0co2 )
111     &                  * dexp (deltaf))
112
113           
114          endif
115
116210     continue
117
118          if ( n_ccn_h2oice(i) .lt. 1e-10 ) then
119c           no dust, no need to compute nucleation!
120            nucrate_h2oice(i)=0.
121            goto 200
122          endif
123
124          if (rad_h2oice(i).gt.3000.*rstar) then
125            zefshapeco2 = (2.+mtetalocalh)*(1.-mtetalocalh)*
126     &                (1.-mtetalocalh) / 4.
127          else
128            zefshapeco2 = fshapeco2(mtetalocalh,rad_h2oice(i)/rstar) ! same m for dust/h2o ice
129          endif
130
131          fistar = (4./3.*pi) * sigco2 * (rstar * rstar) *
132     &             zefshapeco2
133          deltaf = (2.*desorpco2-surfdifco2-fistar)/
134     &             (kbz*temp)
135          deltaf = min( max(deltaf, -100.d0), 100.d0)
136
137          if (deltaf.eq.-100.) then
138            nucrate_h2oice(i) = 0.
139          else
140            nucrate_h2oice(i)= dble(sqrt ( fistar /
141     &               (3.*pi*kbz*temp*(gstar*gstar)) )
142     &                  * kbz * temp * rstar
143     &                  * rstar * 4. * pi
144     &                  * ( nco2*rad_h2oice(i) )
145     &                  * ( nco2*rad_h2oice(i) )
146     &                  / ( zefshapeco2 * nusco2 * m0co2 )
147     &                  * dexp (deltaf))
148          endif
149         
150         
151
152200     continue
153
154      else
155
156        do i=1,nbinco2_cld
157          nucrate(i) = 0.
158          nucrate_h2oice(i) = 0.
159        enddo
160
161      endif
162
163      return
164      end
165
166*********************************************************
167      double precision function fshapeco2(cost,rap)
168      implicit none
169*        function computing the f(m,x) factor           *
170* related to energy required to form a critical embryo  *
171*********************************************************
172
173      double precision cost,rap
174      double precision yeah
175
176          !! PHI
177          yeah = sqrt( 1. - 2.*cost*rap + rap*rap )
178          !! FSHAPECO2 = TERM A
179          fshapeco2 = (1.-cost*rap) / yeah
180          fshapeco2 = fshapeco2 * fshapeco2 * fshapeco2
181          fshapeco2 = 1. + fshapeco2
182          !! ... + TERM B
183          yeah = (rap-cost)/yeah
184          fshapeco2 = fshapeco2 +
185     & rap*rap*rap*(2.-3.*yeah+yeah*yeah*yeah)
186          !! ... + TERM C
187          fshapeco2 = fshapeco2 + 3. * cost * rap * rap * (yeah-1.)
188          !! FACTOR 1/2
189          fshapeco2 = 0.5*fshapeco2
190
191      return
192      end
Note: See TracBrowser for help on using the repository browser.