1 | subroutine 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 | do iq=1,nq |
---|
127 | |
---|
128 | call generic_tracer_index(nq,iq,igcm_generic_vap,igcm_generic_ice,call_ice_vap_generic) |
---|
129 | |
---|
130 | if(call_ice_vap_generic) then ! to call only one time the ice/vap pair of a tracer |
---|
131 | |
---|
132 | m=constants_mass(iq) |
---|
133 | delta_vapH=constants_delta_vapH(iq) |
---|
134 | Tref=constants_Tref(iq) |
---|
135 | Pref=constants_Pref(iq) |
---|
136 | epsi_generic=constants_epsi_generic(iq) |
---|
137 | RLVTT_generic=constants_RLVTT_generic(iq) |
---|
138 | |
---|
139 | RCPD = cpp |
---|
140 | |
---|
141 | |
---|
142 | if (firstcall) then |
---|
143 | |
---|
144 | metallicity=0.0 ! default value --- is not used here but necessary to call function Psat_generic |
---|
145 | call getin_p("metallicity",metallicity) ! --- is not used here but necessary to call function Psat_generic |
---|
146 | |
---|
147 | i_vap_generic=igcm_generic_vap |
---|
148 | i_ice_generic=igcm_generic_ice |
---|
149 | |
---|
150 | write(*,*) "rain: i_ice_generic=",i_ice_generic |
---|
151 | write(*,*) " i_vap_generic=",i_vap_generic |
---|
152 | |
---|
153 | write(*,*) "re-evaporate precipitations?" |
---|
154 | evap_prec_generic=.true. ! default value |
---|
155 | call getin_p("evap_prec_generic",evap_prec_generic) |
---|
156 | write(*,*) " evap_prec_generic = ",evap_prec_generic |
---|
157 | |
---|
158 | if (evap_prec_generic) then |
---|
159 | write(*,*) "multiplicative constant in reevaporation" |
---|
160 | evap_coeff_generic=1. ! default value |
---|
161 | call getin_p("evap_coeff_generic",evap_coeff_generic) |
---|
162 | write(*,*) " evap_coeff_generic = ",evap_coeff_generic |
---|
163 | end if |
---|
164 | |
---|
165 | write(*,*) "Precipitation scheme to use?" |
---|
166 | precip_scheme_generic=1 ! default value |
---|
167 | call getin_p("precip_scheme_generic",precip_scheme_generic) |
---|
168 | write(*,*) " precip_scheme_generic = ",precip_scheme_generic |
---|
169 | |
---|
170 | if (precip_scheme_generic.eq.1) then |
---|
171 | write(*,*) "rainthreshold_generic in simple scheme?" |
---|
172 | rainthreshold_generic=0. ! default value |
---|
173 | call getin_p("rainthreshold_generic",rainthreshold_generic) |
---|
174 | write(*,*) " rainthreshold_generic = ",rainthreshold_generic |
---|
175 | |
---|
176 | !else |
---|
177 | ! write(*,*) "precip_scheme_generic = ", precip_scheme_generic |
---|
178 | ! write(*,*) "this precip_scheme_generic has not been implemented yet," |
---|
179 | ! write(*,*) "only precip_scheme_generic = 1 has been implemented" |
---|
180 | |
---|
181 | else if (precip_scheme_generic.eq.2.or.precip_scheme_generic.eq.3) then |
---|
182 | |
---|
183 | write(*,*) "cloud GCS saturation level in non simple scheme?" |
---|
184 | cloud_sat_generic=2.6e-4 ! default value |
---|
185 | call getin_p("cloud_sat_generic",cloud_sat_generic) |
---|
186 | write(*,*) " cloud_sat_generic = ",cloud_sat_generic |
---|
187 | |
---|
188 | write(*,*) "precipitation timescale in non simple scheme?" |
---|
189 | precip_timescale_generic=3600. ! default value |
---|
190 | call getin_p("precip_timescale_generic",precip_timescale_generic) |
---|
191 | write(*,*) " precip_timescale_generic = ",precip_timescale_generic |
---|
192 | |
---|
193 | else if (precip_scheme_generic.eq.4) then |
---|
194 | |
---|
195 | write(*,*) "multiplicative constant in Boucher 95 precip scheme" |
---|
196 | Cboucher_generic=1. ! default value |
---|
197 | call getin_p("Cboucher_generic",Cboucher_generic) |
---|
198 | write(*,*) " Cboucher_generic = ",Cboucher_generic |
---|
199 | |
---|
200 | c1=1.00*1.097/rhowater*Cboucher_generic*Kboucher |
---|
201 | |
---|
202 | endif |
---|
203 | |
---|
204 | PRINT*, 'in rain_generic.F, ninter=', ninter |
---|
205 | |
---|
206 | firstcall = .false. |
---|
207 | |
---|
208 | endif ! of if (firstcall) |
---|
209 | |
---|
210 | ! GCM -----> subroutine variables |
---|
211 | do k = 1, nlayer |
---|
212 | do i = 1, ngrid |
---|
213 | |
---|
214 | zt(i,k) = t(i,k)+pdt(i,k)*ptimestep ! a big fat bug was here |
---|
215 | q(i,k) = pq(i,k,i_vap_generic)+pdq(i,k,i_vap_generic)*ptimestep |
---|
216 | ql(i,k) = pq(i,k,i_ice_generic)+pdq(i,k,i_ice_generic)*ptimestep |
---|
217 | |
---|
218 | !q(i,k) = pq(i,k,i_vap_generic)!+pdq(i,k,i_vap_generic) |
---|
219 | !ql(i,k) = pq(i,k,i_ice_generic)!+pdq(i,k,i_ice_generic) |
---|
220 | |
---|
221 | if(q(i,k).lt.0.)then ! if this is not done, we don't conserve GCS |
---|
222 | q(i,k)=0. ! vap |
---|
223 | endif |
---|
224 | if(ql(i,k).lt.0.)then |
---|
225 | ql(i,k)=0. ! ice |
---|
226 | endif |
---|
227 | enddo |
---|
228 | enddo |
---|
229 | |
---|
230 | ! Initialise the outputs |
---|
231 | d_t(1:ngrid,1:nlayer) = 0.0 |
---|
232 | d_q(1:ngrid,1:nlayer) = 0.0 |
---|
233 | d_ql(1:ngrid,1:nlayer) = 0.0 |
---|
234 | precip_rate(1:ngrid) = 0.0 |
---|
235 | precip_rate_tmp(1:ngrid) = 0.0 |
---|
236 | |
---|
237 | ! calculate saturation mixing ratio |
---|
238 | do k = 1, nlayer |
---|
239 | do i = 1, ngrid |
---|
240 | p_tmp = pplay(i,k) |
---|
241 | call Psat_generic(zt(i,k),p_tmp,metallicity,psat_tmp,zqs(i,k)) ! calculates psat_tmp & zqs(i,k) |
---|
242 | ! metallicity --- is not used here but necessary to call function Psat_generic |
---|
243 | 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) |
---|
244 | call Tsat_generic(p_tmp,Tsat(i,k)) ! calculates Tsat(i,k) |
---|
245 | |
---|
246 | enddo |
---|
247 | enddo |
---|
248 | |
---|
249 | ! get column / layer conversion factor |
---|
250 | do k = 1, nlayer |
---|
251 | do i = 1, ngrid |
---|
252 | dmass(i,k)=(pplev(i,k)-pplev(i,k+1))/g ! mass per m² in each layer |
---|
253 | enddo |
---|
254 | enddo |
---|
255 | |
---|
256 | ! Vertical loop (from top to bottom) |
---|
257 | ! We carry the rain with us and calculate that added by warm/cold precipitation |
---|
258 | ! processes and that subtracted by evaporation at each level. |
---|
259 | ! We go from a layer to the next layer below and make the rain fall |
---|
260 | do k = nlayer, 1, -1 |
---|
261 | |
---|
262 | if (evap_prec_generic) then ! note no rneb dependence! |
---|
263 | do i = 1, ngrid |
---|
264 | if (precip_rate(i) .GT.0.) then ! if rain from upper layers has fallen in the current layer box |
---|
265 | |
---|
266 | if(zt(i,k).gt.Tsat(i,k))then ! if temperature of the layer box is greater than Tsat |
---|
267 | ! treat the case where all liquid water should boil |
---|
268 | zqev=MIN((zt(i,k)-Tsat(i,k))*RCPD*dmass(i,k)/RLVTT_generic/ptimestep,precip_rate(i)) ! on évapore |
---|
269 | precip_rate(i)=MAX(precip_rate(i)-zqev,0.) ! we withdraw from precip_rate the evaporated part |
---|
270 | d_q(i,k)=zqev/dmass(i,k)*ptimestep ! quantité évaporée |
---|
271 | d_t(i,k) = - d_q(i,k) * RLVTT_generic/RCPD |
---|
272 | else |
---|
273 | zqev = MAX (0.0, (zqs(i,k)-q(i,k)))*dmass(i,k)/(ptimestep*(1.d0+dqsat(i,k))) |
---|
274 | !max evaporation to reach saturation implictly accounting for temperature reduction |
---|
275 | zqevt= MAX (0.0, evap_coeff_generic*2.0e-5*(1.0-q(i,k)/zqs(i,k)) & !default was 2.e-5 |
---|
276 | *sqrt(precip_rate(i))*dmass(i,k)/pplay(i,k)*zt(i,k)*R) ! BC modif here, is it R or r/(mu/1000) ? |
---|
277 | zqev = MIN (zqev, zqevt) |
---|
278 | zqev = MAX (zqev, 0.0) ! a priori inutile d'après les précédentes lignes |
---|
279 | precip_rate_tmp(i)= precip_rate(i) - zqev |
---|
280 | precip_rate_tmp(i)= max(precip_rate_tmp(i),0.0) |
---|
281 | |
---|
282 | d_q(i,k) = - (precip_rate_tmp(i)-precip_rate(i))/dmass(i,k)*ptimestep |
---|
283 | d_t(i,k) = - d_q(i,k) * RLVTT_generic/RCPD |
---|
284 | precip_rate(i) = precip_rate_tmp(i) |
---|
285 | end if |
---|
286 | #ifdef MESOSCALE |
---|
287 | d_t(i,k) = d_t(i,k)+(pphi(i,k+1)-pphi(i,k))*precip_rate(i)*ptimestep/(RCPD*dmass(i,k)) |
---|
288 | ! JL22. Accounts for gravitational energy of falling precipitations (probably not to be used in the GCM |
---|
289 | ! where the counterpart is not included in the dynamics.) |
---|
290 | #endif |
---|
291 | endif ! of if (precip_rate(i) .GT.0.) |
---|
292 | enddo |
---|
293 | endif ! of if (evap_prec_generic) |
---|
294 | |
---|
295 | zoliq(1:ngrid) = 0.0 |
---|
296 | |
---|
297 | |
---|
298 | if(precip_scheme_generic.eq.1)then |
---|
299 | |
---|
300 | do i = 1, ngrid |
---|
301 | |
---|
302 | lconvert=rainthreshold_generic |
---|
303 | |
---|
304 | if (ql(i,k).gt.1.e-9) then |
---|
305 | zneb(i) = MAX(rneb(i,k), seuil_neb) ! in mesoscale rneb = 0 or 1 |
---|
306 | if ((ql(i,k)/zneb(i)).gt.lconvert)then ! precipitate! |
---|
307 | d_ql(i,k) = -MAX((ql(i,k)-lconvert*zneb(i)),0.0) |
---|
308 | precip_rate(i) = precip_rate(i) - d_ql(i,k)*dmass(i,k)/ptimestep |
---|
309 | endif |
---|
310 | endif |
---|
311 | enddo |
---|
312 | |
---|
313 | elseif (precip_scheme_generic.ge.2) then |
---|
314 | |
---|
315 | do i = 1, ngrid |
---|
316 | if (rneb(i,k).GT.0.0) then |
---|
317 | zoliq(i) = ql(i,k) |
---|
318 | zrho(i) = pplay(i,k) / ( zt(i,k) * R ) |
---|
319 | zdz(i) = (pplev(i,k)-pplev(i,k+1)) / (zrho(i)*g) |
---|
320 | zfrac(i) = (zt(i,k)-T_h2O_ice_clouds) / (T_h2O_ice_liq-T_h2O_ice_clouds) |
---|
321 | zfrac(i) = MAX(zfrac(i), 0.0) |
---|
322 | zfrac(i) = MIN(zfrac(i), 1.0) |
---|
323 | zneb(i) = MAX(rneb(i,k), seuil_neb) |
---|
324 | endif |
---|
325 | enddo |
---|
326 | |
---|
327 | !recalculate liquid GCS particle radii |
---|
328 | call h2o_cloudrad(ngrid,nlayer,ql,reffh2oliq,reffh2oice) |
---|
329 | |
---|
330 | SELECT CASE(precip_scheme_generic) |
---|
331 | |
---|
332 | !precip scheme from Sundquist 78 |
---|
333 | CASE(2) |
---|
334 | |
---|
335 | do n = 1, ninter |
---|
336 | do i = 1, ngrid |
---|
337 | if (rneb(i,k).GT.0.0) then |
---|
338 | ! this is the ONLY place where zneb, precip_timescale_generic and cloud_sat_generic are used |
---|
339 | |
---|
340 | zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale_generic)) * zoliq(i) & |
---|
341 | * (1.0-EXP(-(zoliq(i)/zneb(i)/cloud_sat_generic)**2)) * zfrac(i) |
---|
342 | zrhol(i) = zrho(i) * zoliq(i) / zneb(i) |
---|
343 | zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i) & |
---|
344 | *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly... |
---|
345 | ztot(i) = zchau(i) + zfroi(i) |
---|
346 | |
---|
347 | if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0 |
---|
348 | ztot(i) = MIN(MAX(ztot(i),0.0),zoliq(i)) |
---|
349 | zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0) |
---|
350 | endif |
---|
351 | enddo |
---|
352 | enddo |
---|
353 | |
---|
354 | !precip scheme modified from Sundquist 78 (in q**3) |
---|
355 | CASE(3) |
---|
356 | do n = 1, ninter |
---|
357 | do i = 1, ngrid |
---|
358 | if (rneb(i,k).GT.0.0) then |
---|
359 | ! this is the ONLY place where zneb, precip_timescale_generic and cloud_sat_generic are used |
---|
360 | |
---|
361 | zchau(i) = (ptimestep/(FLOAT(ninter)*precip_timescale_generic*cloud_sat_generic**2)) * (zoliq(i)/zneb(i))**3 |
---|
362 | zrhol(i) = zrho(i) * zoliq(i) / zneb(i) |
---|
363 | zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i) & |
---|
364 | *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly... |
---|
365 | ztot(i) = zchau(i) + zfroi(i) |
---|
366 | |
---|
367 | if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0 |
---|
368 | ztot(i) = MIN(MAX(ztot(i),0.0),zoliq(i)) |
---|
369 | zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0) |
---|
370 | endif |
---|
371 | enddo |
---|
372 | enddo |
---|
373 | |
---|
374 | !precip scheme modified from Boucher 95 |
---|
375 | CASE(4) |
---|
376 | do n = 1, ninter |
---|
377 | do i = 1, ngrid |
---|
378 | if (rneb(i,k).GT.0.0) then |
---|
379 | ! this is the ONLY place where zneb and c1 are used |
---|
380 | zchau(i) = ptimestep/FLOAT(ninter) *c1* zrho(i) & |
---|
381 | *(zoliq(i)/zneb(i))**2*reffh2oliq(i,k)*zneb(i)* zfrac(i) |
---|
382 | zrhol(i) = zrho(i) * zoliq(i) / zneb(i) |
---|
383 | zfroi(i) = ptimestep/FLOAT(ninter)/zdz(i)*zoliq(i) & |
---|
384 | *fall2v(reffh2oice(i,k)) * (1.0-zfrac(i)) ! zfroi behaves oddly... |
---|
385 | ztot(i) = zchau(i) + zfroi(i) |
---|
386 | |
---|
387 | if (zneb(i).EQ.seuil_neb) ztot(i) = 0.0 |
---|
388 | ztot(i) = MIN(MAX(ztot(i),0.0),zoliq(i)) |
---|
389 | zoliq(i) = MAX(zoliq(i)-ztot(i), 0.0) |
---|
390 | endif |
---|
391 | enddo |
---|
392 | enddo |
---|
393 | |
---|
394 | END SELECT ! precip_scheme_generic |
---|
395 | |
---|
396 | ! Change in cloud density and surface GCS values |
---|
397 | do i = 1, ngrid |
---|
398 | if (rneb(i,k).GT.0.0) then |
---|
399 | d_ql(i,k) = (zoliq(i) - ql(i,k))!/ptimestep |
---|
400 | precip_rate(i) = precip_rate(i)+ MAX(ql(i,k)-zoliq(i),0.0)*dmass(i,k)/ptimestep |
---|
401 | endif |
---|
402 | enddo |
---|
403 | |
---|
404 | endif ! if precip_scheme_generic=1 |
---|
405 | |
---|
406 | enddo ! of do k = nlayer, 1, -1 |
---|
407 | |
---|
408 | ! Rain or snow on the ground |
---|
409 | do i = 1, ngrid |
---|
410 | if(precip_rate(i).lt.0.0)then |
---|
411 | print*,'Droplets of negative rain are falling...' |
---|
412 | call abort |
---|
413 | endif |
---|
414 | if (t(i,1) .LT. T_h2O_ice_liq) then |
---|
415 | dqssnow_generic(i,i_ice_generic) = precip_rate(i) |
---|
416 | dqsrain_generic(i,i_ice_generic) = 0.0 |
---|
417 | else |
---|
418 | dqssnow_generic(i,i_ice_generic) = 0.0 |
---|
419 | dqsrain_generic(i,i_ice_generic) = precip_rate(i) ! liquid water = ice for now |
---|
420 | endif |
---|
421 | |
---|
422 | ! For now we force : |
---|
423 | dqsrain_generic(i,i_ice_generic) = precip_rate(i) |
---|
424 | dqssnow_generic(i,i_ice_generic) = 0.0 |
---|
425 | enddo |
---|
426 | |
---|
427 | ! now subroutine -----> GCM variables |
---|
428 | if (evap_prec_generic) then |
---|
429 | dq_rain_generic_vap(1:ngrid,1:nlayer,i_vap_generic)=d_q(1:ngrid,1:nlayer)/ptimestep |
---|
430 | d_t(1:ngrid,1:nlayer)=d_t(1:ngrid,1:nlayer)/ptimestep |
---|
431 | do i=1,ngrid |
---|
432 | reevap_precip(i,i_vap_generic)=0. |
---|
433 | do k=1,nlayer |
---|
434 | reevap_precip(i,i_vap_generic)=reevap_precip(i,i_vap_generic)+dq_rain_generic_vap(i,k,i_vap_generic)*dmass(i,k) |
---|
435 | enddo |
---|
436 | enddo |
---|
437 | else |
---|
438 | dq_rain_generic_vap(1:ngrid,1:nlayer,i_vap_generic)=0.0 |
---|
439 | d_t(1:ngrid,1:nlayer)=0.0 |
---|
440 | endif |
---|
441 | dq_rain_generic_cld(1:ngrid,1:nlayer,i_ice_generic) = d_ql(1:ngrid,1:nlayer)/ptimestep |
---|
442 | |
---|
443 | end if ! if(call_ice_vap_generic) |
---|
444 | |
---|
445 | end do ! do iq=1,nq loop on tracers |
---|
446 | |
---|
447 | end subroutine rain_generic |
---|