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

Last change on this file since 1617 was 1617, checked in by jaudouard, 8 years ago

Added modifications for CO2 clouds scheme in physiq_mod.F and added several routines and variables for CO2 microphysics. October 2016 Joachim Audouard

File size: 8.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
13! T. Navarro, June 2012
14! CO2 clouds added 09/16 by J. Audouard
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
19real, parameter :: r3icemin = 1.e-30   ! ie ricemin  = 0.0001 microns
20real, parameter :: ricemin  = 1.e-10
21
22real, parameter :: r3icemax = 125.e-12 ! ie ricemax  = 500 microns
23real, parameter :: ricemax  = 500.e-6
24
25real, parameter :: r3iceco2min = 1.e-30
26real, parameter :: riceco2min  = 1.e-10
27
28real, parameter :: r3iceco2max = 125.e-12
29real, parameter :: riceco2max  = 500.e-6
30
31
32real, parameter :: qice_threshold  = 1.e-15 ! 1.e-10
33real, parameter :: qice_co2_threshold  = 1.e-30 ! 1.e-10
34
35real, parameter :: nccn_threshold  =  1.
36real, parameter :: qccn_threshold  =  1.e-30
37
38real, parameter :: r3ccnmin = 1.e-24    ! ie rccnmin = 10n m
39real, parameter :: rccnmin  = 1.e-8
40
41real, parameter :: r3ccnmax = 125.e-18  ! ie rccnmax  = 5 microns
42real, parameter :: rccnmax  = 5.e-6
43
44
45
46real, parameter :: ndust_threshold  = 1.
47real, parameter :: qdust_threshold  = 1.e-20
48
49real, parameter :: r3dustmin = 1.e-24  ! ie rdustmin = 0.01 microns
50real, parameter :: rdustmin  = 1.e-8
51
52real, parameter :: r3dustmax = 125.e-12 ! ie rdustmax  = 500 microns
53real, parameter :: rdustmax  = 500.e-6 
54
55
56real, parameter :: rdust0 = 0.8e-6
57
58
59     
60
61contains
62
63
64!============================================================================
65!============================================================================
66!============================================================================
67! Update ice radius if microphys == true
68subroutine updaterice_micro(qice,qccn,nccn,coeff,rice,rhocloud)
69use tracer_mod, only: rho_dust, rho_ice
70USE comcstfi_h
71implicit none
72
73real, intent(in)  :: qice,qccn,nccn
74real, intent(in)  :: coeff         ! this coeff is tauscaling if microphy = T (possibly ccn_factor^-1 otherwise)
75real, intent(out) :: rice,rhocloud ! rhocloud is needed for sedimentation and is also a good diagnostic variable
76real nccn_true,qccn_true ! radius of the ice crystal core to the power of 3
77
78   
79nccn_true = max(nccn * coeff, 1.e-30)
80qccn_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
87rhocloud = (qice*rho_ice + qccn_true*rho_dust) / (qice + qccn_true)
88rhocloud = min(max(rhocloud,rho_ice),rho_dust)
89
90rice = (qice + qccn_true) * 0.75 / pi / rhocloud / nccn_true
91 
92if (rice .le. r3icemin) then
93  rice = ricemin
94else if (rice .ge. r3icemax) then
95  rice = ricemax
96else
97  rice = rice**(1./3.) ! here rice is always positive
98endif
99
100
101end subroutine updaterice_micro
102!============================================================================
103!============================================================================
104!============================================================================
105
106subroutine updaterice_microco2(qice,qccn,nccn,coeff,rice,rhocloudco2)
107use tracer_mod, only: rho_dust, rho_ice_co2
108USE comcstfi_h, only:  pi
109implicit none
110!By CL and JA 09/16
111
112real, intent(in)  :: qice,qccn,nccn
113real, intent(in)  :: coeff         ! this coeff is tauscaling if microphy = T (possibly ccn_factor^-1 otherwise)
114real, intent(out) :: rice,rhocloudco2 ! rhocloud is needed for sedimentation and is also a good diagnostic variable
115real nccn_true,qccn_true ! nombre et masse de CCN
116   
117nccn_true = max(nccn * coeff, 1.e-30)
118qccn_true = max(qccn * coeff, 1.e-30)
119
120
121  rhocloudco2 = (qice *rho_ice_co2 + qccn_true*rho_dust) / (qice + qccn_true)
122
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
136end subroutine updaterice_microco2
137
138
139
140!============================================================================
141!============================================================================
142!============================================================================
143! Update ice radius from a typical profile if microphys == false
144subroutine updaterice_typ(qice,tau,pzlay,rice)
145use tracer_mod, only: rho_ice
146USE comcstfi_h
147implicit none
148
149real, intent(in)  :: qice
150real, intent(in)  :: tau   ! tau for dust
151real, intent(in)  :: pzlay ! altitude at the middle of the layers
152real, intent(out) :: rice
153real 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     
161if (nccn .le. 1) then
162
163  rice = ricemin
164
165else
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
179endif
180 
181end 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 !!!
198subroutine updaterdust(qdust,ndust,rdust,tauscaling)
199use tracer_mod, only: r3n_q
200USE comcstfi_h
201implicit none
202
203real, intent(in) :: qdust,ndust ! needed if doubleq
204real, intent(in), optional :: tauscaling ! useful for realistic thresholds
205real, intent(out) :: rdust
206real coeff
207
208if (present(tauscaling)) then
209  coeff = tauscaling ! thresholds on realistic values
210else
211  coeff = 1. ! thresholds on virtual values
212endif
213
214if ((ndust .le. ndust_threshold/coeff) .or. (qdust .le. qdust_threshold/coeff)) then
215
216   rdust = rdustmin
217
218else
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
230endif
231     
232
233end 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)
250subroutine updaterccn(qccn,nccn,rccn,tauscaling)
251use tracer_mod, only: rho_dust
252USE comcstfi_h
253implicit none
254
255real, intent(in) :: qccn,nccn ! needed if doubleq
256real, intent(in), optional :: tauscaling ! useful for realistic thresholds
257
258real, intent(out) :: rccn
259
260real coeff
261
262
263if (present(tauscaling)) then
264  coeff = tauscaling ! threshold on realistic values
265else
266  coeff = 1. ! threshold on virtual values
267endif
268
269if ((nccn .le. nccn_threshold/coeff) .or. (qccn .le. qccn_threshold/coeff)) then
270
271   rccn = rccnmin
272
273else
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
285endif
286     
287end subroutine updaterccn
288!============================================================================
289!============================================================================
290!============================================================================
291
292
293
294end module updaterad
295
Note: See TracBrowser for help on using the repository browser.