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

Last change on this file since 1920 was 1816, checked in by jaudouard, 7 years ago

Commit for CO2 clouds microphysics.

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