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