source: trunk/LMDZ.GENERIC/libf/phystd/rain_generic.F90 @ 3058

Last change on this file since 3058 was 3043, checked in by jleconte, 15 months ago

fixed some bugs and errors

File size: 20.1 KB
Line 
1subroutine rain_generic(ngrid,nlayer,nq,ptimestep,pplev,pplay,pphi,t,pdt,pq,pdq,d_t,dq_rain_generic_vap,dq_rain_generic_cld,dqsrain_generic,dqssnow_generic,reevap_precip,rneb)
2
3
4   use ioipsl_getin_p_mod, only: getin_p
5   use generic_cloud_common_h
6   use watercommon_h, only: T_h2O_ice_liq,T_h2O_ice_clouds,rhowater
7   ! T_h2O_ice_clouds,rhowater  are only used for precip_scheme_generic >=2
8   use radii_mod, only: h2o_cloudrad ! only used for precip_scheme_generic >=2
9   use tracer_h
10   use comcstfi_mod, only: g, r, cpp
11   use generic_tracer_index_mod, only: generic_tracer_index
12   implicit none
13 
14   !==================================================================
15   !     
16   !     Purpose
17   !     -------
18   !     Calculates precipitation for generic condensable tracers, using simplified microphysics.
19   !
20   !     GCS : generic condensable specie
21   !     
22   !     Authors
23   !     -------
24   !     Adapted from rain.F90 by Noé Clément (2022)
25   !     
26   !==================================================================
27 
28   ! Arguments
29   integer,intent(in) :: ngrid ! number of atmospheric columns
30   integer,intent(in) :: nlayer ! number of atmospheric layers
31   integer,intent(in) :: nq ! number of tracers
32   real,intent(in) :: ptimestep    ! time interval
33   real,intent(in) :: pplev(ngrid,nlayer+1) ! inter-layer pressure (Pa)
34   real,intent(in) :: pplay(ngrid,nlayer)   ! mid-layer pressure (Pa)
35   real,intent(in) :: pphi(ngrid,nlayer)   ! mid-layer geopotential
36   real,intent(in) :: t(ngrid,nlayer) ! input temperature (K)
37   real,intent(in) :: pdt(ngrid,nlayer) ! input tendency on temperature (K/s)     
38   real,intent(in) :: pq(ngrid,nlayer,nq)  ! tracers (kg/kg)
39   real,intent(in) :: pdq(ngrid,nlayer,nq) ! input tendency on tracers
40   real,intent(out) :: d_t(ngrid,nlayer) ! temperature tendency (K/s)
41   real,intent(out) :: dq_rain_generic_vap(ngrid,nlayer,nq) ! tendency of GCS tracers precipitation (kg/kg.s-1) - vapor
42   real,intent(out) :: dq_rain_generic_cld(ngrid,nlayer,nq) ! tendency of GCS tracers precipitation (kg/kg.s-1) - cloud
43   real,intent(out) :: dqsrain_generic(ngrid,nq)  ! rain flux at the surface (kg.m-2.s-1)
44   real,intent(out) :: dqssnow_generic(ngrid,nq)  ! snow flux at the surface (kg.m-2.s-1)
45   real,intent(out) :: reevap_precip(ngrid,nq)  ! re-evaporation flux of precipitation integrated over the atmospheric column (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 GCS (Kg/Kg)
50   real q(ngrid,nlayer)          ! specific humidity (Kg/Kg)
51   real d_q(ngrid,nlayer)        ! GCS vapor increment
52   real d_ql(ngrid,nlayer)       ! liquid GCS / ice increment
53
54   integer igcm_generic_vap, igcm_generic_ice ! index of the vap and ice fo GCS
55
56   real, save :: RCPD ! equal to cpp
57
58   real, save :: metallicity ! metallicity of planet --- is not used here, but necessary to call function Psat_generic
59   !$OMP THREADPRIVATE(metallicity)
60
61   ! Subroutine options
62   real,parameter :: seuil_neb=0.001  ! Nebulosity threshold
63
64   integer,save :: precip_scheme_generic      ! id number for precipitaion scheme
65   
66   ! for simple scheme  (precip_scheme_generic=1)
67   real,save :: rainthreshold_generic                ! Precipitation threshold in simple scheme
68   
69   ! for sundquist scheme  (precip_scheme_generic=2-3)
70   real,save :: cloud_sat_generic                    ! Precipitation threshold in non simple scheme
71   real,save :: precip_timescale_generic             ! Precipitation timescale
72   
73   ! for Boucher scheme  (precip_scheme_generic=4)
74   real,save :: Cboucher_generic ! Precipitation constant in Boucher 95 scheme
75   real,parameter :: Kboucher=1.19E8
76   real,save :: c1
77   
78   !$OMP THREADPRIVATE(precip_scheme_generic,rainthreshold_generic,cloud_sat_generic,precip_timescale_generic,Cboucher_generic,c1)
79
80   integer,parameter :: ninter=5 ! only used for precip_scheme_generic >=2
81
82   logical,save :: evap_prec_generic ! Does the rain evaporate ?
83   REAL,SAVE :: evap_coeff_generic ! multiplication evaporation constant. 1. gives the model of gregory et al.
84   !$OMP THREADPRIVATE(evap_prec_generic,evap_coeff_generic)
85
86   ! for simple scheme : precip_scheme_generic=1
87   real lconvert
88
89   ! Local variables
90   integer i, k, n, iq
91   REAL zqs(ngrid,nlayer), dqsat(ngrid,nlayer), dlnpsat(ngrid,nlayer), Tsat(ngrid,nlayer), zdelta, zcor
92   REAL precip_rate(ngrid), precip_rate_tmp(ngrid)    ! local precipitation rate in kg of condensed GCS per m^2 per s.
93   REAL zqev, zqevt
94
95   real zoliq(ngrid)
96   real zdz(ngrid),zrho(ngrid),ztot(ngrid), zrhol(ngrid)
97   real zchau(ngrid),zfroi(ngrid),zfrac(ngrid),zneb(ngrid)
98
99   real reffh2oliq(ngrid,nlayer),reffh2oice(ngrid,nlayer)
100
101   real t_tmp, p_tmp, psat_tmp
102   real tnext(ngrid,nlayer)
103
104   real dmass(ngrid,nlayer)   ! mass of air in each grid cell
105   real dWtot
106 
107 
108   ! Indices of GCS vapour and GCS ice tracers
109   integer, save :: i_vap_generic=0  ! GCS vapour
110   integer, save :: i_ice_generic=0  ! GCS ice
111   !$OMP THREADPRIVATE(i_vap_generic,i_ice_generic)
112
113   LOGICAL,save :: firstcall=.true.
114   !$OMP THREADPRIVATE(firstcall)
115
116   ! to call only one time the ice/vap pair of a tracer
117   logical call_ice_vap_generic
118
119   ! Online functions
120   real fallv, fall2v, zzz ! falling speed of ice crystals
121   fallv (zzz) = 3.29 * ((zzz)**0.16)
122   fall2v (zzz) =10.6 * ((zzz)**0.31)  !for use with radii
123
124   ! Let's loop on tracers
125
126
127   dq_rain_generic_vap(1:ngrid,1:nlayer,1:nq) = 0.0
128   dq_rain_generic_cld(1:ngrid,1:nlayer,1:nq) = 0.0
129   dqsrain_generic(:,:) = 0.0
130
131   do iq=1,nq
132
133      call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic)
134
135      if(call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer
136     
137         m=constants_mass(iq)
138         delta_vapH=constants_delta_vapH(iq)
139         Tref=constants_Tref(iq)
140         Pref=constants_Pref(iq)
141         epsi_generic=constants_epsi_generic(iq)
142         RLVTT_generic=constants_RLVTT_generic(iq)
143
144         RCPD = cpp
145         
146
147         if (firstcall) then
148
149            metallicity=0.0 ! default value --- is not used here but necessary to call function Psat_generic
150            call getin_p("metallicity",metallicity) ! --- is not used here but necessary to call function Psat_generic
151
152            i_vap_generic=igcm_generic_vap
153            i_ice_generic=igcm_generic_ice
154           
155            write(*,*) "rain: i_ice_generic=",i_ice_generic
156            write(*,*) "      i_vap_generic=",i_vap_generic
157
158            write(*,*) "re-evaporate precipitations?"
159            evap_prec_generic=.true. ! default value
160            call getin_p("evap_prec_generic",evap_prec_generic)
161            write(*,*) " evap_prec_generic = ",evap_prec_generic
162
163            if (evap_prec_generic) then
164               write(*,*) "multiplicative constant in reevaporation"
165               evap_coeff_generic=1.   ! default value
166               call getin_p("evap_coeff_generic",evap_coeff_generic)
167               write(*,*) " evap_coeff_generic = ",evap_coeff_generic   
168            end if
169
170            write(*,*) "Precipitation scheme to use?"
171            precip_scheme_generic=1 ! default value
172            call getin_p("precip_scheme_generic",precip_scheme_generic)
173            write(*,*) " precip_scheme_generic = ",precip_scheme_generic
174
175            if (precip_scheme_generic.eq.1) then
176               write(*,*) "rainthreshold_generic in simple scheme?"
177               rainthreshold_generic=0. ! default value
178               call getin_p("rainthreshold_generic",rainthreshold_generic)
179               write(*,*) " rainthreshold_generic = ",rainthreshold_generic
180           
181            !else
182            !   write(*,*) "precip_scheme_generic = ", precip_scheme_generic
183            !   write(*,*) "this precip_scheme_generic has not been implemented yet,"
184            !   write(*,*) "only precip_scheme_generic = 1 has been implemented"
185               
186            else if (precip_scheme_generic.eq.2.or.precip_scheme_generic.eq.3) then
187               
188               write(*,*) "cloud GCS saturation level in non simple scheme?"
189               cloud_sat_generic=2.6e-4   ! default value
190               call getin_p("cloud_sat_generic",cloud_sat_generic)
191               write(*,*) " cloud_sat_generic = ",cloud_sat_generic
192           
193               write(*,*) "precipitation timescale in non simple scheme?"
194               precip_timescale_generic=3600.  ! default value
195               call getin_p("precip_timescale_generic",precip_timescale_generic)
196               write(*,*) " precip_timescale_generic = ",precip_timescale_generic
197
198            else if (precip_scheme_generic.eq.4) then
199               
200               write(*,*) "multiplicative constant in Boucher 95 precip scheme"
201               Cboucher_generic=1.   ! default value
202               call getin_p("Cboucher_generic",Cboucher_generic)
203               write(*,*) " Cboucher_generic = ",Cboucher_generic       
204               
205               c1=1.00*1.097/rhowater*Cboucher_generic*Kboucher 
206
207            endif
208
209            PRINT*, 'in rain_generic.F, ninter=', ninter
210
211            firstcall = .false.
212
213         endif ! of if (firstcall)
214
215         ! GCM -----> subroutine variables
216         do k = 1, nlayer
217            do i = 1, ngrid
218
219               zt(i,k)   = t(i,k)+pdt(i,k)*ptimestep ! a big fat bug was here
220               q(i,k)    = pq(i,k,i_vap_generic)+pdq(i,k,i_vap_generic)*ptimestep
221               ql(i,k)   = pq(i,k,i_ice_generic)+pdq(i,k,i_ice_generic)*ptimestep
222
223               !q(i,k)    = pq(i,k,i_vap_generic)!+pdq(i,k,i_vap_generic)
224               !ql(i,k)   = pq(i,k,i_ice_generic)!+pdq(i,k,i_ice_generic)
225
226               if(q(i,k).lt.0.)then ! if this is not done, we don't conserve GCS
227                  q(i,k)=0. ! vap
228               endif
229               if(ql(i,k).lt.0.)then
230                  ql(i,k)=0. ! ice
231               endif
232            enddo
233         enddo
234
235         ! Initialise the outputs
236         d_t(1:ngrid,1:nlayer) = 0.0
237         d_q(1:ngrid,1:nlayer) = 0.0
238         d_ql(1:ngrid,1:nlayer) = 0.0
239         precip_rate(1:ngrid) = 0.0
240         precip_rate_tmp(1:ngrid) = 0.0
241
242         dq_rain_generic_vap(1:ngrid,1:nlayer,1:nq) = 0.0
243         dq_rain_generic_cld(1:ngrid,1:nlayer,1:nq) = 0.0
244         dqsrain_generic(1:ngrid,1:nq) = 0.0
245         dqssnow_generic(1:ngrid,1:nq) = 0.0
246
247         ! calculate saturation mixing ratio
248         do k = 1, nlayer
249            do i = 1, ngrid
250               p_tmp = pplay(i,k)
251               call Psat_generic(zt(i,k),p_tmp,metallicity,psat_tmp,zqs(i,k)) ! calculates psat_tmp & zqs(i,k)
252               ! metallicity --- is not used here but necessary to call function Psat_generic
253               call Lcpdqsat_generic(zt(i,k),p_tmp,psat_tmp,zqs(i,k),dqsat(i,k),dlnpsat(i,k)) ! calculates dqsat(i,k) & dlnpsat(i,k)
254               call Tsat_generic(p_tmp,Tsat(i,k)) ! calculates Tsat(i,k)
255
256            enddo
257         enddo
258
259         ! get column / layer conversion factor
260         do k = 1, nlayer
261            do i = 1, ngrid
262               dmass(i,k)=(pplev(i,k)-pplev(i,k+1))/g ! mass per m² in each layer
263            enddo
264         enddo
265
266         ! Vertical loop (from top to bottom)
267         ! We carry the rain with us and calculate that added by warm/cold precipitation
268         ! processes and that subtracted by evaporation at each level.
269         ! We go from a layer to the next layer below and make the rain fall
270         do k = nlayer, 1, -1
271
272            if (evap_prec_generic) then ! note no rneb dependence!
273               do i = 1, ngrid
274                  if (precip_rate(i) .GT.0.) then ! if rain from upper layers has fallen in the current layer box
275
276                     if(zt(i,k).gt.Tsat(i,k))then ! if temperature of the layer box is greater than Tsat
277                        ! treat the case where all liquid water should boil
278                        zqev=MIN((zt(i,k)-Tsat(i,k))*RCPD*dmass(i,k)/RLVTT_generic/ptimestep,precip_rate(i)) ! on évapore
279                        precip_rate(i)=MAX(precip_rate(i)-zqev,0.) ! we withdraw from precip_rate the evaporated part
280                        d_q(i,k)=zqev/dmass(i,k)*ptimestep ! quantité évaporée
281                        d_t(i,k) = - d_q(i,k) * RLVTT_generic/RCPD
282                     else
283                        zqev = MAX (0.0, (zqs(i,k)-q(i,k)))*dmass(i,k)/(ptimestep*(1.d0+dqsat(i,k)))
284                        !max evaporation to reach saturation implictly accounting for temperature reduction
285                        zqevt= MAX (0.0, evap_coeff_generic*2.0e-5*(1.0-q(i,k)/zqs(i,k))    & !default was 2.e-5
286                              *sqrt(precip_rate(i))*dmass(i,k)/pplay(i,k)*zt(i,k)*R) ! BC modif here, is it R or r/(mu/1000) ?
287                        zqev  = MIN (zqev, zqevt)
288                        zqev  = MAX (zqev, 0.0) ! a priori inutile d'après les précédentes lignes
289                        precip_rate_tmp(i)= precip_rate(i) - zqev
290                        precip_rate_tmp(i)= max(precip_rate_tmp(i),0.0)
291
292                        d_q(i,k) = - (precip_rate_tmp(i)-precip_rate(i))/dmass(i,k)*ptimestep
293                        d_t(i,k) = - d_q(i,k) * RLVTT_generic/RCPD
294                        precip_rate(i)  = precip_rate_tmp(i)
295                     end if
296#ifdef MESOSCALE
297                     d_t(i,k) = d_t(i,k)+(pphi(i,k+1)-pphi(i,k))*precip_rate(i)*ptimestep/(RCPD*dmass(i,k))
298                     ! JL22. Accounts for gravitational energy of falling precipitations (probably not to be used in the GCM
299                     !   where the counterpart is not included in the dynamics.)
300#endif
301                  endif ! of if (precip_rate(i) .GT.0.)
302               enddo
303            endif ! of if (evap_prec_generic)
304
305            zoliq(1:ngrid) = 0.0
306
307
308            if(precip_scheme_generic.eq.1)then
309
310               do i = 1, ngrid
311                 
312                  lconvert=rainthreshold_generic
313
314                  if (ql(i,k).gt.1.e-9) then
315                     zneb(i)  = MAX(rneb(i,k), seuil_neb) ! in mesoscale rneb = 0 or 1
316                     if ((ql(i,k)/zneb(i)).gt.lconvert)then ! precipitate!
317                        d_ql(i,k) = -MAX((ql(i,k)-lconvert*zneb(i)),0.0)
318                        precip_rate(i)   = precip_rate(i) - d_ql(i,k)*dmass(i,k)/ptimestep
319                     endif
320                  endif
321               enddo
322
323            elseif (precip_scheme_generic.ge.2) then
324           
325               do i = 1, ngrid
326                  if (rneb(i,k).GT.0.0) then
327                     zoliq(i) = ql(i,k)
328                     zrho(i)  = pplay(i,k) / ( zt(i,k) * R )
329                     zdz(i)   = (pplev(i,k)-pplev(i,k+1)) / (zrho(i)*g)
330                     zfrac(i) = (zt(i,k)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds)
331                     zfrac(i) = MAX(zfrac(i), 0.0)
332                     zfrac(i) = MIN(zfrac(i), 1.0)
333                     zneb(i)  = MAX(rneb(i,k), seuil_neb)
334                  endif
335               enddo
336
337               !recalculate liquid GCS particle radii
338               call h2o_cloudrad(ngrid,nlayer,ql,reffh2oliq,reffh2oice)
339
340               SELECT CASE(precip_scheme_generic)
341
342               !precip scheme from Sundquist 78
343               CASE(2)
344
345                  do n = 1, ninter
346                     do i = 1, ngrid
347                        if (rneb(i,k).GT.0.0) then
348                           ! this is the ONLY place where zneb, precip_timescale_generic and cloud_sat_generic are used
349
350                           zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale_generic)) * zoliq(i)      &
351                                    * (1.0-EXP(-(zoliq(i)/zneb(i)/cloud_sat_generic)**2)) * zfrac(i)
352                           zrhol(i) = zrho(i) * zoliq(i) / zneb(i)
353                           zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i)    &
354                                    *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly...
355                           ztot(i)  = zchau(i) + zfroi(i)
356
357                           if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0
358                              ztot(i)  = MIN(MAX(ztot(i),0.0),zoliq(i))
359                              zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0)
360                           endif
361                     enddo
362                  enddo           
363
364               !precip scheme modified from Sundquist 78 (in q**3)
365               CASE(3)   
366                  do n = 1, ninter
367                     do i = 1, ngrid
368                        if (rneb(i,k).GT.0.0) then
369                        ! this is the ONLY place where zneb, precip_timescale_generic and cloud_sat_generic are used
370
371                           zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale_generic*cloud_sat_generic**2)) * (zoliq(i)/zneb(i))**3 
372                           zrhol(i) = zrho(i) * zoliq(i) / zneb(i)
373                           zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i)    &
374                                       *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly...
375                           ztot(i)  = zchau(i) + zfroi(i)
376
377                           if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0
378                              ztot(i)  = MIN(MAX(ztot(i),0.0),zoliq(i))
379                              zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0)
380                           endif
381                     enddo
382                  enddo           
383
384               !precip scheme modified from Boucher 95
385               CASE(4)
386                  do n = 1, ninter
387                     do i = 1, ngrid
388                        if (rneb(i,k).GT.0.0) then
389                        ! this is the ONLY place where zneb and c1 are used
390                           zchau(i) = ptimestep/FLOAT(ninter) *c1* zrho(i) &
391                                       *(zoliq(i)/zneb(i))**2*reffh2oliq(i,k)*zneb(i)* zfrac(i)
392                           zrhol(i) = zrho(i) * zoliq(i) / zneb(i)
393                           zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i)    &
394                                       *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly...
395                           ztot(i)  = zchau(i) + zfroi(i)
396
397                           if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0
398                           ztot(i)  = MIN(MAX(ztot(i),0.0),zoliq(i))
399                           zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0)
400                        endif
401                     enddo
402                  enddo           
403
404               END SELECT ! precip_scheme_generic
405
406               ! Change in cloud density and surface GCS values
407               do i = 1, ngrid
408                  if (rneb(i,k).GT.0.0) then
409                     d_ql(i,k) = (zoliq(i) - ql(i,k))!/ptimestep
410                     precip_rate(i)   = precip_rate(i)+ MAX(ql(i,k)-zoliq(i),0.0)*dmass(i,k)/ptimestep
411                  endif
412               enddo
413
414            endif ! if precip_scheme_generic=1
415
416         enddo ! of do k = nlayer, 1, -1
417
418         ! Rain or snow on the ground
419         do i = 1, ngrid
420            if(precip_rate(i).lt.0.0)then
421               print*,'Droplets of negative rain are falling...'
422               call abort
423            endif
424            if (t(i,1) .LT. T_h2O_ice_liq) then
425               dqssnow_generic(i,i_ice_generic) = precip_rate(i)
426               dqsrain_generic(i,i_ice_generic) = 0.0
427            else
428               dqssnow_generic(i,i_ice_generic) = 0.0
429               dqsrain_generic(i,i_ice_generic) = precip_rate(i) ! liquid water = ice for now
430            endif
431
432            ! For now we force :
433            dqsrain_generic(i,i_ice_generic) = precip_rate(i)
434            dqssnow_generic(i,i_ice_generic) = 0.0
435         enddo
436
437         ! now subroutine -----> GCM variables
438         if (evap_prec_generic) then
439            dq_rain_generic_vap(1:ngrid,1:nlayer,i_vap_generic)=d_q(1:ngrid,1:nlayer)/ptimestep
440            d_t(1:ngrid,1:nlayer)=d_t(1:ngrid,1:nlayer)/ptimestep
441            do i=1,ngrid
442               reevap_precip(i,i_vap_generic)=0.
443               do k=1,nlayer
444                  reevap_precip(i,i_vap_generic)=reevap_precip(i,i_vap_generic)+dq_rain_generic_vap(i,k,i_vap_generic)*dmass(i,k)
445               enddo
446            enddo
447         else
448            dq_rain_generic_vap(1:ngrid,1:nlayer,i_vap_generic)=0.0
449            d_t(1:ngrid,1:nlayer)=0.0
450         endif
451         dq_rain_generic_cld(1:ngrid,1:nlayer,i_ice_generic) = d_ql(1:ngrid,1:nlayer)/ptimestep
452
453      end if ! if(call_ice_vap_generic)
454
455   end do ! do iq=1,nq loop on tracers
456
457end subroutine rain_generic
Note: See TracBrowser for help on using the repository browser.