source: trunk/LMDZ.MARS/libf/phymars/updaterad.F90 @ 3118

Last change on this file since 3118 was 2562, checked in by cmathe, 3 years ago

GCM MARS: CO2 clouds microphysics improvements

File size: 10.9 KB
Line 
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! T. Navarro, June 2012
13! CO2 clouds added 09/16 by J. Audouard
14
15! For instance, if R^3 is lower than r3icemin, then R is set to ricemin.
16! So, ricemin should be the cubic root of r3icemin, but not necessarily ...
17
18real, parameter :: r3icemin = 1.e-30   ! ie ricemin  = 0.0001 microns
19real, parameter :: ricemin  = 1.e-10
20real, parameter :: r3icemax = 125.e-12 ! ie ricemax  = 500 microns
21real, parameter :: ricemax  = 500.e-6
22
23double precision, parameter :: r3iceco2min = 1.e-30
24double precision, parameter :: riceco2min  = 1.e-10
25double precision, parameter :: r3iceco2max = 125.e-12
26double precision, parameter :: riceco2max  = 500.e-6
27
28real, parameter :: qice_threshold  = 1.e-15 ! 1.e-10
29real, parameter :: qice_co2_threshold  = 1.e-30 ! 1.e-10
30
31real, parameter :: nccn_threshold  =  1.
32real, parameter :: qccn_threshold  =  1.e-20
33
34real, parameter :: r3ccnmin = 1.e-21    ! ie rccnmin = 0.1 microns
35real, parameter :: rccnmin  = 0.1e-6
36real, parameter :: rccnCO2min  = 1e-9
37
38real, parameter :: r3ccnmax = 125.e-12  ! ie rccnmax  = 500 microns
39real, parameter :: rccnmax  = 500.e-6
40
41
42
43real, parameter :: ndust_threshold  = 1.
44real, parameter :: qdust_threshold  = 1.e-20
45
46real, parameter :: r3dustmin = 1.e-24  ! ie rdustmin = 0.01 microns
47real, parameter :: rdustmin  = 1.e-8
48
49real, parameter :: r3dustmax = 125.e-12 ! ie rdustmax  = 500 microns
50real, parameter :: rdustmax  = 500.e-6 
51
52
53real, parameter :: rdust0 = 0.8e-6
54
55
56     
57
58contains
59
60
61!============================================================================
62!============================================================================
63!============================================================================
64! Update ice radius if microphys == true
65subroutine updaterice_micro(qice,qccn,nccn,coeff,rice,rhocloud)
66use tracer_mod, only: rho_dust, rho_ice
67USE comcstfi_h
68implicit none
69
70real, intent(in)  :: qice,qccn,nccn
71real, intent(in)  :: coeff         ! this coeff is tauscaling if microphy = T (possibly ccn_factor^-1 otherwise)
72real, intent(out) :: rice,rhocloud ! rhocloud is needed for sedimentation and is also a good diagnostic variable
73real nccn_true,qccn_true ! radius of the ice crystal core to the power of 3
74
75   
76nccn_true = max(nccn * coeff, 1.e-30)
77qccn_true = max(qccn * coeff, 1.e-30)
78
79
80!! Nota: It is very dangerous to apply a threshold on qccn or nccn to force rice to be ricemin.
81!! Indeed, one can obtain ricemin for small but non negligible qice values, and therefore hugely opaque clouds.
82 
83
84rhocloud = (qice*rho_ice + qccn_true*rho_dust) / (qice + qccn_true)
85rhocloud = min(max(rhocloud,rho_ice),rho_dust)
86
87rice = (qice + qccn_true) * 0.75 / pi / rhocloud / nccn_true
88 
89if (rice .le. r3icemin) then
90  rice = ricemin
91else if (rice .ge. r3icemax) then
92  rice = ricemax
93else
94  rice = rice**(1./3.) ! here rice is always positive
95endif
96
97
98end subroutine updaterice_micro
99!============================================================================
100!============================================================================
101!============================================================================
102subroutine updaterice_microco2(qice, qccn, nccn, qccn_h2o, nccn_h2o, temperature, coeff, rice, rhocloudco2)
103use tracer_mod, only: rho_dust, rho_ice
104use comcstfi_h, only: pi
105use density_co2_ice_mod, only: density_co2_ice
106
107implicit none
108!CO2 clouds parameter update by CL and JA 09/16
109
110DOUBLE PRECISION, intent(in)  :: qice,qccn,nccn, qccn_h2o, nccn_h2o
111real, intent(in)  :: coeff         ! this coeff is tauscaling if microphy = T (possibly ccn_factor^-1 otherwise)
112real, intent(in)  :: temperature ! temperature effective for density co2_ice computation
113real, intent(out) :: rhocloudco2 ! rhocloud is needed for sedimentation and is also a good diagnostic variable
114double precision, intent(out) :: rice
115real nccn_true,qccn_true, qccn_h2o_true, nccn_h2o_true ! nombre et masse de CCN
116double precision :: rho_ice_co2T ! density co2_ice Temperature-dependent
117
118nccn_true = max(nccn * coeff, 1.e-30)
119qccn_true = max(qccn * coeff, 1.e-30)
120
121nccn_h2o_true = max(nccn_h2o, 1.e-30)
122qccn_h2o_true = max(qccn_h2o, 1.e-30)
123
124call density_co2_ice(dble(temperature), rho_ice_co2T)
125
126rhocloudco2 = (qice * rho_ice_co2T + qccn_true*rho_dust + qccn_h2o_true*rho_ice) / (qice + qccn_true + qccn_h2o_true)
127
128rhocloudco2 = min(max(rhocloudco2,rho_ice_co2T), rho_dust)
129
130rice = (qice + qccn_true + qccn_h2o_true) * 0.75 / pi / rhocloudco2 / (nccn_true + nccn_h2o_true)
131
132rhocloudco2 = (qice * rho_ice_co2T + qccn_true*rho_dust + qccn_h2o_true*rho_ice) / (qice + qccn_true + qccn_h2o_true)
133
134if (rice <= r3iceco2min) then !r3icemin radius power 3 ?
135  rice = riceco2min
136else if (rice >= r3iceco2max) then !r3icemin radius power 3 ?
137  rice = riceco2max
138else
139  rice = rice**(1./3.) ! here rice is always positive
140endif
141
142
143end subroutine updaterice_microco2
144
145
146
147!============================================================================
148!============================================================================
149!============================================================================
150! Update ice radius from a typical profile if microphys == false
151subroutine updaterice_typ(qice,tau,pzlay,rice)
152use tracer_mod, only: rho_ice
153USE comcstfi_h
154implicit none
155
156real, intent(in)  :: qice
157real, intent(in)  :: tau   ! tau for dust
158real, intent(in)  :: pzlay ! altitude at the middle of the layers
159real, intent(out) :: rice
160real rccn,nccn ! radius  and number of ice crystals
161
162!   Typical CCN profile following Montmessin et al. 2004
163!   (N0=2e6 m-3 has been converted into N0=1.3e8 kg-1, otherwise the equation for rice is not homogeneous...)
164     nccn  = 1.3e+8*max(tau,0.001)/0.1*exp(-pzlay/10000.)
165!   The previously used profile was not correct:
166!   Nccn=( epaisseur/masse ) * 2.e+6/0.1*max(tau,0.001)*exp(-pzlay/10000.)
167     
168if (nccn .le. 1) then
169
170  rice = ricemin
171
172else
173
174! Typical dust radius profile:
175  rccn = max(rdust0*exp(-pzlay/18000.),1.e-9)
176  rice  = qice * 0.75 / pi / rho_ice / nccn + rccn*rccn*rccn
177       
178  if (rice .le. r3icemin) then
179    rice = ricemin
180  else if (rice .ge. r3icemax) then
181    rice = ricemax
182  else
183    rice = rice**(1./3.) ! here rice is always positive
184  endif
185
186endif
187 
188end subroutine updaterice_typ
189!============================================================================
190!============================================================================
191!============================================================================
192
193
194
195
196
197
198!============================================================================
199!============================================================================
200!============================================================================
201! This subroutine computes the geometric mean radius(or number median radius)
202! For a lognormal distribution :
203! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
204! To be used with doubleq == true. otherwise, rdust is constant !!!
205subroutine updaterdust(qdust,ndust,rdust,tauscaling)
206use tracer_mod, only: r3n_q
207USE comcstfi_h
208implicit none
209
210real, intent(in) :: qdust,ndust ! needed if doubleq
211real, intent(in), optional :: tauscaling ! useful for realistic thresholds
212real, intent(out) :: rdust
213real coeff
214
215if (present(tauscaling)) then
216  coeff = tauscaling ! thresholds on realistic values
217else
218  coeff = 1. ! thresholds on virtual values
219endif
220
221if ((ndust .le. ndust_threshold/coeff) .or. (qdust .le. qdust_threshold/coeff)) then
222
223   rdust = rdustmin
224
225else
226
227   rdust  = r3n_q * qdust / ndust
228
229   if (rdust .le. r3dustmin) then
230       rdust = rdustmin
231   else if (rdust .ge. r3dustmax) then
232       rdust = rdustmax
233   else
234      rdust = rdust**(1./3.)
235   endif
236
237endif
238     
239
240end subroutine updaterdust
241!============================================================================
242!============================================================================
243!============================================================================
244
245
246
247
248
249
250!============================================================================
251!============================================================================
252!============================================================================
253! This subroutine computes the mass mean radius,
254! used for heterogenous nucleation on CCNs in microphysics.
255! For a lognormal distribution :
256! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
257subroutine updaterccn(qccn,nccn,rccn,tauscaling)
258use tracer_mod, only: rho_dust
259USE comcstfi_h
260implicit none
261
262real, intent(in) :: qccn,nccn ! needed if doubleq
263real, intent(in), optional :: tauscaling ! useful for realistic thresholds
264
265real, intent(out) :: rccn
266
267real coeff
268
269
270if (present(tauscaling)) then
271  coeff = tauscaling ! threshold on realistic values
272else
273  coeff = 1. ! threshold on virtual values
274endif
275
276if ((nccn .le. nccn_threshold/coeff) .or. (qccn .le. qccn_threshold/coeff)) then
277
278   rccn = rccnmin
279
280else
281
282   rccn  =  qccn * 0.75 / pi / rho_dust / nccn
283
284   if (rccn .le. r3ccnmin) then
285       rccn = rccnmin
286   else if (rccn .ge. r3ccnmax) then
287       rccn = rccnmax
288   else
289       rccn = rccn**(1./3.)
290   endif
291
292endif
293     
294end subroutine updaterccn
295!============================================================================
296!============================================================================
297!============================================================================
298
299
300
301
302
303!============================================================================
304!============================================================================
305!============================================================================
306! This subroutine computes the mass mean radius,
307! used for heterogenous nucleation on CCNs in microphysics of CO2.
308! For a lognormal distribution :
309! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
310subroutine updaterccnCO2(qccn,nccn,rccn,tauscaling)
311use tracer_mod, only: rho_dust
312USE comcstfi_h
313implicit none
314
315real, intent(in) :: qccn,nccn ! needed if doubleq
316real, intent(in), optional :: tauscaling ! useful for realistic thresholds
317
318real, intent(out) :: rccn
319
320real coeff
321
322
323if (present(tauscaling)) then
324  coeff = tauscaling ! threshold on realistic values
325else
326  coeff = 1. ! threshold on virtual values
327endif
328
329if ((nccn .le. nccn_threshold/coeff) .or. (qccn .le. qccn_threshold/coeff)) then
330
331   rccn = rccnCO2min
332
333else
334
335   rccn  =  qccn * 0.75 / pi / rho_dust / nccn
336   rccn = rccn**(1./3.)
337
338endif
339
340rccn=min(5.E-4,rccn)
341     
342end subroutine updaterccnCO2
343!============================================================================
344!============================================================================
345!============================================================================
346
347
348
349end module updaterad
350
Note: See TracBrowser for help on using the repository browser.