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