source: codes/icosagcm/branches/SATURN_DYNAMICO/LMDZ.COMMON/libf/phystd/rain.F90

Last change on this file was 313, checked in by ymipsl, 10 years ago
  • implement splitting of XIOS file for lmdz physics
  • Termination is done properly in parallel by calling MPI_ABORT instead of abort or stop

YM

File size: 14.1 KB
Line 
1subroutine rain(ngrid,nlayer,nq,ptimestep,pplev,pplay,t,pdt,pq,pdq,d_t,dqrain,dqsrain,dqssnow,rneb)
2
3
4! to use  'getin'
5!  use ioipsl_getincom
6  use ioipsl_getincom_p
7  use watercommon_h, only: T_h2O_ice_liq,T_h2O_ice_clouds, RLVTT, RCPD, RCPV, RV, RVTMP2,Psat_water,Tsat_water,rhowater
8  use radii_mod, only: h2o_cloudrad
9  USE tracer_h, only: igcm_h2o_vap, igcm_h2o_ice
10  implicit none
11
12!==================================================================
13!     
14!     Purpose
15!     -------
16!     Calculates H2O precipitation using simplified microphysics.
17!     
18!     Authors
19!     -------
20!     Adapted from the LMDTERRE code by R. Wordsworth (2009)
21!     Added rain vaporization in case of T>Tsat
22!     Original author Z. X. Li (1993)
23!     
24!==================================================================
25
26!  include "dimensions.h"
27!  include "dimphys.h"
28  include "comcstfi.h"
29  include "callkeys.h"
30
31!     Arguments
32      integer,intent(in) :: ngrid ! number of atmospheric columns
33      integer,intent(in) :: nlayer ! number of atmospheric layers
34      integer,intent(in) :: nq ! number of tracers
35      real,intent(in) :: ptimestep    ! time interval
36      real,intent(in) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa)
37      real,intent(in) :: pplay(ngrid,nlayer)   ! mid-layer pressure (Pa)
38      real,intent(in) :: t(ngrid,nlayer) ! input temperature (K)
39      real,intent(in) :: pdt(ngrid,nlayer) ! input tendency on temperature (K/s)     
40      real,intent(in) :: pq(ngrid,nlayer,nq)  ! tracers (kg/kg)
41      real,intent(in) :: pdq(ngrid,nlayer,nq) ! input tendency on tracers
42      real,intent(out) :: d_t(ngrid,nlayer) ! temperature tendency (K/s)
43      real,intent(out) :: dqrain(ngrid,nlayer,nq) ! tendency of H2O precipitation (kg/kg.s-1)
44      real,intent(out) :: dqsrain(ngrid)  ! rain flux at the surface (kg.m-2.s-1)
45      real,intent(out) :: dqssnow(ngrid)  ! snow flux at the surface (kg.m-2.s-1)
46      real,intent(in) :: rneb(ngrid,nlayer) ! cloud fraction
47
48      REAL zt(ngrid,nlayer)         ! working temperature (K)
49      REAL ql(ngrid,nlayer)         ! liquid water (Kg/Kg)
50      REAL q(ngrid,nlayer)          ! specific humidity (Kg/Kg)
51      REAL d_q(ngrid,nlayer)        ! water vapor increment
52      REAL d_ql(ngrid,nlayer)       ! liquid water / ice increment
53
54!     Subroutine options
55      REAL,PARAMETER :: seuil_neb=0.001  ! Nebulosity threshold
56
57      INTEGER,save :: precip_scheme      ! id number for precipitaion scheme
58!     for simple scheme  (precip_scheme=1)
59      REAL,SAVE :: rainthreshold                ! Precipitation threshold in simple scheme
60!     for sundquist scheme  (precip_scheme=2-3)
61      REAL,SAVE :: cloud_sat                    ! Precipitation threshold in non simple scheme
62      REAL,SAVE :: precip_timescale             ! Precipitation timescale
63!     for Boucher scheme  (precip_scheme=4)
64      REAL,SAVE :: Cboucher ! Precipitation constant in Boucher 95 scheme
65      REAL,PARAMETER :: Kboucher=1.19E8
66      REAL,SAVE :: c1
67!$OMP THREADPRIVATE(precip_scheme,rainthreshold,cloud_sat,precip_timescale,Cboucher,c1)
68
69      INTEGER,PARAMETER :: ninter=5
70
71      logical,save :: evap_prec ! Does the rain evaporate?
72!$OMP THREADPRIVATE(evap_prec)
73
74!     for simple scheme
75      real,parameter :: t_crit=218.0
76      real lconvert
77
78!     Local variables
79      INTEGER i, k, n
80      REAL zqs(ngrid,nlayer),Tsat(ngrid,nlayer), zdelta, zcor
81      REAL zrfl(ngrid), zrfln(ngrid), zqev, zqevt
82
83      REAL zoliq(ngrid)
84      REAL zdz(ngrid),zrho(ngrid),ztot(ngrid), zrhol(ngrid)
85      REAL zchau(ngrid),zfroi(ngrid),zfrac(ngrid),zneb(ngrid)
86
87      real reffh2oliq(ngrid,nlayer),reffh2oice(ngrid,nlayer)
88 
89      real ttemp, ptemp, psat_tmp
90      real tnext(ngrid,nlayer)
91
92      real l2c(ngrid,nlayer)
93      real dWtot
94
95
96!     Indices of water vapour and water ice tracers
97      INTEGER, SAVE :: i_vap=0  ! water vapour
98      INTEGER, SAVE :: i_ice=0  ! water ice
99!$OMP THREADPRIVATE(i_vap,i_ice)
100
101      LOGICAL,SAVE :: firstcall=.true.
102!$OMP THREADPRIVATE(firstcall)
103
104!     Online functions
105      REAL fallv, fall2v, zzz ! falling speed of ice crystals
106      fallv (zzz) = 3.29 * ((zzz)**0.16)
107      fall2v (zzz) =10.6 * ((zzz)**0.31)  !for use with radii
108
109
110      IF (firstcall) THEN
111
112         i_vap=igcm_h2o_vap
113         i_ice=igcm_h2o_ice
114       
115         write(*,*) "rain: i_ice=",i_ice
116         write(*,*) "      i_vap=",i_vap
117
118         PRINT*, 'in rain.F, ninter=', ninter
119         PRINT*, 'in rain.F, evap_prec=', evap_prec
120
121         write(*,*) "Precipitation scheme to use?"
122         precip_scheme=1 ! default value
123         call getin_p("precip_scheme",precip_scheme)
124         write(*,*) " precip_scheme = ",precip_scheme
125
126         if (precip_scheme.eq.1) then
127            write(*,*) "rainthreshold in simple scheme?"
128            rainthreshold=0. ! default value
129            call getin_p("rainthreshold",rainthreshold)
130            write(*,*) " rainthreshold = ",rainthreshold
131
132         else if (precip_scheme.eq.2.or.precip_scheme.eq.3) then
133            write(*,*) "cloud water saturation level in non simple scheme?"
134            cloud_sat=2.6e-4   ! default value
135            call getin_p("cloud_sat",cloud_sat)
136            write(*,*) " cloud_sat = ",cloud_sat
137            write(*,*) "precipitation timescale in non simple scheme?"
138            precip_timescale=3600.  ! default value
139            call getin_p("precip_timescale",precip_timescale)
140            write(*,*) " precip_timescale = ",precip_timescale
141
142         else if (precip_scheme.eq.4) then
143            write(*,*) "multiplicative constant in Boucher 95 precip scheme"
144            Cboucher=1.   ! default value
145            call getin_p("Cboucher",Cboucher)
146            write(*,*) " Cboucher = ",Cboucher 
147            c1=1.00*1.097/rhowater*Cboucher*Kboucher 
148
149         endif
150
151         write(*,*) "re-evaporate precipitations?"
152         evap_prec=.true. ! default value
153         call getin_p("evap_prec",evap_prec)
154         write(*,*) " evap_prec = ",evap_prec
155
156         firstcall = .false.
157      ENDIF ! of IF (firstcall)
158
159!     GCM -----> subroutine variables
160      DO k = 1, nlayer
161      DO i = 1, ngrid
162
163         zt(i,k)   = t(i,k)+pdt(i,k)*ptimestep ! a big fat bug was here
164         q(i,k)    = pq(i,k,i_vap)+pdq(i,k,i_vap)*ptimestep
165         ql(i,k)   = pq(i,k,i_ice)+pdq(i,k,i_ice)*ptimestep
166
167         !q(i,k)    = pq(i,k,i_vap)!+pdq(i,k,i_vap)
168         !ql(i,k)   = pq(i,k,i_ice)!+pdq(i,k,i_ice)
169
170         if(q(i,k).lt.0.)then ! if this is not done, we don't conserve water
171            q(i,k)=0.
172         endif
173         if(ql(i,k).lt.0.)then
174            ql(i,k)=0.
175         endif
176
177      ENDDO
178      ENDDO
179
180!     Initialise the outputs
181      d_t(1:ngrid,1:nlayer) = 0.0
182      d_q(1:ngrid,1:nlayer) = 0.0
183      d_ql(1:ngrid,1:nlayer) = 0.0
184      zrfl(1:ngrid) = 0.0
185      zrfln(1:ngrid) = 0.0
186
187      ! calculate saturation mixing ratio
188      DO k = 1, nlayer
189         DO i = 1, ngrid
190            ttemp = zt(i,k) 
191            ptemp = pplay(i,k)
192!            call watersat(ttemp,ptemp,zqs(i,k))
193            call Psat_water(ttemp,ptemp,psat_tmp,zqs(i,k))
194            call Tsat_water(ptemp,Tsat(i,k))
195         ENDDO
196      ENDDO
197
198      ! get column / layer conversion factor
199      DO k = 1, nlayer
200         DO i = 1, ngrid
201            l2c(i,k)=(pplev(i,k)-pplev(i,k+1))/g
202         ENDDO
203      ENDDO
204
205      ! Vertical loop (from top to bottom)
206      ! We carry the rain with us and calculate that added by warm/cold precipitation
207      ! processes and that subtracted by evaporation at each level.
208      DO k = nlayer, 1, -1
209
210         IF (evap_prec) THEN ! note no rneb dependence!
211            DO i = 1, ngrid
212               IF (zrfl(i) .GT.0.) THEN
213
214                  if(zt(i,k).gt.Tsat(i,k))then
215!!                   treat the case where all liquid water should boil
216                     zqev=MIN((zt(i,k)-Tsat(i,k))*RCPD*l2c(i,k)/RLVTT/ptimestep,zrfl(i))
217                     zrfl(i)=MAX(zrfl(i)-zqev,0.)
218                     d_q(i,k)=zqev/l2c(i,k)*ptimestep
219                     d_t(i,k) = - d_q(i,k) * RLVTT/RCPD
220                  else
221                     zqev = MAX (0.0, (zqs(i,k)-q(i,k)))*l2c(i,k)/ptimestep !there was a bug here
222                     zqevt= 2.0e-5*(1.0-q(i,k)/zqs(i,k))    & !default was 2.e-5
223                        *sqrt(zrfl(i))*l2c(i,k)/pplay(i,k)*zt(i,k)*R ! BC modif here
224                     zqevt = MAX (zqevt, 0.0)
225                     zqev  = MIN (zqev, zqevt)
226                     zqev  = MAX (zqev, 0.0)
227                     zrfln(i)= zrfl(i) - zqev
228                     zrfln(i)= max(zrfln(i),0.0)
229
230                     d_q(i,k) = - (zrfln(i)-zrfl(i))/l2c(i,k)*ptimestep
231                     !d_t(i,k) = d_q(i,k) * RLVTT/RCPD!/(1.0+RVTMP2*q(i,k)) ! double BC modif here
232                     d_t(i,k) = - d_q(i,k) * RLVTT/RCPD ! was bugged!
233                     zrfl(i)  = zrfln(i)
234                  end if
235                     
236
237               ENDIF ! of IF (zrfl(i) .GT.0.)
238            ENDDO
239         ENDIF ! of IF (evap_prec)
240
241         zoliq(1:ngrid) = 0.0
242
243
244         if(precip_scheme.eq.1)then
245
246            DO i = 1, ngrid
247               ttemp = zt(i,k)
248               IF (ttemp .ge. T_h2O_ice_liq) THEN
249                  lconvert=rainthreshold
250               ELSEIF (ttemp .gt. t_crit) THEN
251                  lconvert=rainthreshold*(1.- t_crit/ttemp)
252                  lconvert=MAX(0.0,lconvert)             
253               ELSE
254                  lconvert=0.
255               ENDIF
256
257
258               IF (ql(i,k).gt.1.e-9) then
259                  zneb(i)  = MAX(rneb(i,k), seuil_neb)
260                  IF ((ql(i,k)/zneb(i)).gt.lconvert)THEN ! precipitate!
261                     d_ql(i,k) = -MAX((ql(i,k)-lconvert*zneb(i)),0.0)
262                     zrfl(i)   = zrfl(i) - d_ql(i,k)*l2c(i,k)/ptimestep
263                  ENDIF
264               ENDIF
265            ENDDO
266
267         elseif (precip_scheme.ge.2) then
268         
269           DO i = 1, ngrid
270               IF (rneb(i,k).GT.0.0) THEN
271                  zoliq(i) = ql(i,k)
272                  zrho(i)  = pplay(i,k) / ( zt(i,k) * R )
273                  zdz(i)   = (pplev(i,k)-pplev(i,k+1)) / (zrho(i)*g)
274                  zfrac(i) = (zt(i,k)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds)
275                  zfrac(i) = MAX(zfrac(i), 0.0)
276                  zfrac(i) = MIN(zfrac(i), 1.0)
277                  zneb(i)  = MAX(rneb(i,k), seuil_neb)
278               ENDIF
279           ENDDO
280
281 !recalculate liquid water particle radii
282           call h2o_cloudrad(ngrid,nlayer,ql,reffh2oliq,reffh2oice)
283
284           SELECT CASE(precip_scheme)
285 !precip scheme from Sundquist 78
286           CASE(2)
287
288            DO n = 1, ninter
289               DO i = 1, ngrid
290                  IF (rneb(i,k).GT.0.0) THEN
291                     ! this is the ONLY place where zneb, precip_timescale and cloud_sat are used
292
293                     zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale)) * zoliq(i)      &
294                          * (1.0-EXP(-(zoliq(i)/zneb(i)/cloud_sat)**2)) * zfrac(i)
295                     zrhol(i) = zrho(i) * zoliq(i) / zneb(i)
296                     zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i)    &
297                          *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly...
298                     ztot(i)  = zchau(i) + zfroi(i)
299
300                     IF (zneb(i).EQ.seuil_neb) ztot(i) = 0.0
301                     ztot(i)  = MIN(MAX(ztot(i),0.0),zoliq(i))
302                     zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0)
303
304                  ENDIF
305               ENDDO
306            ENDDO           
307
308 !precip scheme modified from Sundquist 78 (in q**3)
309           CASE(3)         
310           
311            DO n = 1, ninter
312               DO i = 1, ngrid
313                  IF (rneb(i,k).GT.0.0) THEN
314                     ! this is the ONLY place where zneb, precip_timescale and cloud_sat are used
315
316                     zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale*cloud_sat**2)) * (zoliq(i)/zneb(i))**3 
317                     zrhol(i) = zrho(i) * zoliq(i) / zneb(i)
318                     zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i)    &
319                          *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly...
320                     ztot(i)  = zchau(i) + zfroi(i)
321
322                     IF (zneb(i).EQ.seuil_neb) ztot(i) = 0.0
323                     ztot(i)  = MIN(MAX(ztot(i),0.0),zoliq(i))
324                     zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0)
325
326                  ENDIF
327               ENDDO
328            ENDDO           
329
330 !precip scheme modified from Boucher 95
331           CASE(4)
332
333            DO n = 1, ninter
334               DO i = 1, ngrid
335                  IF (rneb(i,k).GT.0.0) THEN
336                     ! this is the ONLY place where zneb and c1 are used
337
338                     zchau(i) = ptimestep/FLOAT(ninter) *c1* zrho(i) &
339                                    *(zoliq(i)/zneb(i))**2*reffh2oliq(i,k)*zneb(i)* zfrac(i) 
340                     zrhol(i) = zrho(i) * zoliq(i) / zneb(i)
341                     zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i)    &
342                          *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly...
343                     ztot(i)  = zchau(i) + zfroi(i)
344
345                     IF (zneb(i).EQ.seuil_neb) ztot(i) = 0.0
346                     ztot(i)  = MIN(MAX(ztot(i),0.0),zoliq(i))
347                     zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0)
348
349                  ENDIF
350               ENDDO
351            ENDDO           
352
353           END SELECT ! precip_scheme
354
355            ! Change in cloud density and surface H2O values
356            DO i = 1, ngrid
357               IF (rneb(i,k).GT.0.0) THEN
358                  d_ql(i,k) = (zoliq(i) - ql(i,k))!/ptimestep
359                  zrfl(i)   = zrfl(i)+ MAX(ql(i,k)-zoliq(i),0.0)*l2c(i,k)/ptimestep
360               ENDIF
361            ENDDO
362
363
364         endif ! if precip_scheme=1
365
366      ENDDO ! of DO k = nlayer, 1, -1
367
368!     Rain or snow on the ground
369      DO i = 1, ngrid
370         if(zrfl(i).lt.0.0)then
371            print*,'Droplets of negative rain are falling...'
372            call abort_physiq
373         endif
374         IF (t(i,1) .LT. T_h2O_ice_liq) THEN
375            dqssnow(i) = zrfl(i)
376            dqsrain(i) = 0.0
377         ELSE
378            dqssnow(i) = 0.0
379            dqsrain(i) = zrfl(i) ! liquid water = ice for now
380         ENDIF
381      ENDDO
382
383!     now subroutine -----> GCM variables
384      if (evap_prec) then
385        dqrain(1:ngrid,1:nlayer,i_vap)=d_q(1:ngrid,1:nlayer)/ptimestep
386        d_t(1:ngrid,1:nlayer)=d_t(1:ngrid,1:nlayer)/ptimestep
387      else
388        dqrain(1:ngrid,1:nlayer,i_vap)=0.0
389        d_t(1:ngrid,1:nlayer)=0.0
390      endif
391      dqrain(1:ngrid,1:nlayer,i_ice) = d_ql(1:ngrid,1:nlayer)/ptimestep
392
393    end subroutine rain
Note: See TracBrowser for help on using the repository browser.