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

Last change on this file since 1911 was 1818, checked in by jaudouard, 7 years ago

Final update from J.A. about the CO2 clouds scheme for the LMDZ.MARS GCM

  • Property svn:executable set to *
File size: 6.1 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*   CO2 nucleation routine dev. by Constantino        *
19*     Listowski and Joachim Audouard (2016-2017),     *
20*     adapted from the water ice nucleation     
21* It computes two different nucleation rates : one
22* on the dust CCN distribution and the other one on
23* the water ice particles distribution
24*******************************************************
25 ! nucrate = output
26 ! nucrate_h2o en sortie aussi :
27!nucleation sur dust et h2o separement ici
28#include "microphys.h"
29
30c     Inputs
31      DOUBLE PRECISION pco2,sat,vo2co2
32      DOUBLE PRECISION n_ccn(nbinco2_cld), n_ccn_h2oice(nbinco2_cld)
33      REAL temp !temperature
34c     Output
35      DOUBLE PRECISION nucrate(nbinco2_cld)
36      DOUBLE PRECISION nucrate_h2oice(nbinco2_cld) ! h2o as substrate
37      double precision rad_h2oice(nbinco2_cld) ! h2o ice grid (as substrate)
38
39c     Local variables
40      DOUBLE PRECISION nco2
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 fshapeco2   ! function defined at the end of the file
45      DOUBLE PRECISION deltaf     
46      double precision mtetalocal,mtetalocalh ! local mteta in double precision
47      double precision fshapeco2simple,zefshapeco2
48      integer i
49c     *************************************************
50
51      mtetalocal = dble(mtetaco2)  !! use mtetalocal for better performance
52      mtetalocalh=dble(mteta)
53
54
55      if (sat .gt. 1.) then    ! minimum condition to activate nucleation
56
57        nco2   = pco2 / kbz / temp
58        rstar  = 2. * sigco2 * vo2co2 / (kbz*temp*dlog(sat))
59        gstar  = 4. * pi * (rstar * rstar * rstar) / (3.*vo2co2)
60       
61       fshapeco2simple = (2.+mtetalocal)*(1.-mtetalocal)*(1.-mtetalocal)
62     &                   / 4.
63
64c       Loop over size bins
65        do 200 i=1,nbinco2_cld
66c            write(*,*) "IN NUCLEA, i, RAD_CLDCO2(i) = ",i, rad_cldco2(i),
67c     &          n_ccn(i)
68
69          if ( n_ccn(i) .lt. 1e-10 ) then
70c           no dust, no need to compute nucleation!
71            nucrate(i)=0.
72            goto 210
73          endif
74
75          if (rad_cldco2(i).gt.3000.*rstar) then
76            zefshapeco2 = fshapeco2simple
77          else
78            zefshapeco2 = fshapeco2(mtetalocal,rad_cldco2(i)/rstar)
79          endif
80
81          fistar = (4./3.*pi) * sigco2 * (rstar * rstar) *
82     &             zefshapeco2
83          deltaf = (2.*desorpco2-surfdifco2-fistar)/
84     &             (kbz*temp)
85          deltaf = min( max(deltaf, -100.d0), 100.d0)
86
87          if (deltaf.eq.-100.) then
88            nucrate(i) = 0.
89          else
90            nucrate(i)= dble(sqrt ( fistar /
91     &               (3.*pi*kbz*temp*(gstar*gstar)) )
92     &                  * kbz * temp * rstar
93     &                  * rstar * 4. * pi
94     &                  * ( nco2*rad_cldco2(i) )
95     &                  * ( nco2*rad_cldco2(i) )
96     &                  / ( zefshapeco2 * nusco2 * m0co2 )
97     &                  * dexp (deltaf))
98
99           
100          endif
101
102210     continue
103
104          if ( n_ccn_h2oice(i) .lt. 1e-10 ) then
105c           no dust, no need to compute nucleation!
106            nucrate_h2oice(i)=0.
107            goto 200
108          endif
109
110          if (rad_h2oice(i).gt.3000.*rstar) then
111            zefshapeco2 = (2.+mtetalocalh)*(1.-mtetalocalh)*
112     &                (1.-mtetalocalh) / 4.
113          else
114            zefshapeco2 = fshapeco2(mtetalocalh,rad_h2oice(i)/rstar) ! same m for dust/h2o ice
115          endif
116
117          fistar = (4./3.*pi) * sigco2 * (rstar * rstar) *
118     &             zefshapeco2
119          deltaf = (2.*desorpco2-surfdifco2-fistar)/
120     &             (kbz*temp)
121          deltaf = min( max(deltaf, -100.d0), 100.d0)
122
123          if (deltaf.eq.-100.) then
124            nucrate_h2oice(i) = 0.
125          else
126            nucrate_h2oice(i)= dble(sqrt ( fistar /
127     &               (3.*pi*kbz*temp*(gstar*gstar)) )
128     &                  * kbz * temp * rstar
129     &                  * rstar * 4. * pi
130     &                  * ( nco2*rad_h2oice(i) )
131     &                  * ( nco2*rad_h2oice(i) )
132     &                  / ( zefshapeco2 * nusco2 * m0co2 )
133     &                  * dexp (deltaf))
134          endif
135         
136         
137
138200     continue
139
140      else
141
142        do i=1,nbinco2_cld
143          nucrate(i) = 0.
144          nucrate_h2oice(i) = 0.
145        enddo
146
147      endif
148
149      return
150      end
151
152*********************************************************
153      double precision function fshapeco2(cost,rap)
154      implicit none
155*        function computing the f(m,x) factor           *
156* related to energy required to form a critical embryo  *
157*********************************************************
158
159      double precision cost,rap
160      double precision yeah
161
162          !! PHI
163          yeah = sqrt( 1. - 2.*cost*rap + rap*rap )
164          !! FSHAPECO2 = TERM A
165          fshapeco2 = (1.-cost*rap) / yeah
166          fshapeco2 = fshapeco2 * fshapeco2 * fshapeco2
167          fshapeco2 = 1. + fshapeco2
168          !! ... + TERM B
169          yeah = (rap-cost)/yeah
170          fshapeco2 = fshapeco2 +
171     & rap*rap*rap*(2.-3.*yeah+yeah*yeah*yeah)
172          !! ... + TERM C
173          fshapeco2 = fshapeco2 + 3. * cost * rap * rap * (yeah-1.)
174          !! FACTOR 1/2
175          fshapeco2 = 0.5*fshapeco2
176
177      return
178      end
Note: See TracBrowser for help on using the repository browser.