source: trunk/LMDZ.GENERIC/libf/phystd/rain.F90 @ 1704

Last change on this file since 1704 was 1521, checked in by emillour, 9 years ago

All GCMs: Updates to make planetary codes (+Earth) setups converge.

  • Made a "phy_common" directory in libf, to contain routines common (wrt structural nature of underlying code/grid) to all LMDZ-related physics packages.
  • moved all "mod_phys_*" and "mod_grid_phy_lmdz" files from dynlonlat_phylonlat to "phy_common"
  • moved "ioipsl_getincom_p.F90 from "misc" to "phy_common" and modified it to match Earth GCM version and renamed it ioipsl_getin_p_mod.F90
  • added an "abort_physics" (as in Earth GCM) in "phy_common"
  • added a "print_control_mod.F90 (as in Earth GCM) in phy_common
  • made similar changes in LMDZ.GENERIC and LMDZ.MARS

EM

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