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

Last change on this file since 1711 was 1685, checked in by jaudouard, 8 years ago

Further modifications on CO2 clouds scheme. Water ice clouds can now serve as CCN for CO2 clouds

  • Property svn:executable set to *
File size: 7.5 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
68cccccccccccccccccccccccccccccccccccccccccccccccccc
69ccccccccccc ESSAIS TN MTETA = F (T) cccccccccccccc
70c      if (temp .gt. 200) then
71c         mtetalocal = mtetalocal
72c      else if (temp .lt. 190) then
73c         mtetalocal = mtetalocal-0.05
74c      else
75c         mtetalocal = mtetalocal - (190-temp)*0.005
76c      endif
77c----------------exp law, see Trainer 2008, J. Phys. Chem. C 2009, 113, 2036\u20132040
78       !mtetalocal = max(mtetalocal - 6005*exp(-0.065*temp),0.1)
79       !mtetalocal = max(mtetalocal - 6005*exp(-0.068*temp),0.1)
80               !print*, mtetalocal, temp
81cccccccccccccccccccccccccccccccccccccccccccccccccc
82cccccccccccccccccccccccccccccccccccccccccccccccccc
83c      IF (firstcall) THEN
84c          print*, ' ' 
85c          print*, 'dear user, please keep in mind that'
86c          print*, 'contact parameter IS constant'
87          !print*, 'contact parameter IS NOT constant:'
88          !print*, 'max(mteta - 6005*exp(-0.065*temp),0.1)'
89          !print*, 'max(mteta - 6005*exp(-0.068*temp),0.1)'
90c          print*, ' ' 
91c         firstcall=.false.
92c     END IF
93cccccccccccccccccccccccccccccccccccccccccccccccccc
94cccccccccccccccccccccccccccccccccccccccccccccccccc
95   
96c      write(*,*) "IN nuc, SAT = ",sat
97c      write(*,*) "IN nuc, mtetalocal = ",mtetalocal
98
99
100      if (sat .gt. 1.) then    ! minimum condition to activate nucleation
101
102        nco2   = pco2 / kbz / temp
103        rstar  = 2. * sigco2 * vo2co2 / (kbz*temp*dlog(sat))
104        gstar  = 4. * pi * (rstar * rstar * rstar) / (3.*vo2co2)
105       
106       fshapeco2simple = (2.+mtetalocal)*(1.-mtetalocal)*(1.-mtetalocal)
107     &                   / 4.
108
109c       Loop over size bins
110        do 200 i=1,nbinco2_cld
111c            write(*,*) "IN NUCLEA, i, RAD_CLDCO2(i) = ",i, rad_cldco2(i),
112c     &          n_ccn(i)
113
114          if ( n_ccn(i) .lt. 1e-10 ) then
115c           no dust, no need to compute nucleation!
116            nucrate(i)=0.
117            goto 210
118          endif
119
120          if (rad_cldco2(i).gt.3000.*rstar) then
121            zefshapeco2 = fshapeco2simple
122          else
123            zefshapeco2 = fshapeco2(mtetalocal,rad_cldco2(i)/rstar)
124          endif
125
126          fistar = (4./3.*pi) * sigco2 * (rstar * rstar) *
127     &             zefshapeco2
128          deltaf = (2.*desorpco2-surfdifco2-fistar)/
129     &             (kbz*temp)
130          deltaf = min( max(deltaf, -100.d0), 100.d0)
131
132          if (deltaf.eq.-100.) then
133            nucrate(i) = 0.
134          else
135            nucrate(i)= dble(sqrt ( fistar /
136     &               (3.*pi*kbz*temp*(gstar*gstar)) )
137     &                  * kbz * temp * rstar
138     &                  * rstar * 4. * pi
139     &                  * ( nco2*rad_cldco2(i) )
140     &                  * ( nco2*rad_cldco2(i) )
141     &                  / ( zefshapeco2 * nusco2 * m0co2 )
142     &                  * dexp (deltaf))
143
144           
145          endif
146
147210     continue
148
149          if ( n_ccn_h2oice(i) .lt. 1e-10 ) then
150c           no dust, no need to compute nucleation!
151            nucrate_h2oice(i)=0.
152            goto 200
153          endif
154
155          if (rad_h2oice(i).gt.3000.*rstar) then
156            zefshapeco2 = (2.+mtetalocalh)*(1.-mtetalocalh)*
157     &                (1.-mtetalocalh) / 4.
158          else
159            zefshapeco2 = fshapeco2(mtetalocalh,rad_h2oice(i)/rstar) ! same m for dust/h2o ice
160          endif
161
162          fistar = (4./3.*pi) * sigco2 * (rstar * rstar) *
163     &             zefshapeco2
164          deltaf = (2.*desorpco2-surfdifco2-fistar)/
165     &             (kbz*temp)
166          deltaf = min( max(deltaf, -100.d0), 100.d0)
167
168          if (deltaf.eq.-100.) then
169            nucrate_h2oice(i) = 0.
170          else
171            nucrate_h2oice(i)= dble(sqrt ( fistar /
172     &               (3.*pi*kbz*temp*(gstar*gstar)) )
173     &                  * kbz * temp * rstar
174     &                  * rstar * 4. * pi
175     &                  * ( nco2*rad_h2oice(i) )
176     &                  * ( nco2*rad_h2oice(i) )
177     &                  / ( zefshapeco2 * nusco2 * m0co2 )
178     &                  * dexp (deltaf))
179          endif
180         
181         
182
183200     continue
184
185      else
186
187        do i=1,nbinco2_cld
188          nucrate(i) = 0.
189          nucrate_h2oice(i) = 0.
190        enddo
191
192      endif
193
194      return
195      end
196
197*********************************************************
198      double precision function fshapeco2(cost,rap)
199      implicit none
200*        function computing the f(m,x) factor           *
201* related to energy required to form a critical embryo  *
202*********************************************************
203
204      double precision cost,rap
205      double precision yeah
206
207          !! PHI
208          yeah = sqrt( 1. - 2.*cost*rap + rap*rap )
209          !! FSHAPECO2 = TERM A
210          fshapeco2 = (1.-cost*rap) / yeah
211          fshapeco2 = fshapeco2 * fshapeco2 * fshapeco2
212          fshapeco2 = 1. + fshapeco2
213          !! ... + TERM B
214          yeah = (rap-cost)/yeah
215          fshapeco2 = fshapeco2 +
216     & rap*rap*rap*(2.-3.*yeah+yeah*yeah*yeah)
217          !! ... + TERM C
218          fshapeco2 = fshapeco2 + 3. * cost * rap * rap * (yeah-1.)
219          !! FACTOR 1/2
220          fshapeco2 = 0.5*fshapeco2
221
222      return
223      end
Note: See TracBrowser for help on using the repository browser.