source: LMDZ6/branches/Amaury_dev/libf/phylmd/lmdz_lscp_tools.F90 @ 5111

Last change on this file since 5111 was 5111, checked in by abarral, 4 months ago

Put abort_physic into a module
Remove -g option from makelmdz_fcm, since that option is linked to a header file that isn't included anywhere.
(lint) light lint on traversed files

File size: 29.6 KB
Line 
1MODULE lmdz_lscp_tools
2
3    IMPLICIT NONE
4
5CONTAINS
6
7!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
8SUBROUTINE FALLICE_VELOCITY(klon,iwc,temp,rho,pres,ptconv,velo)
9
10    ! Ref:
11    ! Stubenrauch, C. J., Bonazzola, M.,
12    ! Protopapadaki, S. E., & Musat, I. (2019).
13    ! New cloud system metrics to assess bulk
14    ! ice cloud schemes in a GCM. Journal of
15    ! Advances in Modeling Earth Systems, 11,
16    ! 3212–3234. https://doi.org/10.1029/2019MS001642
17   
18    use lmdz_lscp_ini, only: iflag_vice, ffallv_con, ffallv_lsc
19    use lmdz_lscp_ini, only: cice_velo, dice_velo
20
21    IMPLICIT NONE
22
23    INTEGER, INTENT(IN) :: klon
24    REAL, INTENT(IN), DIMENSION(klon) :: iwc       ! specific ice water content [kg/m3]
25    REAL, INTENT(IN), DIMENSION(klon) :: temp      ! temperature [K]
26    REAL, INTENT(IN), DIMENSION(klon) :: rho       ! dry air density [kg/m3]
27    REAL, INTENT(IN), DIMENSION(klon) :: pres      ! air pressure [Pa]
28    LOGICAL, INTENT(IN), DIMENSION(klon) :: ptconv    ! convective point  [-]
29
30    REAL, INTENT(OUT), DIMENSION(klon) :: velo    ! fallspeed velocity of crystals [m/s]
31
32
33    INTEGER i
34    REAL logvm,iwcg,tempc,phpa,fallv_tun
35    REAL m2ice, m2snow, vmice, vmsnow
36    REAL aice, bice, asnow, bsnow
37   
38
39    DO i=1,klon
40
41        IF (ptconv(i)) THEN
42            fallv_tun=ffallv_con
43        ELSE
44            fallv_tun=ffallv_lsc
45        ENDIF
46
47        tempc=temp(i)-273.15 ! celcius temp
48        iwcg=MAX(iwc(i)*1000.,1E-3) ! iwc in g/m3. We set a minimum value to prevent from division by 0
49        phpa=pres(i)/100.    ! pressure in hPa
50
51    IF (iflag_vice == 1) THEN
52        ! so-called 'empirical parameterization' in Stubenrauch et al. 2019
53        if (tempc >= -60.0) then
54
55            logvm= -0.0000414122*tempc*tempc*log(iwcg)-0.00538922*tempc*log(iwcg) &
56                    -0.0516344*log(iwcg)+0.00216078*tempc + 1.9714   
57            velo(i)=exp(logvm)
58        else
59            velo(i)=65.0*(iwcg**0.2)*(150./phpa)**0.15
60        endif
61       
62        velo(i)=fallv_tun*velo(i)/100.0 ! from cm/s to m/s
63
64    ELSE IF (iflag_vice == 2) THEN
65        ! so called  PSDM empirical coherent bulk ice scheme in Stubenrauch et al. 2019
66        aice=0.587
67        bice=2.45
68        asnow=0.0444
69        bsnow=2.1
70       
71        m2ice=((iwcg*0.001/aice)/(exp(13.6-bice*7.76+0.479*bice**2)* &
72                exp((-0.0361+bice*0.0151+0.00149*bice**2)*tempc)))   &
73                **(1./(0.807+bice*0.00581+0.0457*bice**2))
74
75        vmice=100.*1042.4*exp(13.6-(bice+1)*7.76+0.479*(bice+1.)**2)*exp((-0.0361+ &
76                 (bice+1.)*0.0151+0.00149*(bice+1.)**2)*tempc) &
77                 *(m2ice**(0.807+(bice+1.)*0.00581+0.0457*(bice+1.)**2))/(iwcg*0.001/aice)
78
79       
80        vmice=vmice*((1000./phpa)**0.2)
81     
82        m2snow=((iwcg*0.001/asnow)/(exp(13.6-bsnow*7.76+0.479*bsnow**2)* &
83               exp((-0.0361+bsnow*0.0151+0.00149*bsnow**2)*tempc)))         &
84               **(1./(0.807+bsnow*0.00581+0.0457*bsnow**2))
85
86
87        vmsnow=100.*14.3*exp(13.6-(bsnow+.416)*7.76+0.479*(bsnow+.416)**2)&
88                  *exp((-0.0361+(bsnow+.416)*0.0151+0.00149*(bsnow+.416)**2)*tempc)&
89                  *(m2snow**(0.807+(bsnow+.416)*0.00581+0.0457*(bsnow+.416)**2))/(iwcg*0.001/asnow)
90
91        vmsnow=vmsnow*((1000./phpa)**0.35)
92        velo(i)=fallv_tun*min(vmsnow,vmice)/100. ! to m/s
93
94    ELSE
95        ! By default, fallspeed velocity of ice crystals according to Heymsfield & Donner 1990
96        velo(i) = fallv_tun*cice_velo*((iwcg/1000.)**dice_velo)
97    ENDIF
98    ENDDO
99
100END SUBROUTINE FALLICE_VELOCITY
101!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
102
103!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
104SUBROUTINE ICEFRAC_LSCP(klon, temp, iflag_ice_thermo, distcltop, temp_cltop, icefrac, dicefracdT)
105!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
106 
107  ! Compute the ice fraction 1-xliq (see e.g.
108  ! Doutriaux-Boucher & Quaas 2004, section 2.2.)
109  ! as a function of temperature
110  ! see also Fig 3 of Madeleine et al. 2020, JAMES
111!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
112
113
114    USE print_control_mod, ONLY: lunout, prt_level
115    USE lmdz_lscp_ini, ONLY: t_glace_min, t_glace_max, exposant_glace, iflag_t_glace
116    USE lmdz_lscp_ini, ONLY: RTT, dist_liq, temp_nowater
117    USE lmdz_abort_physic, ONLY: abort_physic
118
119    IMPLICIT NONE
120
121
122    INTEGER, INTENT(IN)                 :: klon              ! number of horizontal grid points
123    REAL, INTENT(IN), DIMENSION(klon)   :: temp              ! temperature
124    REAL, INTENT(IN), DIMENSION(klon)   :: distcltop         ! distance to cloud top
125    REAL, INTENT(IN), DIMENSION(klon)   :: temp_cltop        ! temperature of cloud top
126    INTEGER, INTENT(IN)                 :: iflag_ice_thermo
127    REAL, INTENT(OUT), DIMENSION(klon)  :: icefrac
128    REAL, INTENT(OUT), DIMENSION(klon)  :: dicefracdT
129
130
131    INTEGER i
132    REAL    liqfrac_tmp, dicefrac_tmp
133    REAL    Dv, denomdep,beta,qsi,dqsidt
134    LOGICAL ice_thermo
135
136    CHARACTER (len = 20) :: modname = 'lscp_tools'
137    CHARACTER (len = 80) :: abort_message
138
139    IF ((iflag_t_glace<2)) THEN !.OR. (iflag_t_glace.GT.6)) THEN
140       abort_message = 'lscp cannot be used if iflag_t_glace<2 or >6'
141       CALL abort_physic(modname,abort_message,1)
142    ENDIF
143
144    IF (.NOT.((iflag_ice_thermo == 1).OR.(iflag_ice_thermo >= 3))) THEN
145       abort_message = 'lscp cannot be used without ice thermodynamics'
146       CALL abort_physic(modname,abort_message,1)
147    ENDIF
148
149
150    DO i=1,klon
151 
152        ! old function with sole dependence upon temperature
153        IF (iflag_t_glace == 2) THEN
154            liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
155            liqfrac_tmp = MIN(MAX(liqfrac_tmp,0.0),1.0)
156            icefrac(i) = (1.0-liqfrac_tmp)**exposant_glace
157            IF (icefrac(i) >0.) THEN
158                 dicefracdT(i)= exposant_glace * (icefrac(i)**(exposant_glace-1.)) &
159                           / (t_glace_min - t_glace_max)
160            ENDIF
161
162            IF ((icefrac(i)==0).OR.(icefrac(i)==1)) THEN
163                 dicefracdT(i)=0.
164            ENDIF
165
166        ENDIF
167
168        ! function of temperature used in CMIP6 physics
169        IF (iflag_t_glace == 3) THEN
170            liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
171            liqfrac_tmp = MIN(MAX(liqfrac_tmp,0.0),1.0)
172            icefrac(i) = 1.0-liqfrac_tmp**exposant_glace
173            IF ((icefrac(i) >0.) .AND. (liqfrac_tmp > 0.)) THEN
174                dicefracdT(i)= exposant_glace * ((liqfrac_tmp)**(exposant_glace-1.)) &
175                          / (t_glace_min - t_glace_max)
176            ELSE
177                dicefracdT(i)=0.
178            ENDIF
179        ENDIF
180
181        ! for iflag_t_glace .GE. 4, the liquid fraction depends upon temperature at cloud top
182        ! and then decreases with decreasing height
183
184        !with linear function of temperature at cloud top
185        IF (iflag_t_glace == 4) THEN
186                liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
187                liqfrac_tmp = MIN(MAX(liqfrac_tmp,0.0),1.0)
188                icefrac(i)    =  MAX(MIN(1.,1.0 - liqfrac_tmp*exp(-distcltop(i)/dist_liq)),0.)
189                dicefrac_tmp = - temp(i)/(t_glace_max-t_glace_min)
190                dicefracdT(i) = dicefrac_tmp*exp(-distcltop(i)/dist_liq)
191                IF ((liqfrac_tmp <=0) .OR. (liqfrac_tmp >= 1)) THEN
192                        dicefracdT(i) = 0.
193                ENDIF
194        ENDIF
195
196        ! with CMIP6 function of temperature at cloud top
197        IF ((iflag_t_glace == 5) .OR. (iflag_t_glace == 7)) THEN
198                liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
199                liqfrac_tmp =  MIN(MAX(liqfrac_tmp,0.0),1.0)
200                liqfrac_tmp = liqfrac_tmp**exposant_glace
201                icefrac(i)  =  MAX(MIN(1.,1.0 - liqfrac_tmp*exp(-distcltop(i)/dist_liq)),0.)
202                IF ((liqfrac_tmp <=0) .OR. (liqfrac_tmp >= 1)) THEN
203                        dicefracdT(i) = 0.
204                ELSE
205                        dicefracdT(i) = exposant_glace*((liqfrac_tmp)**(exposant_glace-1.))/(t_glace_min- t_glace_max) &
206                                        *exp(-distcltop(i)/dist_liq)
207                ENDIF
208        ENDIF
209
210        ! with modified function of temperature at cloud top
211        ! to get largere values around 260 K, works well with t_glace_min = 241K
212        IF (iflag_t_glace == 6) THEN
213                IF (temp(i) > t_glace_max) THEN
214                        liqfrac_tmp = 1.
215                ELSE
216                        liqfrac_tmp = -((temp(i)-t_glace_max) / (t_glace_max-t_glace_min))**2+1.
217                ENDIF
218                liqfrac_tmp  = MIN(MAX(liqfrac_tmp,0.0),1.0)
219                icefrac(i)   = MAX(MIN(1.,1.0 - liqfrac_tmp*exp(-distcltop(i)/dist_liq)),0.)       
220                IF ((liqfrac_tmp <=0) .OR. (liqfrac_tmp >= 1)) THEN
221                        dicefracdT(i) = 0.
222                ELSE
223                        dicefracdT(i) = 2*((temp(i)-t_glace_max) / (t_glace_max-t_glace_min))/(t_glace_max-t_glace_min) &
224                                  *exp(-distcltop(i)/dist_liq)
225                ENDIF
226        ENDIF
227
228        ! if temperature of cloud top <-40°C,
229        IF (iflag_t_glace >= 4) THEN
230                IF ((temp_cltop(i) <= temp_nowater) .AND. (temp(i) <= t_glace_max)) THEN
231                        icefrac(i) = 1.
232                        dicefracdT(i) = 0.
233                ENDIF
234        ENDIF
235     
236
237     ENDDO ! klon
238
239 
240END SUBROUTINE ICEFRAC_LSCP
241!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
242
243SUBROUTINE ICEFRAC_LSCP_TURB(klon, dtime, temp, pplay, paprsdn, paprsup, qice_ini, snowcld, qtot_incl, cldfra, tke, tke_dissip, qliq, qvap_cld, qice, icefrac, dicefracdT, cldfraliq, sigma2_icefracturb, mean_icefracturb)
244!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
245  ! Compute the liquid, ice and vapour content (+ice fraction) based
246  ! on turbulence (see Fields 2014, Furtado 2016)
247!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
248
249
250   USE lmdz_lscp_ini, ONLY: prt_level, lunout
251   USE lmdz_lscp_ini, ONLY: RCPD, RLSTT, RLVTT, RLMLT, RVTMP2, RTT, RD, RG, RV, RPI
252   USE lmdz_lscp_ini, ONLY: seuil_neb, temp_nowater
253   USE lmdz_lscp_ini, ONLY: tau_mixenv, lmix_mpc, naero5, gamma_snwretro, gamma_taud, capa_crystal
254   USE lmdz_lscp_ini, ONLY: eps
255
256   IMPLICIT NONE
257
258   INTEGER,   INTENT(IN)                           :: klon              !--number of horizontal grid points
259   REAL,      INTENT(IN)                           :: dtime             !--time step [s]
260
261   REAL,      INTENT(IN),       DIMENSION(klon)    :: temp              !--temperature
262   REAL,      INTENT(IN),       DIMENSION(klon)    :: pplay             !--pressure in the middle of the layer       [Pa]
263   REAL,      INTENT(IN),       DIMENSION(klon)    :: paprsdn           !--pressure at the bottom interface of the layer [Pa]
264   REAL,      INTENT(IN),       DIMENSION(klon)    :: paprsup           !--pressure at the top interface of the layer [Pa]
265   REAL,      INTENT(IN),       DIMENSION(klon)    :: qtot_incl         !--specific total cloud water content, in-cloud content [kg/kg]
266   REAL,      INTENT(IN),       DIMENSION(klon)    :: cldfra            !--cloud fraction in gridbox [-]
267   REAL,      INTENT(IN),       DIMENSION(klon)    :: tke               !--turbulent kinetic energy [m2/s2]
268   REAL,      INTENT(IN),       DIMENSION(klon)    :: tke_dissip        !--TKE dissipation [m2/s3]
269
270   REAL,      INTENT(IN),       DIMENSION(klon)    :: qice_ini          !--initial specific ice content gridbox-mean [kg/kg]
271   REAL,      INTENT(IN),       DIMENSION(klon)    :: snowcld
272   REAL,      INTENT(OUT),      DIMENSION(klon)    :: qliq              !--specific liquid content gridbox-mean [kg/kg]
273   REAL,      INTENT(OUT),      DIMENSION(klon)    :: qvap_cld          !--specific cloud vapor content, gridbox-mean [kg/kg]
274   REAL,      INTENT(OUT),      DIMENSION(klon)    :: qice              !--specific ice content gridbox-mean [kg/kg]
275   REAL,      INTENT(OUT),      DIMENSION(klon)    :: icefrac           !--fraction of ice in condensed water [-]
276   REAL,      INTENT(OUT),      DIMENSION(klon)    :: dicefracdT
277
278   REAL,      INTENT(OUT),      DIMENSION(klon)    :: cldfraliq         !--fraction of cldfra which is liquid only
279   REAL,      INTENT(OUT),      DIMENSION(klon)    :: sigma2_icefracturb     !--Temporary
280   REAL,      INTENT(OUT),      DIMENSION(klon)    :: mean_icefracturb      !--Temporary
281
282   REAL, DIMENSION(klon) :: qzero, qsatl, dqsatl, qsati, dqsati         !--specific humidity saturation values
283   INTEGER :: i
284
285   REAL :: qvap_incl, qice_incl, qliq_incl, qiceini_incl                !--In-cloud specific quantities [kg/kg]
286   REAL :: qsnowcld_incl
287   !REAL :: capa_crystal                                                 !--Capacitance of ice crystals  [-]
288   REAL :: water_vapor_diff                                             !--Water-vapour diffusion coefficient in air [m2/s] (function of T&P)
289   REAL :: air_thermal_conduct                                          !--Thermal conductivity of air [J/m/K/s] (function of T)
290   REAL :: C0                                                           !--Lagrangian structure function [-]
291   REAL :: tau_mixingenv
292   REAL :: tau_dissipturb
293   REAL :: invtau_phaserelax
294   REAL :: sigma2_pdf, mean_pdf
295   REAL :: ai, bi, B0
296   REAL :: sursat_iceliq
297   REAL :: sursat_env
298   REAL :: liqfra_max
299   REAL :: sursat_iceext
300   REAL :: nb_crystals                                                  !--number concentration of ice crystals [#/m3]
301   REAL :: moment1_PSD                                                  !--1st moment of ice PSD
302   REAL :: N0_PSD, lambda_PSD                                           !--parameters of the exponential PSD
303
304   REAL :: rho_ice                                                      !--ice density [kg/m3]
305   REAL :: cldfra1D
306   REAL :: deltaz, rho_air
307   REAL :: psati                                                        !--saturation vapor pressure wrt i [Pa]
308   
309   C0            = 10.                                                  !--value assumed in Field2014           
310   rho_ice       = 950.
311   sursat_iceext = -0.1
312   !capa_crystal  = 1. !r_ice                                       
313   qzero(:)      = 0.
314   cldfraliq(:)  = 0.
315   icefrac(:)    = 0.
316   dicefracdT(:) = 0.
317
318   sigma2_icefracturb(:) = 0.
319   mean_icefracturb(:)  = 0.
320
321   !--wrt liquid water
322   CALL calc_qsat_ecmwf(klon,temp(:),qzero(:),pplay(:),RTT,1,.FALSE.,qsatl(:),dqsatl(:))
323   !--wrt ice
324   CALL calc_qsat_ecmwf(klon,temp(:),qzero(:),pplay(:),RTT,2,.FALSE.,qsati(:),dqsati(:))
325
326
327    DO i=1,klon
328
329
330     rho_air  = pplay(i) / temp(i) / RD
331     !deltaz   = ( paprsdn(i) - paprsup(i) ) / RG / rho_air(i)
332     ! because cldfra is intent in, but can be locally modified due to test
333     cldfra1D = cldfra(i)
334     IF (cldfra(i) <= 0.) THEN
335        qvap_cld(i)   = 0.
336        qliq(i)       = 0.
337        qice(i)       = 0.
338        cldfraliq(i)  = 0.
339        icefrac(i)    = 0.
340        dicefracdT(i) = 0.
341
342     ! If there is a cloud
343     ELSE
344        IF (cldfra(i) >= 1.0) THEN
345           cldfra1D = 1.0
346        END IF
347       
348        ! T>0°C, no ice allowed
349        IF ( temp(i) >= RTT ) THEN
350           qvap_cld(i)   = qsatl(i) * cldfra1D
351           qliq(i)       = MAX(0.0,qtot_incl(i)-qsatl(i))  * cldfra1D
352           qice(i)       = 0.
353           cldfraliq(i)  = 1.
354           icefrac(i)    = 0.
355           dicefracdT(i) = 0.
356       
357        ! T<-38°C, no liquid allowed
358        ELSE IF ( temp(i) <= temp_nowater) THEN
359           qvap_cld(i)   = qsati(i) * cldfra1D
360           qliq(i)       = 0.
361           qice(i)       = MAX(0.0,qtot_incl(i)-qsati(i)) * cldfra1D
362           cldfraliq(i)  = 0.
363           icefrac(i)    = 1.
364           dicefracdT(i) = 0.
365
366        ! MPC temperature
367        ELSE
368           ! Not enough TKE     
369           IF ( tke_dissip(i) <= eps )  THEN
370              qvap_cld(i)   = qsati(i) * cldfra1D
371              qliq(i)       = 0.
372              qice(i)       = MAX(0.,qtot_incl(i)-qsati(i)) * cldfra1D   
373              cldfraliq(i)  = 0.
374              icefrac(i)    = 1.
375              dicefracdT(i) = 0.
376           
377           ! Enough TKE   
378           ELSE   
379              !---------------------------------------------------------
380              !--               ICE SUPERSATURATION PDF   
381              !---------------------------------------------------------
382              !--If -38°C< T <0°C and there is enough turbulence,
383              !--we compute the cloud liquid properties with a Gaussian PDF
384              !--of ice supersaturation F(Si) (Field2014, Furtado2016).
385              !--Parameters of the PDF are function of turbulence and
386              !--microphysics/existing ice.
387
388              sursat_iceliq = qsatl(i)/qsati(i) - 1.
389              psati         = qsati(i) * pplay(i) / (RD/RV)
390
391              !-------------- MICROPHYSICAL TERMS --------------
392              !--We assume an exponential ice PSD whose parameters
393              !--are computed following Morrison&Gettelman 2008
394              !--Ice number density is assumed equals to INP density
395              !--which is a function of temperature (DeMott 2010) 
396              !--bi and B0 are microphysical function characterizing
397              !--vapor/ice interactions
398              !--tau_phase_relax is the typical time of vapor deposition
399              !--onto ice crystals
400             
401              qiceini_incl  = qice_ini(i) / cldfra1D
402              qsnowcld_incl = snowcld(i) * RG * dtime / ( paprsdn(i) - paprsup(i) ) / cldfra1D
403              sursat_env    = max(0., (qtot_incl(i) - qiceini_incl)/qsati(i) - 1.)
404              IF ( qiceini_incl > eps ) THEN
405                nb_crystals = 1.e3 * 5.94e-5 * ( RTT - temp(i) )**3.33 * naero5**(0.0264*(RTT-temp(i))+0.0033)
406                lambda_PSD  = ( (RPI*rho_ice*nb_crystals*24.) / (6.*(qiceini_incl + gamma_snwretro * qsnowcld_incl)) ) ** (1./3.)
407                N0_PSD      = nb_crystals * lambda_PSD
408                moment1_PSD = N0_PSD/2./lambda_PSD**2
409              ELSE
410                moment1_PSD = 0.
411              ENDIF
412
413              !--Formulae for air thermal conductivity and water vapor diffusivity
414              !--comes respectively from Beard and Pruppacher (1971)
415              !--and  Hall and Pruppacher (1976)
416
417              air_thermal_conduct = ( 5.69 + 0.017 * ( temp(i) - RTT ) ) * 1.e-3 * 4.184
418              water_vapor_diff    = 2.11*1e-5 * ( temp(i) / RTT )**1.94 * ( 101325 / pplay(i) )
419             
420              bi = 1./((qsati(i)+qsatl(i))/2.) + RLSTT**2 / RCPD / RV / temp(i)**2
421              B0 = 4. * RPI * capa_crystal * 1. / (  RLSTT**2 / air_thermal_conduct / RV / temp(i)**2  &
422                                                  +  RV * temp(i) / psati / water_vapor_diff  )
423
424              invtau_phaserelax  = (bi * B0 * moment1_PSD )
425
426!             Old way of estimating moment1 : spherical crystals + monodisperse PSD             
427!             nb_crystals = rho_air * qiceini_incl / ( 4. / 3. * RPI * r_ice**3. * rho_ice )
428!             moment1_PSD = nb_crystals * r_ice
429
430              !----------------- TURBULENT SOURCE/SINK TERMS -----------------
431              !--Tau_mixingenv is the time needed to homogeneize the parcel
432              !--with its environment by turbulent diffusion over the parcel
433              !--length scale
434              !--if lmix_mpc <0, tau_mixigenv value is prescribed
435              !--else tau_mixigenv value is derived from tke_dissip and lmix_mpc
436              !--Tau_dissipturb is the time needed turbulence to decay due to
437              !--viscosity
438
439              ai = RG / RD / temp(i) * ( RD * RLSTT / RCPD / RV / temp(i) - 1. )
440              IF ( lmix_mpc > 0 ) THEN
441                 tau_mixingenv = ( lmix_mpc**2. / tke_dissip(i) )**(1./3.)
442              ELSE
443                 tau_mixingenv = tau_mixenv
444              ENDIF
445             
446              tau_dissipturb = gamma_taud * 2. * 2./3. * tke(i) / tke_dissip(i) / C0
447
448              !--------------------- PDF COMPUTATIONS ---------------------
449              !--Formulae for sigma2_pdf (variance), mean of PDF in Furtado2016
450              !--cloud liquid fraction and in-cloud liquid content are given
451              !--by integrating resp. F(Si) and Si*F(Si)
452              !--Liquid is limited by the available water vapor trough a
453              !--maximal liquid fraction
454
455              liqfra_max = MAX(0., (MIN (1.,( qtot_incl(i) - qiceini_incl - qsati(i) * (1 + sursat_iceext ) ) / ( qsatl(i) - qsati(i) ) ) ) )
456              sigma2_pdf = 1./2. * ( ai**2 ) *  2./3. * tke(i) * tau_dissipturb / ( invtau_phaserelax + 1./tau_mixingenv )
457              mean_pdf   = sursat_env * 1./tau_mixingenv / ( invtau_phaserelax + 1./tau_mixingenv )
458              cldfraliq(i) = 0.5 * (1. - erf( ( sursat_iceliq - mean_pdf) / (SQRT(2.* sigma2_pdf) ) ) )
459              !IF (cldfraliq(i) .GT. liqfra_max) THEN
460              !    cldfraliq(i) = liqfra_max
461              !ENDIF
462             
463              qliq_incl = qsati(i) * SQRT(sigma2_pdf) / SQRT(2.*RPI) * EXP( -1.*(sursat_iceliq - mean_pdf)**2. / (2.*sigma2_pdf) )  &
464                        - qsati(i) * cldfraliq(i) * (sursat_iceliq - mean_pdf )
465             
466              sigma2_icefracturb(i)= sigma2_pdf
467              mean_icefracturb(i)  = mean_pdf     
468              !------------ ICE AMOUNT AND WATER CONSERVATION  ------------
469
470              IF ( (qliq_incl <= eps) .OR. (cldfraliq(i) <= eps) ) THEN
471                  qliq_incl    = 0.
472                  cldfraliq(i) = 0.
473              END IF
474             
475              !--Choice for in-cloud vapor :
476              !--1.Weighted mean between qvap in MPC parts and in ice-only parts
477              !--2.Always at ice saturation
478              qvap_incl = MAX(qsati(i), ( 1. - cldfraliq(i) ) * (sursat_iceext + 1.) * qsati(i) + cldfraliq(i) * qsatl(i) )
479               
480              IF ( qvap_incl  >= qtot_incl(i) ) THEN
481                 qvap_incl = qsati(i)
482                 qliq_incl = qtot_incl(i) - qvap_incl
483                 qice_incl = 0.
484
485              ELSEIF ( (qvap_incl + qliq_incl) >= qtot_incl(i) ) THEN
486                 qliq_incl = MAX(0.0,qtot_incl(i) - qvap_incl)
487                 qice_incl = 0.
488              ELSE
489                 qice_incl = qtot_incl(i) - qvap_incl - qliq_incl
490              END IF
491
492              qvap_cld(i)   = qvap_incl * cldfra1D
493              qliq(i)       = qliq_incl * cldfra1D
494              qice(i)       = qice_incl * cldfra1D
495              icefrac(i)    = qice(i) / ( qice(i) + qliq(i) )
496              dicefracdT(i) = 0.
497              !PRINT*,'MPC turb'
498
499           END IF ! Enough TKE
500
501        END IF ! ! MPC temperature
502
503     END IF ! cldfra
504   
505   ENDDO ! klon
506END SUBROUTINE ICEFRAC_LSCP_TURB
507!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
508
509
510SUBROUTINE CALC_QSAT_ECMWF(klon,temp,qtot,pressure,tref,phase,flagth,qs,dqs)
511!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
512    ! Calculate qsat following ECMWF method
513!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
514
515
516    IMPLICIT NONE
517
518    include "YOMCST.h"
519    include "YOETHF.h"
520    include "FCTTRE.h"
521
522    INTEGER, INTENT(IN) :: klon  ! number of horizontal grid points
523    REAL, INTENT(IN), DIMENSION(klon) :: temp     ! temperature in K
524    REAL, INTENT(IN), DIMENSION(klon) :: qtot     ! total specific water in kg/kg
525    REAL, INTENT(IN), DIMENSION(klon) :: pressure ! pressure in Pa
526    REAL, INTENT(IN)                  :: tref     ! reference temperature in K
527    LOGICAL, INTENT(IN) :: flagth     ! flag for qsat calculation for thermals
528    INTEGER, INTENT(IN) :: phase
529    ! phase: 0=depend on temperature sign (temp>tref -> liquid, temp<tref, solid)
530    !        1=liquid
531    !        2=solid
532
533    REAL, INTENT(OUT), DIMENSION(klon) :: qs      ! saturation specific humidity [kg/kg]
534    REAL, INTENT(OUT), DIMENSION(klon) :: dqs     ! derivation of saturation specific humidity wrt T
535
536    REAL delta, cor, cvm5
537    INTEGER i
538
539    DO i=1,klon
540
541    IF (phase == 1) THEN
542        delta=0.
543    ELSEIF (phase == 2) THEN
544        delta=1.
545    ELSE
546        delta=MAX(0.,SIGN(1.,tref-temp(i)))
547    ENDIF
548
549    IF (flagth) THEN
550    cvm5=R5LES*(1.-delta) + R5IES*delta
551    ELSE
552    cvm5 = R5LES*RLVTT*(1.-delta) + R5IES*RLSTT*delta
553    cvm5 = cvm5 /RCPD/(1.0+RVTMP2*(qtot(i)))
554    ENDIF
555
556    qs(i)= R2ES*FOEEW(temp(i),delta)/pressure(i)
557    qs(i)=MIN(0.5,qs(i))
558    cor=1./(1.-RETV*qs(i))
559    qs(i)=qs(i)*cor
560    dqs(i)= FOEDE(temp(i),delta,cvm5,qs(i),cor)
561
562    END DO
563
564END SUBROUTINE CALC_QSAT_ECMWF
565!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
566
567
568!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
569SUBROUTINE CALC_GAMMASAT(klon,temp,qtot,pressure,gammasat,dgammasatdt)
570
571!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
572    ! programme that calculates the gammasat parameter that determines the
573    ! homogeneous condensation thresholds for cold (<0oC) clouds
574    ! condensation at q>gammasat*qsat
575    ! Etienne Vignon, March 2021
576!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
577
578    use lmdz_lscp_ini, only: iflag_gammasat, t_glace_min, RTT
579
580    IMPLICIT NONE
581
582
583    INTEGER, INTENT(IN) :: klon                       ! number of horizontal grid points
584    REAL, INTENT(IN), DIMENSION(klon) :: temp         ! temperature in K
585    REAL, INTENT(IN), DIMENSION(klon) :: qtot         ! total specific water in kg/kg
586
587    REAL, INTENT(IN), DIMENSION(klon) :: pressure     ! pressure in Pa
588
589    REAL, INTENT(OUT), DIMENSION(klon) :: gammasat    ! coefficient to multiply qsat with to calculate saturation
590    REAL, INTENT(OUT), DIMENSION(klon) :: dgammasatdt ! derivative of gammasat wrt temperature
591
592    REAL, DIMENSION(klon) ::  qsi,qsl,dqsl,dqsi
593    REAL  fcirrus, fac
594    REAL, PARAMETER :: acirrus=2.349
595    REAL, PARAMETER :: bcirrus=259.0
596
597    INTEGER i
598   
599        CALL CALC_QSAT_ECMWF(klon,temp,qtot,pressure,RTT,1,.FALSE.,qsl,dqsl)
600        CALL CALC_QSAT_ECMWF(klon,temp,qtot,pressure,RTT,2,.FALSE.,qsi,dqsi)
601
602    DO i=1,klon
603
604        IF (temp(i) >= RTT) THEN
605            ! warm clouds: condensation at saturation wrt liquid
606            gammasat(i)=1.
607            dgammasatdt(i)=0.
608
609        ELSEIF ((temp(i) < RTT) .AND. (temp(i) > t_glace_min)) THEN
610           
611            IF (iflag_gammasat >= 2) THEN
612               gammasat(i)=qsl(i)/qsi(i)
613               dgammasatdt(i)=(dqsl(i)*qsi(i)-dqsi(i)*qsl(i))/qsi(i)/qsi(i)
614            ELSE
615               gammasat(i)=1.
616               dgammasatdt(i)=0.
617            ENDIF
618
619        ELSE
620
621            IF (iflag_gammasat >=1) THEN
622            ! homogeneous freezing of aerosols, according to
623            ! Koop, 2000 and Karcher 2008, QJRMS
624            ! 'Cirrus regime'
625               fcirrus=acirrus-temp(i)/bcirrus
626               IF (fcirrus > qsl(i)/qsi(i)) THEN
627                  gammasat(i)=qsl(i)/qsi(i)
628                  dgammasatdt(i)=(dqsl(i)*qsi(i)-dqsi(i)*qsl(i))/qsi(i)/qsi(i)
629               ELSE
630                  gammasat(i)=fcirrus
631                  dgammasatdt(i)=-1.0/bcirrus
632               ENDIF
633           
634            ELSE
635
636               gammasat(i)=1.
637               dgammasatdt(i)=0.
638
639            ENDIF
640
641        ENDIF
642   
643    END DO
644
645
646END SUBROUTINE CALC_GAMMASAT
647!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
648
649
650!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
651SUBROUTINE DISTANCE_TO_CLOUD_TOP(klon,klev,k,temp,pplay,paprs,rneb,distcltop1D,temp_cltop)
652!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
653   
654   USE lmdz_lscp_ini, ONLY: rd,rg,tresh_cl
655
656   IMPLICIT NONE
657   
658   INTEGER, INTENT(IN) :: klon,klev                !number of horizontal and vertical grid points
659   INTEGER, INTENT(IN) :: k                        ! vertical index
660   REAL, INTENT(IN), DIMENSION(klon,klev) :: temp  ! temperature in K
661   REAL, INTENT(IN), DIMENSION(klon,klev) :: pplay ! pressure middle layer in Pa
662   REAL, INTENT(IN), DIMENSION(klon,klev+1) :: paprs ! pressure interfaces in Pa
663   REAL, INTENT(IN), DIMENSION(klon,klev) :: rneb  ! cloud fraction
664
665   REAL, INTENT(OUT), DIMENSION(klon) :: distcltop1D  ! distance from cloud top
666   REAL, INTENT(OUT), DIMENSION(klon) :: temp_cltop     ! temperature of cloud top
667   
668   REAL dzlay(klon,klev)
669   REAL zlay(klon,klev)
670   REAL dzinterf
671   INTEGER i,k_top, kvert
672   LOGICAL bool_cl
673
674
675   DO i=1,klon
676         ! Initialization height middle of first layer
677          dzlay(i,1) = Rd * temp(i,1) / rg * log(paprs(i,1)/paprs(i,2))
678          zlay(i,1) = dzlay(i,1)/2
679
680          DO kvert=2,klev
681                 IF (kvert==klev) THEN
682                       dzlay(i,kvert) = 2*(rd * temp(i,kvert) / rg * log(paprs(i,kvert)/pplay(i,kvert)))
683                 ELSE
684                       dzlay(i,kvert) = rd * temp(i,kvert) / rg * log(paprs(i,kvert)/paprs(i,kvert+1))
685                 ENDIF
686                       dzinterf       = rd * temp(i,kvert) / rg * log(pplay(i,kvert-1)/pplay(i,kvert))
687                       zlay(i,kvert)  = zlay(i,kvert-1) + dzinterf
688           ENDDO
689   ENDDO
690   
691   DO i=1,klon
692          k_top = k
693          IF (rneb(i,k) <= tresh_cl) THEN
694                 bool_cl = .FALSE.
695          ELSE
696                 bool_cl = .TRUE.
697          ENDIF
698
699          DO WHILE ((bool_cl) .AND. (k_top <= klev))
700          ! find cloud top
701                IF (rneb(i,k_top) > tresh_cl) THEN
702                      k_top = k_top + 1
703                ELSE
704                      bool_cl = .FALSE.
705                      k_top   = k_top - 1
706                ENDIF
707          ENDDO
708          k_top=min(k_top,klev)
709
710          !dist to top is dist between current layer and layer of cloud top (from middle to middle) + dist middle to
711          !interf for layer of cloud top
712          distcltop1D(i) = zlay(i,k_top) - zlay(i,k) + dzlay(i,k_top)/2
713          temp_cltop(i)  = temp(i,k_top)
714   ENDDO ! klon
715
716END SUBROUTINE DISTANCE_TO_CLOUD_TOP
717!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
718
719END MODULE lmdz_lscp_tools
720
721
Note: See TracBrowser for help on using the repository browser.