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

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

Added modifications for CO2 clouds scheme in physiq_mod.F and added several routines and variables for CO2 microphysics. October 2016 Joachim Audouard

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