source: lmdz_wrf/WRFV3/lmdz/clouds_gno.F90 @ 1

Last change on this file since 1 was 1, checked in by lfita, 10 years ago
  • -- --- Opening of the WRF+LMDZ coupling repository --- -- -

WRF: version v3.3
LMDZ: version v1818

More details in:

File size: 8.6 KB
Line 
1!
2! $Header$
3!
4!C
5!C================================================================================
6!C
7      SUBROUTINE CLOUDS_GNO(klon,ND,R,RS,QSUB,PTCONV,RATQSC,CLDF)
8      IMPLICIT NONE
9!C     
10!C--------------------------------------------------------------------------------
11!C
12!C Inputs:
13!C
14!C  ND----------: Number of vertical levels
15!C  R--------ND-: Domain-averaged mixing ratio of total water
16!C  RS-------ND-: Mean saturation humidity mixing ratio within the gridbox
17!C  QSUB-----ND-: Mixing ratio of condensed water within clouds associated
18!C                with SUBGRID-SCALE condensation processes (here, it is
19!C                predicted by the convection scheme)
20!C Outputs:
21!C
22!C  PTCONV-----ND-: Point convectif = TRUE
23!C  RATQSC-----ND-: Largeur normalisee de la distribution
24!C  CLDF-----ND-: Fraction nuageuse
25!C
26!C--------------------------------------------------------------------------------
27
28
29      INTEGER klon,ND
30      REAL R(klon,ND),  RS(klon,ND), QSUB(klon,ND)
31      LOGICAL PTCONV(klon,ND)
32      REAL RATQSC(klon,ND)
33      REAL CLDF(klon,ND)
34
35!c -- parameters controlling the iteration:
36!c --    nmax    : maximum nb of iterations (hopefully never reached)
37!c --    epsilon : accuracy of the numerical resolution
38!c --    vmax    : v-value above which we use an asymptotic expression for ERF(v)
39
40      INTEGER nmax
41      PARAMETER ( nmax = 10)
42      REAL epsilon, vmax0, vmax(klon)
43      PARAMETER ( epsilon = 0.02, vmax0 = 2.0 )
44
45      REAL min_mu, min_Q
46      PARAMETER ( min_mu =  1.e-12, min_Q=1.e-12 )
47      INTEGER i,K, n, m
48      REAL mu(klon), qsat, delta(klon), beta(klon)
49      real zu2,zv2
50      REAL xx(klon), aux(klon), coeff, block
51      REAL  dist, fprime, det
52      REAL pi, u, v, erfcu, erfcv
53      REAL  xx1, xx2
54      real erf,hsqrtlog_2,v2
55      real sqrtpi,sqrt2,zx1,zx2,exdel
56!c lconv = true si le calcul a converge (entre autre si qsub < min_q)
57       LOGICAL lconv(klon)
58
59!cdir arraycomb
60      cldf  (1:klon,1:ND)=0.0        ! cym
61      ratqsc(1:klon,1:ND)=0.0
62      ptconv(1:klon,1:ND)=.false.
63!cdir end arraycomb
64     
65      pi = ACOS(-1.)
66      sqrtpi=sqrt(pi)
67      sqrt2=sqrt(2.)
68      hsqrtlog_2=0.5*SQRT(log(2.))
69
70      DO 500 K = 1, ND
71
72                                    do i=1,klon ! vector
73      mu(i) = R(i,K)
74      mu(i) = MAX(mu(i),min_mu)
75      qsat = RS(i,K)
76      qsat = MAX(qsat,min_mu)
77      delta(i) = log(mu(i)/qsat)
78!c                                   enddo ! vector
79
80!C
81!C ***          There is no subgrid-scale condensation;        ***
82!C ***   the scheme becomes equivalent to an "all-or-nothing"  ***
83!C ***             large-scale condensation scheme.            ***
84!C
85
86!C
87!C ***     Some condensation is produced at the subgrid-scale       ***
88!C ***                                                              ***
89!C ***       PDF = generalized log-normal distribution (GNO)        ***
90!C ***   (k<0 because a lower bound is considered for the PDF)      ***
91!C ***                                                              ***
92!C ***  -> Determine x (the parameter k of the GNO PDF) such        ***
93!C ***  that the contribution of subgrid-scale processes to         ***
94!C ***  the in-cloud water content is equal to QSUB(K)              ***
95!C ***  (equations (13), (14), (15) + Appendix B of the paper)      ***
96!C ***                                                              ***
97!C ***    Here, an iterative method is used for this purpose        ***
98!C ***    (other numerical methods might be more efficient)         ***
99!C ***                                                              ***
100!C ***          NB: the "error function" is called ERF              ***
101!C ***                 (ERF in double precision)                   ***
102!C
103
104!c  On commence par eliminer les cas pour lesquels on n'a pas
105!c  suffisamment d'eau nuageuse.
106
107!c                                   do i=1,klon ! vector
108
109      IF ( QSUB(i,K) .lt. min_Q ) THEN
110        ptconv(i,k)=.false.
111        ratqsc(i,k)=0.
112        lconv(i)  = .true.
113
114!c   Rien on a deja initialise
115
116      ELSE
117
118        lconv(i)  = .FALSE.
119        vmax(i) = vmax0
120
121        beta(i) = QSUB(i,K)/mu(i) + EXP( -MIN(0.0,delta(i)) )
122
123!c --  roots of equation v > vmax:
124
125        det = delta(i) + vmax(i)*vmax(i)
126        if (det.LE.0.0) vmax(i) = vmax0 + 1.0
127        det = delta(i) + vmax(i)*vmax(i)
128
129        if (det.LE.0.) then
130          xx(i) = -0.0001
131        else
132         zx1=-sqrt2*vmax(i)
133         zx2=SQRT(1.0+delta(i)/(vmax(i)*vmax(i)))
134         xx1=zx1*(1.0-zx2)
135         xx2=zx1*(1.0+zx2)
136         xx(i) = 1.01 * xx1
137         if ( xx1 .GE. 0.0 ) xx(i) = 0.5*xx2
138        endif
139        if (delta(i).LT.0.) xx(i) = -hsqrtlog_2
140
141      ENDIF
142
143                                    enddo       ! vector
144
145!c----------------------------------------------------------------------
146!c   Debut des nmax iterations pour trouver la solution.
147!c----------------------------------------------------------------------
148
149      DO n = 1, nmax
150
151                                    do i=1,klon ! vector
152        if (.not.lconv(i)) then
153
154          u = delta(i)/(xx(i)*sqrt2) + xx(i)/(2.*sqrt2)
155          v = delta(i)/(xx(i)*sqrt2) - xx(i)/(2.*sqrt2)
156          v2 = v*v
157
158          IF ( v .GT. vmax(i) ) THEN
159
160            IF (     ABS(u)  .GT. vmax(i)                                            &
161       &          .AND.  delta(i) .LT. 0. ) THEN
162
163!c -- use asymptotic expression of erf for u and v large:
164!c ( -> analytic solution for xx )
165             exdel=beta(i)*EXP(delta(i))
166             aux(i) = 2.0*delta(i)*(1.-exdel)                                        &
167       &                       /(1.+exdel)
168             if (aux(i).lt.0.) then
169!c                print*,'AUX(',i,',',k,')<0',aux(i),delta(i),beta(i)
170                aux(i)=0.
171             endif
172             xx(i) = -SQRT(aux(i))
173             block = EXP(-v*v) / v / sqrtpi
174             dist = 0.0
175             fprime = 1.0
176
177            ELSE
178
179!c -- erfv -> 1.0, use an asymptotic expression of erfv for v large:
180
181             erfcu = 1.0-ERF(u)
182!c  !!! ATTENTION : rajout d'un seuil pour l'exponentiel
183             aux(i) = sqrtpi*erfcu*EXP(min(v2,100.))
184             coeff = 1.0 - 0.5/(v2) + 0.75/(v2*v2)
185             block = coeff * EXP(-v2) / v / sqrtpi
186             dist = v * aux(i) / coeff - beta(i)
187             fprime = 2.0 / xx(i) * (v2)                                             &
188       &           * ( EXP(-delta(i)) - u * aux(i) / coeff )                         &
189       &           / coeff
190           
191            ENDIF ! ABS(u)
192
193          ELSE
194
195!c -- general case:
196
197           erfcu = 1.0-ERF(u)
198           erfcv = 1.0-ERF(v)
199           block = erfcv
200           dist = erfcu / erfcv - beta(i)
201           zu2=u*u
202           zv2=v2
203           if(zu2.gt.20..or. zv2.gt.20.) then
204!c              print*,'ATTENTION !!! xx(',i,') =', xx(i)
205!c           print*,'ATTENTION !!! klon,ND,R,RS,QSUB,PTCONV,RATQSC,CLDF',
206!c     .klon,ND,R(i,k),RS(i,k),QSUB(i,k),PTCONV(i,k),RATQSC(i,k),
207!c     .CLDF(i,k)
208!c              print*,'ATTENTION !!! zu2 zv2 =',zu2(i),zv2(i)
209              zu2=20.
210              zv2=20.
211             fprime = 0.
212           else
213             fprime = 2. /sqrtpi /xx(i) /(erfcv*erfcv)                               &
214       &           * (   erfcv*v*EXP(-zu2)                                           &
215       &               - erfcu*u*EXP(-zv2) )
216           endif
217          ENDIF ! x
218
219!c -- test numerical convergence:
220
221!          if (beta(i).lt.1.e-10) then
222!              print*,'avant test ',i,k,lconv(i),u(i),v(i),beta(i)
223!              stop
224!          endif
225          if (abs(fprime).lt.1.e-11) then
226!              print*,'avant test fprime<.e-11 '
227!     s        ,i,k,lconv(i),u(i),v(i),beta(i),fprime(i)
228!              print*,'klon,ND,R,RS,QSUB',
229!     s        klon,ND,R(i,k),rs(i,k),qsub(i,k)
230              fprime=sign(1.e-11,fprime)
231          endif
232
233
234          if ( ABS(dist/beta(i)) .LT. epsilon ) then
235!c           print*,'v-u **2',(v(i)-u(i))**2
236!c           print*,'exp v-u **2',exp((v(i)-u(i))**2)
237            ptconv(i,K) = .TRUE.
238            lconv(i)=.true.
239!c  borne pour l'exponentielle
240            ratqsc(i,k)=min(2.*(v-u)*(v-u),20.)
241            ratqsc(i,k)=sqrt(exp(ratqsc(i,k))-1.)
242            CLDF(i,K) = 0.5 * block
243          else
244            xx(i) = xx(i) - dist/fprime
245          endif
246!c         print*,'apres test ',i,k,lconv(i)
247
248        endif ! lconv
249                                    enddo       ! vector
250
251!c----------------------------------------------------------------------
252!c   Fin des nmax iterations pour trouver la solution.
253        ENDDO ! n
254!c----------------------------------------------------------------------
255
256500   CONTINUE  ! K
257
258       RETURN
259       END SUBROUTINE CLOUDS_GNO
260
261 
262
Note: See TracBrowser for help on using the repository browser.