[1617] | 1 | module updaterad |
---|
[740] | 2 | |
---|
| 3 | |
---|
| 4 | |
---|
| 5 | ! This module intents to group together all ice and dust radius computations done in the GCM, |
---|
| 6 | ! so that it stays coherent through the code and make (numerous) radius bugs easier to debug. |
---|
| 7 | ! All different thresholds values are defined below. |
---|
| 8 | |
---|
| 9 | ! Radius computation do not always occur on the whole grid (cf. improvedcloud). |
---|
| 10 | ! So, subroutines are designed for scalar values instead of tables |
---|
| 11 | |
---|
| 12 | |
---|
| 13 | ! T. Navarro, June 2012 |
---|
[1617] | 14 | ! CO2 clouds added 09/16 by J. Audouard |
---|
[740] | 15 | |
---|
| 16 | ! For instance, if R^3 is lower than r3icemin, then R is set to ricemin. |
---|
| 17 | ! So, ricemin should be the cubic root of r3icemin, but not necessarily ... |
---|
| 18 | |
---|
| 19 | real, parameter :: r3icemin = 1.e-30 ! ie ricemin = 0.0001 microns |
---|
| 20 | real, parameter :: ricemin = 1.e-10 |
---|
| 21 | |
---|
| 22 | real, parameter :: r3icemax = 125.e-12 ! ie ricemax = 500 microns |
---|
| 23 | real, parameter :: ricemax = 500.e-6 |
---|
| 24 | |
---|
[1617] | 25 | real, parameter :: r3iceco2min = 1.e-30 |
---|
| 26 | real, parameter :: riceco2min = 1.e-10 |
---|
[740] | 27 | |
---|
[1617] | 28 | real, parameter :: r3iceco2max = 125.e-12 |
---|
| 29 | real, parameter :: riceco2max = 500.e-6 |
---|
| 30 | |
---|
| 31 | |
---|
| 32 | real, parameter :: qice_threshold = 1.e-15 ! 1.e-10 |
---|
| 33 | real, parameter :: qice_co2_threshold = 1.e-30 ! 1.e-10 |
---|
| 34 | |
---|
[740] | 35 | real, parameter :: nccn_threshold = 1. |
---|
[1619] | 36 | real, parameter :: qccn_threshold = 1.e-20 |
---|
[740] | 37 | |
---|
[1619] | 38 | real, parameter :: r3ccnmin = 1.e-21 ! ie rccnmin = 0.1 microns |
---|
| 39 | real, parameter :: rccnmin = 0.1e-6 |
---|
[740] | 40 | |
---|
[1619] | 41 | real, parameter :: r3ccnmax = 125.e-12 ! ie rccnmax = 500 microns |
---|
| 42 | real, parameter :: rccnmax = 500.e-6 |
---|
[740] | 43 | |
---|
| 44 | |
---|
| 45 | |
---|
| 46 | real, parameter :: ndust_threshold = 1. |
---|
| 47 | real, parameter :: qdust_threshold = 1.e-20 |
---|
| 48 | |
---|
[1617] | 49 | real, parameter :: r3dustmin = 1.e-24 ! ie rdustmin = 0.01 microns |
---|
| 50 | real, parameter :: rdustmin = 1.e-8 |
---|
[740] | 51 | |
---|
| 52 | real, parameter :: r3dustmax = 125.e-12 ! ie rdustmax = 500 microns |
---|
| 53 | real, parameter :: rdustmax = 500.e-6 |
---|
| 54 | |
---|
| 55 | |
---|
[745] | 56 | real, parameter :: rdust0 = 0.8e-6 |
---|
[740] | 57 | |
---|
| 58 | |
---|
| 59 | |
---|
| 60 | |
---|
| 61 | contains |
---|
| 62 | |
---|
| 63 | |
---|
| 64 | !============================================================================ |
---|
| 65 | !============================================================================ |
---|
| 66 | !============================================================================ |
---|
| 67 | ! Update ice radius if microphys == true |
---|
| 68 | subroutine updaterice_micro(qice,qccn,nccn,coeff,rice,rhocloud) |
---|
[1036] | 69 | use tracer_mod, only: rho_dust, rho_ice |
---|
[1226] | 70 | USE comcstfi_h |
---|
[740] | 71 | implicit none |
---|
| 72 | |
---|
| 73 | real, intent(in) :: qice,qccn,nccn |
---|
| 74 | real, intent(in) :: coeff ! this coeff is tauscaling if microphy = T (possibly ccn_factor^-1 otherwise) |
---|
| 75 | real, intent(out) :: rice,rhocloud ! rhocloud is needed for sedimentation and is also a good diagnostic variable |
---|
| 76 | real nccn_true,qccn_true ! radius of the ice crystal core to the power of 3 |
---|
| 77 | |
---|
| 78 | |
---|
| 79 | nccn_true = max(nccn * coeff, 1.e-30) |
---|
| 80 | qccn_true = max(qccn * coeff, 1.e-30) |
---|
| 81 | |
---|
| 82 | |
---|
| 83 | !! Nota: It is very dangerous to apply a threshold on qccn or nccn to force rice to be ricemin. |
---|
| 84 | !! Indeed, one can obtain ricemin for small but non negligible qice values, and therefore hugely opaque clouds. |
---|
| 85 | |
---|
| 86 | |
---|
[1307] | 87 | rhocloud = (qice*rho_ice + qccn_true*rho_dust) / (qice + qccn_true) |
---|
| 88 | rhocloud = min(max(rhocloud,rho_ice),rho_dust) |
---|
[740] | 89 | |
---|
[1307] | 90 | rice = (qice + qccn_true) * 0.75 / pi / rhocloud / nccn_true |
---|
[740] | 91 | |
---|
[1307] | 92 | if (rice .le. r3icemin) then |
---|
| 93 | rice = ricemin |
---|
| 94 | else if (rice .ge. r3icemax) then |
---|
| 95 | rice = ricemax |
---|
[740] | 96 | else |
---|
[1307] | 97 | rice = rice**(1./3.) ! here rice is always positive |
---|
| 98 | endif |
---|
[740] | 99 | |
---|
| 100 | |
---|
| 101 | end subroutine updaterice_micro |
---|
| 102 | !============================================================================ |
---|
| 103 | !============================================================================ |
---|
| 104 | !============================================================================ |
---|
| 105 | |
---|
[1617] | 106 | subroutine updaterice_microco2(qice,qccn,nccn,coeff,rice,rhocloudco2) |
---|
| 107 | use tracer_mod, only: rho_dust, rho_ice_co2 |
---|
| 108 | USE comcstfi_h, only: pi |
---|
| 109 | implicit none |
---|
| 110 | !By CL and JA 09/16 |
---|
[740] | 111 | |
---|
[1617] | 112 | real, intent(in) :: qice,qccn,nccn |
---|
| 113 | real, intent(in) :: coeff ! this coeff is tauscaling if microphy = T (possibly ccn_factor^-1 otherwise) |
---|
| 114 | real, intent(out) :: rice,rhocloudco2 ! rhocloud is needed for sedimentation and is also a good diagnostic variable |
---|
| 115 | real nccn_true,qccn_true ! nombre et masse de CCN |
---|
| 116 | |
---|
| 117 | nccn_true = max(nccn * coeff, 1.e-30) |
---|
| 118 | qccn_true = max(qccn * coeff, 1.e-30) |
---|
[740] | 119 | |
---|
| 120 | |
---|
[1617] | 121 | rhocloudco2 = (qice *rho_ice_co2 + qccn_true*rho_dust) / (qice + qccn_true) |
---|
[740] | 122 | |
---|
[1617] | 123 | rhocloudco2 = min(max(rhocloudco2,rho_ice_co2),rho_dust) |
---|
| 124 | |
---|
| 125 | rice = (qice + qccn_true) * 0.75 / pi / rhocloudco2 / nccn_true |
---|
| 126 | |
---|
| 127 | if (rice .le. r3iceco2min) then !r3icemin radius power 3 ? |
---|
| 128 | rice = riceco2min |
---|
| 129 | else if (rice .ge. r3iceco2max) then !r3icemin radius power 3 ? |
---|
| 130 | rice = riceco2max |
---|
| 131 | else |
---|
| 132 | rice = rice**(1./3.) ! here rice is always positive |
---|
| 133 | endif |
---|
| 134 | |
---|
| 135 | |
---|
| 136 | end subroutine updaterice_microco2 |
---|
| 137 | |
---|
| 138 | |
---|
| 139 | |
---|
[740] | 140 | !============================================================================ |
---|
| 141 | !============================================================================ |
---|
| 142 | !============================================================================ |
---|
| 143 | ! Update ice radius from a typical profile if microphys == false |
---|
| 144 | subroutine updaterice_typ(qice,tau,pzlay,rice) |
---|
[1036] | 145 | use tracer_mod, only: rho_ice |
---|
[1226] | 146 | USE comcstfi_h |
---|
[740] | 147 | implicit none |
---|
| 148 | |
---|
| 149 | real, intent(in) :: qice |
---|
| 150 | real, intent(in) :: tau ! tau for dust |
---|
| 151 | real, intent(in) :: pzlay ! altitude at the middle of the layers |
---|
| 152 | real, intent(out) :: rice |
---|
| 153 | real rccn,nccn ! radius and number of ice crystals |
---|
| 154 | |
---|
| 155 | ! Typical CCN profile following Montmessin et al. 2004 |
---|
| 156 | ! (N0=2e6 m-3 has been converted into N0=1.3e8 kg-1, otherwise the equation for rice is not homogeneous...) |
---|
| 157 | nccn = 1.3e+8*max(tau,0.001)/0.1*exp(-pzlay/10000.) |
---|
| 158 | ! The previously used profile was not correct: |
---|
| 159 | ! Nccn=( epaisseur/masse ) * 2.e+6/0.1*max(tau,0.001)*exp(-pzlay/10000.) |
---|
| 160 | |
---|
| 161 | if (nccn .le. 1) then |
---|
| 162 | |
---|
| 163 | rice = ricemin |
---|
| 164 | |
---|
| 165 | else |
---|
| 166 | |
---|
| 167 | ! Typical dust radius profile: |
---|
| 168 | rccn = max(rdust0*exp(-pzlay/18000.),1.e-9) |
---|
| 169 | rice = qice * 0.75 / pi / rho_ice / nccn + rccn*rccn*rccn |
---|
| 170 | |
---|
| 171 | if (rice .le. r3icemin) then |
---|
| 172 | rice = ricemin |
---|
| 173 | else if (rice .ge. r3icemax) then |
---|
| 174 | rice = ricemax |
---|
| 175 | else |
---|
| 176 | rice = rice**(1./3.) ! here rice is always positive |
---|
| 177 | endif |
---|
| 178 | |
---|
| 179 | endif |
---|
| 180 | |
---|
| 181 | end subroutine updaterice_typ |
---|
| 182 | !============================================================================ |
---|
| 183 | !============================================================================ |
---|
| 184 | !============================================================================ |
---|
| 185 | |
---|
| 186 | |
---|
| 187 | |
---|
| 188 | |
---|
| 189 | |
---|
| 190 | |
---|
| 191 | !============================================================================ |
---|
| 192 | !============================================================================ |
---|
| 193 | !============================================================================ |
---|
| 194 | ! This subroutine computes the geometric mean radius(or number median radius) |
---|
| 195 | ! For a lognormal distribution : |
---|
| 196 | ! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2) |
---|
| 197 | ! To be used with doubleq == true. otherwise, rdust is constant !!! |
---|
| 198 | subroutine updaterdust(qdust,ndust,rdust,tauscaling) |
---|
[1036] | 199 | use tracer_mod, only: r3n_q |
---|
[1226] | 200 | USE comcstfi_h |
---|
[740] | 201 | implicit none |
---|
| 202 | |
---|
| 203 | real, intent(in) :: qdust,ndust ! needed if doubleq |
---|
| 204 | real, intent(in), optional :: tauscaling ! useful for realistic thresholds |
---|
| 205 | real, intent(out) :: rdust |
---|
| 206 | real coeff |
---|
| 207 | |
---|
| 208 | if (present(tauscaling)) then |
---|
| 209 | coeff = tauscaling ! thresholds on realistic values |
---|
| 210 | else |
---|
| 211 | coeff = 1. ! thresholds on virtual values |
---|
| 212 | endif |
---|
| 213 | |
---|
| 214 | if ((ndust .le. ndust_threshold/coeff) .or. (qdust .le. qdust_threshold/coeff)) then |
---|
| 215 | |
---|
| 216 | rdust = rdustmin |
---|
| 217 | |
---|
| 218 | else |
---|
| 219 | |
---|
| 220 | rdust = r3n_q * qdust / ndust |
---|
| 221 | |
---|
| 222 | if (rdust .le. r3dustmin) then |
---|
| 223 | rdust = rdustmin |
---|
| 224 | else if (rdust .ge. r3dustmax) then |
---|
| 225 | rdust = rdustmax |
---|
| 226 | else |
---|
| 227 | rdust = rdust**(1./3.) |
---|
| 228 | endif |
---|
| 229 | |
---|
| 230 | endif |
---|
| 231 | |
---|
| 232 | |
---|
| 233 | end subroutine updaterdust |
---|
| 234 | !============================================================================ |
---|
| 235 | !============================================================================ |
---|
| 236 | !============================================================================ |
---|
| 237 | |
---|
| 238 | |
---|
| 239 | |
---|
| 240 | |
---|
| 241 | |
---|
| 242 | |
---|
| 243 | !============================================================================ |
---|
| 244 | !============================================================================ |
---|
| 245 | !============================================================================ |
---|
| 246 | ! This subroutine computes the mass mean radius, |
---|
| 247 | ! used for heterogenous nucleation on CCNs in microphysics. |
---|
| 248 | ! For a lognormal distribution : |
---|
| 249 | ! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2) |
---|
| 250 | subroutine updaterccn(qccn,nccn,rccn,tauscaling) |
---|
[1036] | 251 | use tracer_mod, only: rho_dust |
---|
[1226] | 252 | USE comcstfi_h |
---|
[740] | 253 | implicit none |
---|
| 254 | |
---|
| 255 | real, intent(in) :: qccn,nccn ! needed if doubleq |
---|
| 256 | real, intent(in), optional :: tauscaling ! useful for realistic thresholds |
---|
| 257 | |
---|
| 258 | real, intent(out) :: rccn |
---|
| 259 | |
---|
| 260 | real coeff |
---|
| 261 | |
---|
| 262 | |
---|
| 263 | if (present(tauscaling)) then |
---|
| 264 | coeff = tauscaling ! threshold on realistic values |
---|
| 265 | else |
---|
| 266 | coeff = 1. ! threshold on virtual values |
---|
| 267 | endif |
---|
| 268 | |
---|
| 269 | if ((nccn .le. nccn_threshold/coeff) .or. (qccn .le. qccn_threshold/coeff)) then |
---|
| 270 | |
---|
| 271 | rccn = rccnmin |
---|
| 272 | |
---|
| 273 | else |
---|
| 274 | |
---|
| 275 | rccn = qccn * 0.75 / pi / rho_dust / nccn |
---|
| 276 | |
---|
| 277 | if (rccn .le. r3ccnmin) then |
---|
| 278 | rccn = rccnmin |
---|
| 279 | else if (rccn .ge. r3ccnmax) then |
---|
| 280 | rccn = rccnmax |
---|
| 281 | else |
---|
| 282 | rccn = rccn**(1./3.) |
---|
| 283 | endif |
---|
| 284 | |
---|
| 285 | endif |
---|
| 286 | |
---|
| 287 | end subroutine updaterccn |
---|
| 288 | !============================================================================ |
---|
| 289 | !============================================================================ |
---|
| 290 | !============================================================================ |
---|
| 291 | |
---|
| 292 | |
---|
| 293 | |
---|
| 294 | end module updaterad |
---|
| 295 | |
---|