source: LMDZ6/trunk/libf/phylmd/lmdz_lscp_tools.f90 @ 5649

Last change on this file since 5649 was 5633, checked in by evignon, 2 months ago

small changes related to mixed phase clouds treatment
Etienne and Lea

File size: 39.0 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 .EQ. 1) THEN
52        ! so-called 'empirical parameterization' in Stubenrauch et al. 2019
53        if (tempc .GE. -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 .EQ. 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
118    IMPLICIT NONE
119
120
121    INTEGER, INTENT(IN)                 :: klon              ! number of horizontal grid points
122    REAL, INTENT(IN), DIMENSION(klon)   :: temp              ! temperature
123    REAL, INTENT(IN), DIMENSION(klon)   :: distcltop         ! distance to cloud top
124    REAL, INTENT(IN), DIMENSION(klon)   :: temp_cltop        ! temperature of cloud top
125    INTEGER, INTENT(IN)                 :: iflag_ice_thermo
126    REAL, INTENT(OUT), DIMENSION(klon)  :: icefrac
127    REAL, INTENT(OUT), DIMENSION(klon)  :: dicefracdT
128
129
130    INTEGER i
131    REAL    liqfrac_tmp, dicefrac_tmp
132    REAL    Dv, denomdep,beta,qsi,dqsidt
133    LOGICAL ice_thermo
134
135    CHARACTER (len = 20) :: modname = 'lscp_tools'
136    CHARACTER (len = 80) :: abort_message
137
138    IF ((iflag_t_glace.LT.2)) THEN !.OR. (iflag_t_glace.GT.6)) THEN
139       abort_message = 'lscp cannot be used if iflag_t_glace<2 or >6'
140       CALL abort_physic(modname,abort_message,1)
141    ENDIF
142
143    IF (.NOT.((iflag_ice_thermo .EQ. 1).OR.(iflag_ice_thermo .GE. 3))) THEN
144       abort_message = 'lscp cannot be used without ice thermodynamics'
145       CALL abort_physic(modname,abort_message,1)
146    ENDIF
147
148
149    DO i=1,klon
150 
151        ! old function with sole dependence upon temperature
152        IF (iflag_t_glace .EQ. 2) THEN
153            liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
154            liqfrac_tmp = MIN(MAX(liqfrac_tmp,0.0),1.0)
155            icefrac(i) = (1.0-liqfrac_tmp)**exposant_glace
156            IF (icefrac(i) .GT.0.) THEN
157                 dicefracdT(i)= exposant_glace * (icefrac(i)**(exposant_glace-1.)) &
158                           / (t_glace_min - t_glace_max)
159            ENDIF
160
161            IF ((icefrac(i).EQ.0).OR.(icefrac(i).EQ.1)) THEN
162                 dicefracdT(i)=0.
163            ENDIF
164
165        ENDIF
166
167        ! function of temperature used in CMIP6 physics
168        IF (iflag_t_glace .EQ. 3) THEN
169            liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
170            liqfrac_tmp = MIN(MAX(liqfrac_tmp,0.0),1.0)
171            icefrac(i) = 1.0-liqfrac_tmp**exposant_glace
172            IF ((icefrac(i) .GT.0.) .AND. (liqfrac_tmp .GT. 0.)) THEN
173                dicefracdT(i)= exposant_glace * ((liqfrac_tmp)**(exposant_glace-1.)) &
174                          / (t_glace_min - t_glace_max)
175            ELSE
176                dicefracdT(i)=0.
177            ENDIF
178        ENDIF
179
180        ! for iflag_t_glace .GE. 4, the liquid fraction depends upon temperature at cloud top
181        ! and then decreases with decreasing height
182
183        !with linear function of temperature at cloud top
184        IF (iflag_t_glace .EQ. 4) THEN 
185                liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
186                liqfrac_tmp = MIN(MAX(liqfrac_tmp,0.0),1.0)
187                icefrac(i)    =  MAX(MIN(1.,1.0 - liqfrac_tmp*exp(-distcltop(i)/dist_liq)),0.)
188                dicefrac_tmp = - temp(i)/(t_glace_max-t_glace_min)
189                dicefracdT(i) = dicefrac_tmp*exp(-distcltop(i)/dist_liq)
190                IF ((liqfrac_tmp .LE.0) .OR. (liqfrac_tmp .GE. 1)) THEN
191                        dicefracdT(i) = 0.
192                ENDIF
193        ENDIF
194
195        ! with CMIP6 function of temperature at cloud top
196        IF ((iflag_t_glace .EQ. 5) .OR. (iflag_t_glace .EQ. 7)) THEN
197                liqfrac_tmp = (temp(i)-t_glace_min) / (t_glace_max-t_glace_min)
198                liqfrac_tmp =  MIN(MAX(liqfrac_tmp,0.0),1.0)
199                liqfrac_tmp = liqfrac_tmp**exposant_glace
200                icefrac(i)  =  MAX(MIN(1.,1.0 - liqfrac_tmp*exp(-distcltop(i)/dist_liq)),0.)
201                IF ((liqfrac_tmp .LE.0) .OR. (liqfrac_tmp .GE. 1)) THEN
202                        dicefracdT(i) = 0.
203                ELSE
204                        dicefracdT(i) = exposant_glace*((liqfrac_tmp)**(exposant_glace-1.))/(t_glace_min- t_glace_max) &
205                                        *exp(-distcltop(i)/dist_liq)
206                ENDIF
207        ENDIF
208
209        ! with modified function of temperature at cloud top
210        ! to get largere values around 260 K, works well with t_glace_min = 241K
211        IF (iflag_t_glace .EQ. 6) THEN
212                IF (temp(i) .GT. t_glace_max) THEN
213                        liqfrac_tmp = 1.
214                ELSE
215                        liqfrac_tmp = -((temp(i)-t_glace_max) / (t_glace_max-t_glace_min))**2+1.
216                ENDIF
217                liqfrac_tmp  = MIN(MAX(liqfrac_tmp,0.0),1.0)
218                icefrac(i)   = MAX(MIN(1.,1.0 - liqfrac_tmp*exp(-distcltop(i)/dist_liq)),0.)       
219                IF ((liqfrac_tmp .LE.0) .OR. (liqfrac_tmp .GE. 1)) THEN
220                        dicefracdT(i) = 0.
221                ELSE
222                        dicefracdT(i) = 2*((temp(i)-t_glace_max) / (t_glace_max-t_glace_min))/(t_glace_max-t_glace_min) &
223                                  *exp(-distcltop(i)/dist_liq)
224                ENDIF
225        ENDIF
226
227        ! if temperature or temperature of cloud top <-40°C,
228        IF (iflag_t_glace .GE. 4) THEN
229                IF ((temp_cltop(i) .LE. temp_nowater) .AND. (temp(i) .LE. t_glace_max)) THEN
230                        icefrac(i) = 1.
231                        dicefracdT(i) = 0.
232                ENDIF
233        ENDIF
234     
235
236     ENDDO ! klon
237     RETURN
238 
239END SUBROUTINE ICEFRAC_LSCP
240!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
241
242
243SUBROUTINE ICEFRAC_LSCP_TURB(klon, dtime, pticefracturb, temp, pplay, paprsdn, paprsup, wvel, qice_ini, snowcld, snowfracld, qtot_incl, cldfra, tke,   &
244                             tke_dissip, sursat_e, invtau_e, qliq, qvap_cld, qice, icefrac, dicefracdT, cldfraliq, sigma2_icefracturb, mean_icefracturb)
245!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
246  ! Compute the liquid, ice and vapour content (+ice fraction) based
247  ! on turbulence (see Fields 2014, Furtado 2016, Raillard 2025)
248  ! L.Raillard (23/09/24)
249  ! E.Vignon (03/2025) : additional elements for treatment of convective
250  !                      boundary layer clouds
251!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
252
253
254   USE lmdz_lscp_ini, ONLY : prt_level, lunout
255   USE lmdz_lscp_ini, ONLY : RCPD, RLSTT, RLVTT, RLMLT, RVTMP2, RTT, RD, RG, RV, RPI
256   USE lmdz_lscp_ini, ONLY : seuil_neb, temp_nowater
257   USE lmdz_lscp_ini, ONLY : naero5, gamma_snwretro, gamma_taud, capa_crystal, rho_ice
258   USE lmdz_lscp_ini, ONLY : eps, snow_fallspeed
259
260   IMPLICIT NONE
261
262   INTEGER,   INTENT(IN)                           :: klon                !--number of horizontal grid points
263   REAL,      INTENT(IN)                           :: dtime               !--time step [s]
264   LOGICAL,   INTENT(IN),       DIMENSION(klon)    :: pticefracturb       !--grid points concerned by this routine 
265   REAL,      INTENT(IN),       DIMENSION(klon)    :: temp                !--temperature
266   REAL,      INTENT(IN),       DIMENSION(klon)    :: pplay               !--pressure in the middle of the layer           [Pa]
267   REAL,      INTENT(IN),       DIMENSION(klon)    :: paprsdn             !--pressure at the bottom interface of the layer [Pa]
268   REAL,      INTENT(IN),       DIMENSION(klon)    :: paprsup             !--pressure at the top interface of the layer    [Pa]
269   REAL,      INTENT(IN),       DIMENSION(klon)    :: wvel                !--vertical velocity                             [m/s]
270   REAL,      INTENT(IN),       DIMENSION(klon)    :: qtot_incl           !--specific total cloud water in-cloud content   [kg/kg]
271   REAL,      INTENT(IN),       DIMENSION(klon)    :: cldfra              !--cloud fraction in gridbox                     [-]
272   REAL,      INTENT(IN),       DIMENSION(klon)    :: tke                 !--turbulent kinetic energy                      [m2/s2]
273   REAL,      INTENT(IN),       DIMENSION(klon)    :: tke_dissip          !--TKE dissipation                               [m2/s3]
274
275   REAL,      INTENT(IN),       DIMENSION(klon)    :: qice_ini            !--initial specific ice content gridbox-mean     [kg/kg]
276   REAL,      INTENT(IN),       DIMENSION(klon)    :: snowcld             !--in-cloud snowfall flux                        [kg/m2/s]
277   REAL,      INTENT(IN),       DIMENSION(klon)    :: snowfracld          !--cloudy precip fraction                        [-]
278   REAL,      INTENT(IN),       DIMENSION(klon)    :: sursat_e            !--environment supersaturation                   [-]
279   REAL,      INTENT(IN),       DIMENSION(klon)    :: invtau_e            !--inverse time-scale of mixing with environment [s-1]
280   REAL,      INTENT(OUT),      DIMENSION(klon)    :: qliq                !--specific liquid content gridbox-mean          [kg/kg]
281   REAL,      INTENT(OUT),      DIMENSION(klon)    :: qvap_cld            !--specific cloud vapor content, gridbox-mean    [kg/kg]
282   REAL,      INTENT(OUT),      DIMENSION(klon)    :: qice                !--specific ice content gridbox-mean             [kg/kg]
283
284   REAL,      INTENT(INOUT),    DIMENSION(klon)    :: icefrac             !--fraction of ice in condensed water            [-]
285   REAL,      INTENT(INOUT),    DIMENSION(klon)    :: dicefracdT
286
287   REAL,      INTENT(OUT),      DIMENSION(klon)    :: cldfraliq           !--fraction of cldfra where liquid               [-]
288   REAL,      INTENT(OUT),      DIMENSION(klon)    :: sigma2_icefracturb  !--Sigma2 of the ice supersaturation PDF         [-]
289   REAL,      INTENT(OUT),      DIMENSION(klon)    :: mean_icefracturb    !--Mean of the ice supersaturation PDF           [-]
290
291   REAL, DIMENSION(klon) :: qzero, qsatl, dqsatl, qsati, dqsati           !--specific humidity saturation values
292   INTEGER :: i
293
294   REAL :: qvap_incl, qice_incl, qliq_incl, qiceini_incl                  !--In-cloud specific quantities                  [kg/kg]
295   REAL :: water_vapor_diff                                               !--Water-vapour diffusion coeff in air (f(T,P))  [m2/s]
296   REAL :: air_thermal_conduct                                            !--Thermal conductivity of air (f(T))            [J/m/K/s]
297   REAL :: C0                                                             !--Lagrangian structure function                 [-]
298   REAL :: tau_dissipturb
299   REAL :: invtau_phaserelax
300   REAL :: sigma2_pdf
301   REAL :: ai, bi, B0
302   REAL :: sursat_iceliq
303   REAL :: sursat_equ
304   REAL :: liqfra_max
305   REAL :: sursat_iceext
306   REAL :: nb_crystals                                                    !--number concentration of ice crystals          [#/m3]
307   REAL :: moment1_PSD                                                    !--1st moment of ice PSD
308   REAL :: N0_PSD, lambda_PSD                                             !--parameters of the exponential PSD
309
310   REAL :: cldfra1D
311   REAL :: rho_air
312   REAL :: psati                                                          !--saturation vapor pressure wrt ice             [Pa]
313
314                                                                       
315    REAL :: tempvig1, tempvig2
316
317   tempvig1    = -21.06 + RTT
318   tempvig2    = -30.35 + RTT
319   C0            = 10.                                                    !--value assumed in Field2014           
320   sursat_iceext = -0.1
321   qzero(:)      = 0.
322   cldfraliq(:)  = 0.
323   qliq(:)       = 0.
324   qice(:)       = 0.
325   qvap_cld(:)   = 0.
326   sigma2_icefracturb(:) = 0.
327   mean_icefracturb(:)   = 0.
328
329   !--wrt liquid
330   CALL calc_qsat_ecmwf(klon,temp(:),qzero(:),pplay(:),RTT,1,.false.,qsatl(:),dqsatl(:))
331   !--wrt ice
332   CALL calc_qsat_ecmwf(klon,temp(:),qzero(:),pplay(:),RTT,2,.false.,qsati(:),dqsati(:))
333
334
335    DO i=1,klon
336     rho_air  = pplay(i) / temp(i) / RD
337     ! because cldfra is intent in, but can be locally modified due to test
338     cldfra1D = cldfra(i)
339     ! activate param for concerned grid points and for cloudy conditions
340     IF ((pticefracturb(i)) .AND. (cldfra(i) .GT. 0.)) THEN
341        IF (cldfra(i) .GE. 1.0) THEN
342           cldfra1D = 1.0
343        END IF
344       
345        ! T>0°C, no ice allowed
346        IF ( temp(i) .GE. RTT ) THEN
347           qvap_cld(i)   = qsatl(i) * cldfra1D
348           qliq(i)       = MAX(0.0,qtot_incl(i)-qsatl(i))  * cldfra1D
349           qice(i)       = 0.
350           cldfraliq(i)  = 1.
351           icefrac(i)    = 0.
352           dicefracdT(i) = 0.
353       
354        ! T<-38°C, no liquid allowed
355        ELSE IF ( temp(i) .LE. temp_nowater) THEN
356           qvap_cld(i)   = qsati(i) * cldfra1D
357           qliq(i)       = 0.
358           qice(i)       = MAX(0.0,qtot_incl(i)-qsati(i)) * cldfra1D
359           cldfraliq(i)  = 0.
360           icefrac(i)    = 1.
361           dicefracdT(i) = 0.
362
363
364        !---------------------------------------------------------
365        !--             MIXED PHASE TEMPERATURE REGIME         
366        !---------------------------------------------------------
367        !--In the mixed phase regime (-38°C< T <0°C) we distinguish
368        !--3 possible subcases.
369        !--1.  No pre-existing ice 
370        !--2A. Pre-existing ice and no turbulence
371        !--2B. Pre-existing ice and turbulence
372        ELSE
373
374           ! gamma_snwretro controls the contribution of snowflakes to the negative feedback
375           ! note that for reasons related to inetarctions with the condensation iteration in lscp_main
376           ! we consider here the mean snowflake concentration in the mesh (not the in-cloud concentration)
377           ! when poprecip is active, it will be worth testing considering the incloud fraction, dividing
378           ! by snowfracld     
379           ! qiceini_incl  = qice_ini(i) / cldfra1D + &
380           !              gamma_snwretro * snowcld(i) * RG * dtime / ( paprsdn(i) - paprsup(i) )
381           ! assuming constant snowfall velocity
382           qiceini_incl  = qice_ini(i) / cldfra1D + gamma_snwretro * snowcld(i) / pplay(i) * RD * temp(i) / snow_fallspeed
383
384           !--1. No preexisting ice and no mixing with environment: if vertical motion, fully liquid
385           !--cloud else fully iced cloud
386           IF ( (qiceini_incl .LT. eps) .AND. (invtau_e(i) .LT. eps) ) THEN
387              IF ( (wvel(i) .GT. eps) .OR. (tke(i) .GT. eps) ) THEN
388                 qvap_cld(i)   = qsatl(i) * cldfra1D
389                 qliq(i)       = MAX(0.,qtot_incl(i)-qsatl(i)) * cldfra1D
390                 qice(i)       = 0.
391                 cldfraliq(i)  = 1.
392                 icefrac(i)    = 0.
393                 dicefracdT(i) = 0.
394              ELSE
395                 qvap_cld(i)   = qsati(i) * cldfra1D
396                 qliq(i)       = 0.
397                 qice(i)       = MAX(0.,qtot_incl(i)-qsati(i)) * cldfra1D
398                 cldfraliq(i)  = 0.
399                 icefrac(i)    = 1.
400                 dicefracdT(i) = 0.
401              ENDIF
402           
403
404           !--2. Pre-existing ice and/or mixing with environment:computation of ice properties for
405           !--feedback
406           ELSE
407
408              sursat_iceliq = qsatl(i)/qsati(i) - 1.
409              psati         = qsati(i) * pplay(i) / (RD/RV)
410             
411              !--We assume an exponential ice PSD whose parameters
412              !--are computed following Morrison&Gettelman 2008
413              !--Ice number density is assumed equals to INP density
414              !--which is for naero5>0 a function of temperature (DeMott 2010)   
415              !--bi and B0 are microphysical function characterizing
416              !--vapor/ice interactions
417              !--tau_phase_relax is the typical time of vapor deposition
418              !--onto ice crystals
419
420              !--For naero5<=0 INP density is derived from the empirical fit
421              !--from MARCUS campaign from Vignon 2021
422              !--/!\ Note that option is very specific and should be use for
423              !--the Southern Ocean and the Antarctic
424
425              IF (naero5 .LE. 0) THEN
426                 IF ( temp(i) .GT. tempvig1 ) THEN
427                      nb_crystals = 1.e3 * 10**(-0.14*(temp(i)-tempvig1) - 2.88)
428                 ELSE IF ( temp(i) .GT. tempvig2 ) THEN
429                      nb_crystals = 1.e3 * 10**(-0.31*(temp(i)-tempvig1) - 2.88)
430                 ELSE
431                      nb_crystals = 1.e3 * 10**(0.)
432                 ENDIF
433              ELSE
434                 nb_crystals = 1.e3 * 5.94e-5 * ( RTT - temp(i) )**3.33 * naero5**(0.0264*(RTT-temp(i))+0.0033)
435              ENDIF
436              lambda_PSD  = ( (RPI*rho_ice*nb_crystals) / (rho_air * MAX(qiceini_incl , eps) ) ) ** (1./3.)
437              N0_PSD      = nb_crystals * lambda_PSD
438              moment1_PSD = N0_PSD/lambda_PSD**2
439
440              !--Formulae for air thermal conductivity and water vapor diffusivity
441              !--comes respectively from Beard and Pruppacher (1971)
442              !--and  Hall and Pruppacher (1976)
443
444              air_thermal_conduct = ( 5.69 + 0.017 * ( temp(i) - RTT ) ) * 1.e-3 * 4.184
445              water_vapor_diff    = 2.11*1e-5 * ( temp(i) / RTT )**1.94 * ( 101325 / pplay(i) )
446             
447              bi = 1./((qsati(i)+qsatl(i))/2.) + RLSTT**2 / RCPD / RV / temp(i)**2
448              B0 = 4. * RPI * capa_crystal * 1. / (  RLSTT**2 / air_thermal_conduct / RV / temp(i)**2  &
449                                                  +  RV * temp(i) / psati / water_vapor_diff  )
450              invtau_phaserelax = bi * B0 * moment1_PSD
451             
452              ai = RG / RD / temp(i) * ( RD * RLSTT / RCPD / RV / temp(i) - 1. )
453              sursat_equ    = (ai * wvel(i) + sursat_e(i)*invtau_e(i)) / (invtau_phaserelax + invtau_e(i))
454              ! as sursaturation is by definition lower than -1 and
455              ! because local supersaturation > 1 are never found in the atmosphere
456
457              !--2A. No TKE : stationnary binary solution depending on vertical velocity and mixing with env.
458              ! If Sequ > Siw liquid cloud, else ice cloud
459              IF ( tke_dissip(i) .LE. eps )  THEN
460                 sigma2_icefracturb(i)= 0.
461                 mean_icefracturb(i)  = sursat_equ
462                 IF (sursat_equ .GT. sursat_iceliq) THEN
463                    qvap_cld(i)   = qsatl(i) * cldfra1D
464                    qliq(i)       = MAX(0.,qtot_incl(i)-qsatl(i)) * cldfra1D
465                    qice(i)       = 0.
466                    cldfraliq(i)  = 1.
467                    icefrac(i)    = 0.
468                    dicefracdT(i) = 0.
469                 ELSE
470                    qvap_cld(i)   = qsati(i) * cldfra1D
471                    qliq(i)       = 0.
472                    qice(i)       = MAX(0.,qtot_incl(i)-qsati(i)) * cldfra1D
473                    cldfraliq(i)  = 0.
474                    icefrac(i)    = 1.
475                    dicefracdT(i) = 0.
476                 ENDIF
477                 
478              !--2B. TKE and ice : ice supersaturation PDF
479              !--we compute the cloud liquid properties with a Gaussian PDF
480              !--of ice supersaturation F(Si) (Field2014, Furtado2016).
481              !--Parameters of the PDF are function of turbulence and
482              !--microphysics/existing ice.
483              ELSE 
484                     
485                 !--Tau_dissipturb is the time needed for turbulence to decay
486                 !--due to viscosity
487                 tau_dissipturb = gamma_taud * 2. * 2./3. * tke(i) / tke_dissip(i) / C0
488
489                 !--------------------- PDF COMPUTATIONS ---------------------
490                 !--Formulae for sigma2_pdf (variance), mean of PDF in Raillard2025
491                 !--cloud liquid fraction and in-cloud liquid content are given
492                 !--by integrating resp. F(Si) and Si*F(Si)
493                 !--Liquid is limited by the available water vapor trough a
494                 !--maximal liquid fraction
495                 !--qice_ini(i) / cldfra1D = qiceincld without precip
496
497                 liqfra_max = MAX(0., (MIN (1.,( qtot_incl(i) - (qice_ini(i) / cldfra1D) - qsati(i) * (1 + sursat_iceext ) ) / ( qsatl(i) - qsati(i) ) ) ) )
498                 sigma2_pdf = 1./2. * ( ai**2 ) *  2./3. * tke(i) * tau_dissipturb / (invtau_phaserelax + invtau_e(i))
499                 ! sursat ranges between -1 and 1, so we prevent sigma2 so exceed 1
500                 cldfraliq(i) = 0.5 * (1. - erf( ( sursat_iceliq - sursat_equ) / (SQRT(2.* sigma2_pdf) ) ) )
501                 IF (cldfraliq(i) .GT. liqfra_max) THEN
502                     cldfraliq(i) = liqfra_max
503                 ENDIF
504                 
505                 qliq_incl = qsati(i) * SQRT(sigma2_pdf) / SQRT(2.*RPI) * EXP( -1.*(sursat_iceliq - sursat_equ)**2. / (2.*sigma2_pdf) )  &
506                           - qsati(i) * cldfraliq(i) * (sursat_iceliq - sursat_equ )
507                 
508                 sigma2_icefracturb(i)= sigma2_pdf
509                 mean_icefracturb(i)  = sursat_equ
510     
511                 !------------ SPECIFIC VAPOR CONTENT AND WATER CONSERVATION  ------------
512
513                 IF ( (qliq_incl .LE. eps) .OR. (cldfraliq(i) .LE. eps) ) THEN
514                     qliq_incl    = 0.
515                     cldfraliq(i) = 0.
516                 END IF
517                 
518                 !--Specific humidity is the max between qsati and the weighted mean between
519                 !--qv in MPC patches and qv in ice-only parts. We assume that MPC parts are
520                 !--always at qsatl and ice-only parts slightly subsaturated (qsati*sursat_iceext+1)
521                 !--The whole cloud can therefore be supersaturated but never subsaturated.
522
523                 qvap_incl = MAX(qsati(i), ( 1. - cldfraliq(i) ) * (sursat_iceext + 1.) * qsati(i) + cldfraliq(i) * qsatl(i) )
524
525                 IF ( qvap_incl  .GE. qtot_incl(i) ) THEN
526                    qvap_incl = qsati(i)
527                    qliq_incl = MAX(0.0,qtot_incl(i) - qvap_incl)
528                    qice_incl = 0.
529
530                 ELSEIF ( (qvap_incl + qliq_incl) .GE. qtot_incl(i) ) THEN
531                    qliq_incl = MAX(0.0,qtot_incl(i) - qvap_incl)
532                    qice_incl = 0.
533                 ELSE
534                    qice_incl = qtot_incl(i) - qvap_incl - qliq_incl
535                 END IF
536
537                 qvap_cld(i)   = qvap_incl * cldfra1D
538                 qliq(i)       = qliq_incl * cldfra1D
539                 qice(i)       = qice_incl * cldfra1D
540                 IF ((qice(i)+qliq(i)) .GT. 0.) THEN
541                    icefrac(i)    = qice(i) / ( qice(i) + qliq(i) )
542                 ELSE
543                    icefrac(i)    = 1. ! to keep computation of qsat wrt ice in condensation loop in lmdz_lscp_main
544                 ENDIF
545                 dicefracdT(i) = 0.
546
547              END IF ! Enough TKE
548           
549           END IF ! End qini
550
551        END IF ! ! MPC temperature
552
553     END IF ! pticefracturb and cldfra
554   
555   ENDDO ! klon
556END SUBROUTINE ICEFRAC_LSCP_TURB
557!
558!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
559
560
561SUBROUTINE CALC_QSAT_ECMWF(klon,temp,qtot,pressure,tref,phase,flagth,qs,dqs)
562!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
563    ! Calculate qsat following ECMWF method
564!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
565
566
567USE yoethf_mod_h
568        USE yomcst_mod_h
569IMPLICIT NONE
570
571
572    include "FCTTRE.h"
573
574    INTEGER, INTENT(IN) :: klon  ! number of horizontal grid points
575    REAL, INTENT(IN), DIMENSION(klon) :: temp     ! temperature in K
576    REAL, INTENT(IN), DIMENSION(klon) :: qtot     ! total specific water in kg/kg
577    REAL, INTENT(IN), DIMENSION(klon) :: pressure ! pressure in Pa
578    REAL, INTENT(IN)                  :: tref     ! reference temperature in K
579    LOGICAL, INTENT(IN) :: flagth     ! flag for qsat calculation for thermals
580    INTEGER, INTENT(IN) :: phase
581    ! phase: 0=depend on temperature sign (temp>tref -> liquid, temp<tref, solid)
582    !        1=liquid
583    !        2=solid
584
585    REAL, INTENT(OUT), DIMENSION(klon) :: qs      ! saturation specific humidity [kg/kg]
586    REAL, INTENT(OUT), DIMENSION(klon) :: dqs     ! derivation of saturation specific humidity wrt T
587
588    REAL delta, cor, cvm5
589    INTEGER i
590
591    DO i=1,klon
592
593    IF (phase .EQ. 1) THEN
594        delta=0.
595    ELSEIF (phase .EQ. 2) THEN
596        delta=1.
597    ELSE
598        delta=MAX(0.,SIGN(1.,tref-temp(i)))
599    ENDIF
600
601    IF (flagth) THEN
602    cvm5=R5LES*(1.-delta) + R5IES*delta
603    ELSE
604    cvm5 = R5LES*RLVTT*(1.-delta) + R5IES*RLSTT*delta
605    cvm5 = cvm5 /RCPD/(1.0+RVTMP2*(qtot(i)))
606    ENDIF
607
608    qs(i)= R2ES*FOEEW(temp(i),delta)/pressure(i)
609    qs(i)=MIN(0.5,qs(i))
610    cor=1./(1.-RETV*qs(i))
611    qs(i)=qs(i)*cor
612    dqs(i)= FOEDE(temp(i),delta,cvm5,qs(i),cor)
613
614    END DO
615
616END SUBROUTINE CALC_QSAT_ECMWF
617!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
618
619
620!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
621SUBROUTINE CALC_GAMMASAT(klon,temp,qtot,pressure,gammasat,dgammasatdt)
622
623!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
624    ! programme that calculates the gammasat parameter that determines the
625    ! homogeneous condensation thresholds for cold (<0oC) clouds
626    ! condensation at q>gammasat*qsat
627    ! Etienne Vignon, March 2021
628!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
629
630    use lmdz_lscp_ini, only: iflag_gammasat, temp_nowater, RTT
631    use lmdz_lscp_ini, only: a_homofreez, b_homofreez, delta_hetfreez
632
633    IMPLICIT NONE
634
635
636    INTEGER, INTENT(IN) :: klon                       ! number of horizontal grid points
637    REAL, INTENT(IN), DIMENSION(klon) :: temp         ! temperature in K
638    REAL, INTENT(IN), DIMENSION(klon) :: qtot         ! total specific water in kg/kg
639
640    REAL, INTENT(IN), DIMENSION(klon) :: pressure     ! pressure in Pa
641
642    REAL, INTENT(OUT), DIMENSION(klon) :: gammasat    ! coefficient to multiply qsat with to calculate saturation
643    REAL, INTENT(OUT), DIMENSION(klon) :: dgammasatdt ! derivative of gammasat wrt temperature
644
645    REAL, DIMENSION(klon) ::  qsi,qsl,dqsl,dqsi
646    REAL  f_homofreez, fac
647
648    INTEGER i
649   
650        CALL CALC_QSAT_ECMWF(klon,temp,qtot,pressure,RTT,1,.false.,qsl,dqsl)
651        CALL CALC_QSAT_ECMWF(klon,temp,qtot,pressure,RTT,2,.false.,qsi,dqsi)
652
653    DO i = 1, klon
654
655        IF ( temp(i) .GE. RTT ) THEN
656            ! warm clouds: condensation at saturation wrt liquid
657            gammasat(i) = 1.
658            dgammasatdt(i) = 0.
659
660        ELSE
661            ! cold clouds: qsi > qsl
662           
663            ! homogeneous freezing of aerosols, according to
664            ! Koop, 2000 and Ren and MacKenzie, 2005 (QJRMS)
665            ! 'Cirrus regime'
666            ! if f_homofreez > qsl / qsi, liquid nucleation
667            ! if f_homofreez < qsl / qsi, homogeneous freezing of aerosols
668            ! Note: f_homofreez = qsl / qsi for temp ~= -38degC
669            f_homofreez = a_homofreez - temp(i) / b_homofreez
670           
671            IF ( iflag_gammasat .GE. 3 ) THEN
672              ! condensation at homogeneous freezing threshold for temp < -38 degC
673              ! condensation at liquid saturation for temp > -38 degC
674              IF ( f_homofreez .LE. qsl(i) / qsi(i) ) THEN
675                gammasat(i) = f_homofreez
676                dgammasatdt(i) = - 1. / b_homofreez
677              ELSE
678                gammasat(i) = qsl(i) / qsi(i)
679                dgammasatdt(i) = ( dqsl(i) * qsi(i) - dqsi(i) * qsl(i) ) / qsi(i) / qsi(i)
680              ENDIF
681
682            ELSEIF ( iflag_gammasat .EQ. 2 ) THEN
683              ! condensation at homogeneous freezing threshold for temp < -38 degC
684              ! condensation at a threshold linearly decreasing between homogeneous
685              ! freezing and ice saturation for -38 degC < temp < temp_nowater
686              ! condensation at ice saturation for temp > temp_nowater
687              ! If temp_nowater = 235.15 K, this is equivalent to iflag_gammasat = 1
688              IF ( f_homofreez .LE. qsl(i) / qsi(i) ) THEN
689                gammasat(i) = f_homofreez
690                dgammasatdt(i) = - 1. / b_homofreez
691              ELSEIF ( temp(i) .LE. temp_nowater ) THEN
692                ! Here, we assume that f_homofreez = qsl / qsi for temp = -38 degC = 235.15 K
693                dgammasatdt(i) = ( a_homofreez - 235.15 / b_homofreez - 1. ) &
694                               / ( 235.15 - temp_nowater )
695                gammasat(i) = dgammasatdt(i) * ( temp(i) - temp_nowater ) + 1.
696              ELSE
697                gammasat(i) = 1.
698                dgammasatdt(i) = 0.
699              ENDIF
700
701            ELSEIF ( iflag_gammasat .EQ. 1 ) THEN
702              ! condensation at homogeneous freezing threshold for temp < -38 degC
703              ! condensation at ice saturation for temp > -38 degC
704              IF ( f_homofreez .LE. qsl(i) / qsi(i) ) THEN
705                gammasat(i) = f_homofreez
706                dgammasatdt(i) = - 1. / b_homofreez
707              ELSE
708                gammasat(i) = 1.
709                dgammasatdt(i) = 0.
710              ENDIF
711
712            ELSE
713              ! condensation at ice saturation for temp < -38 degC
714              ! condensation at ice saturation for temp > -38 degC
715              gammasat(i) = 1.
716              dgammasatdt(i) = 0.
717
718            ENDIF
719
720            ! Note that the delta_hetfreez parameter allows to linearly decrease the
721            ! condensation threshold between the calculated threshold and the ice saturation
722            ! for delta_hetfreez = 1, the threshold is the calculated condensation threshold
723            ! for delta_hetfreez = 0, the threshold is the ice saturation
724            gammasat(i) = ( 1. - delta_hetfreez ) + delta_hetfreez * gammasat(i)
725            dgammasatdt(i) = delta_hetfreez * dgammasatdt(i)
726
727        ENDIF
728   
729    END DO
730
731
732END SUBROUTINE CALC_GAMMASAT
733!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
734
735
736!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
737SUBROUTINE DISTANCE_TO_CLOUD_TOP(klon,klev,k,temp,pplay,paprs,rneb,distcltop1D,temp_cltop)
738!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
739   
740   USE lmdz_lscp_ini, ONLY : rd,rg,tresh_cl
741
742   IMPLICIT NONE
743   
744   INTEGER, INTENT(IN) :: klon,klev                !number of horizontal and vertical grid points
745   INTEGER, INTENT(IN) :: k                        ! vertical index
746   REAL, INTENT(IN), DIMENSION(klon,klev) :: temp  ! temperature in K
747   REAL, INTENT(IN), DIMENSION(klon,klev) :: pplay ! pressure middle layer in Pa
748   REAL, INTENT(IN), DIMENSION(klon,klev+1) :: paprs ! pressure interfaces in Pa
749   REAL, INTENT(IN), DIMENSION(klon,klev) :: rneb  ! cloud fraction
750
751   REAL, INTENT(OUT), DIMENSION(klon) :: distcltop1D  ! distance from cloud top
752   REAL, INTENT(OUT), DIMENSION(klon) :: temp_cltop     ! temperature of cloud top
753   
754   REAL dzlay(klon,klev)
755   REAL zlay(klon,klev)
756   REAL dzinterf
757   INTEGER i,k_top, kvert
758   LOGICAL bool_cl
759
760
761   DO i=1,klon
762         ! Initialization height middle of first layer
763          dzlay(i,1) = Rd * temp(i,1) / rg * log(paprs(i,1)/paprs(i,2))
764          zlay(i,1) = dzlay(i,1)/2
765
766          DO kvert=2,klev
767                 IF (kvert.EQ.klev) THEN
768                       dzlay(i,kvert) = 2*(rd * temp(i,kvert) / rg * log(paprs(i,kvert)/pplay(i,kvert)))
769                 ELSE
770                       dzlay(i,kvert) = rd * temp(i,kvert) / rg * log(paprs(i,kvert)/paprs(i,kvert+1))
771                 ENDIF
772                       dzinterf       = rd * temp(i,kvert) / rg * log(pplay(i,kvert-1)/pplay(i,kvert))
773                       zlay(i,kvert)  = zlay(i,kvert-1) + dzinterf
774           ENDDO
775   ENDDO
776   
777   DO i=1,klon
778          k_top = k
779          IF (rneb(i,k) .LE. tresh_cl) THEN
780                 bool_cl = .FALSE.
781          ELSE
782                 bool_cl = .TRUE.
783          ENDIF
784
785          DO WHILE ((bool_cl) .AND. (k_top .LE. klev))
786          ! find cloud top
787                IF (rneb(i,k_top) .GT. tresh_cl) THEN
788                      k_top = k_top + 1
789                ELSE
790                      bool_cl = .FALSE.
791                      k_top   = k_top - 1
792                ENDIF
793          ENDDO
794          k_top=min(k_top,klev)
795
796          !dist to top is dist between current layer and layer of cloud top (from middle to middle) + dist middle to
797          !interf for layer of cloud top
798          distcltop1D(i) = zlay(i,k_top) - zlay(i,k) + dzlay(i,k_top)/2
799          temp_cltop(i)  = temp(i,k_top)
800   ENDDO ! klon
801
802END SUBROUTINE DISTANCE_TO_CLOUD_TOP
803!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
804
805!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
806FUNCTION GAMMAINC ( p, x )
807
808!*****************************************************************************80
809!
810!! GAMMAINC computes the regularized lower incomplete Gamma Integral
811!
812!  Modified:
813!
814!    20 January 2008
815!
816!  Author:
817!
818!    Original FORTRAN77 version by B Shea.
819!    FORTRAN90 version by John Burkardt.
820!
821!  Reference:
822!
823!    B Shea,
824!    Algorithm AS 239:
825!    Chi-squared and Incomplete Gamma Integral,
826!    Applied Statistics,
827!    Volume 37, Number 3, 1988, pages 466-473.
828!
829!  Parameters:
830!
831!    Input, real X, P, the parameters of the incomplete
832!    gamma ratio.  0 <= X, and 0 < P.
833!
834!    Output, real GAMMAINC, the value of the incomplete
835!    Gamma integral.
836!
837  IMPLICIT NONE
838
839  REAL A
840  REAL AN
841  REAL ARG
842  REAL B
843  REAL C
844  REAL, PARAMETER :: ELIMIT = - 88.0E+00
845  REAL GAMMAINC
846  REAL, PARAMETER :: OFLO = 1.0E+37
847  REAL P
848  REAL, PARAMETER :: PLIMIT = 1000.0E+00
849  REAL PN1
850  REAL PN2
851  REAL PN3
852  REAL PN4
853  REAL PN5
854  REAL PN6
855  REAL RN
856  REAL, PARAMETER :: TOL = 1.0E-14
857  REAL X
858  REAL, PARAMETER :: XBIG = 1.0E+08
859
860  GAMMAINC = 0.0E+00
861
862  IF ( X == 0.0E+00 ) THEN
863    GAMMAINC = 0.0E+00
864    RETURN
865  END IF
866!
867!  IF P IS LARGE, USE A NORMAL APPROXIMATION.
868!
869  IF ( PLIMIT < P ) THEN
870
871    PN1 = 3.0E+00 * SQRT ( P ) * ( ( X / P )**( 1.0E+00 / 3.0E+00 ) &
872    + 1.0E+00 / ( 9.0E+00 * P ) - 1.0E+00 )
873
874    GAMMAINC = 0.5E+00 * ( 1. + ERF ( PN1 ) )
875    RETURN
876
877  END IF
878!
879!  IF X IS LARGE SET GAMMAD = 1.
880!
881  IF ( XBIG < X ) THEN
882    GAMMAINC = 1.0E+00
883    RETURN
884  END IF
885!
886!  USE PEARSON'S SERIES EXPANSION.
887!  (NOTE THAT P IS NOT LARGE ENOUGH TO FORCE OVERFLOW IN ALOGAM).
888!
889  IF ( X <= 1.0E+00 .OR. X < P ) THEN
890
891    ARG = P * LOG ( X ) - X - LOG_GAMMA ( P + 1.0E+00 )
892    C = 1.0E+00
893    GAMMAINC = 1.0E+00
894    A = P
895
896    DO
897
898      A = A + 1.0E+00
899      C = C * X / A
900      GAMMAINC = GAMMAINC + C
901
902      IF ( C <= TOL ) THEN
903        EXIT
904      END IF
905
906    END DO
907
908    ARG = ARG + LOG ( GAMMAINC )
909
910    IF ( ELIMIT <= ARG ) THEN
911      GAMMAINC = EXP ( ARG )
912    ELSE
913      GAMMAINC = 0.0E+00
914    END IF
915!
916!  USE A CONTINUED FRACTION EXPANSION.
917!
918  ELSE
919
920    ARG = P * LOG ( X ) - X - LOG_GAMMA ( P )
921    A = 1.0E+00 - P
922    B = A + X + 1.0E+00
923    C = 0.0E+00
924    PN1 = 1.0E+00
925    PN2 = X
926    PN3 = X + 1.0E+00
927    PN4 = X * B
928    GAMMAINC = PN3 / PN4
929
930    DO
931
932      A = A + 1.0E+00
933      B = B + 2.0E+00
934      C = C + 1.0E+00
935      AN = A * C
936      PN5 = B * PN3 - AN * PN1
937      PN6 = B * PN4 - AN * PN2
938
939      IF ( PN6 /= 0.0E+00 ) THEN
940
941        RN = PN5 / PN6
942
943        IF ( ABS ( GAMMAINC - RN ) <= MIN ( TOL, TOL * RN ) ) THEN
944          EXIT
945        END IF
946
947        GAMMAINC = RN
948
949      END IF
950
951      PN1 = PN3
952      PN2 = PN4
953      PN3 = PN5
954      PN4 = PN6
955!
956!  RE-SCALE TERMS IN CONTINUED FRACTION IF TERMS ARE LARGE.
957!
958      IF ( OFLO <= ABS ( PN5 ) ) THEN
959        PN1 = PN1 / OFLO
960        PN2 = PN2 / OFLO
961        PN3 = PN3 / OFLO
962        PN4 = PN4 / OFLO
963      END IF
964
965    END DO
966
967    ARG = ARG + LOG ( GAMMAINC )
968
969    IF ( ELIMIT <= ARG ) THEN
970      GAMMAINC = 1.0E+00 - EXP ( ARG )
971    ELSE
972      GAMMAINC = 1.0E+00
973    END IF
974
975  END IF
976
977  RETURN
978END FUNCTION GAMMAINC
979!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
980
981END MODULE lmdz_lscp_tools
982
983
Note: See TracBrowser for help on using the repository browser.