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

Last change on this file since 1754 was 1632, checked in by aslmd, 8 years ago

compiler complain about type mismatch (double vs. real)

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
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
25double precision, parameter :: r3iceco2min = 1.e-30
26double precision, parameter :: riceco2min  = 1.e-10
27
28double precision, parameter :: r3iceco2max = 125.e-12
29double precision, 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-20
37
38real, parameter :: r3ccnmin = 1.e-21    ! ie rccnmin = 0.1 microns
39real, parameter :: rccnmin  = 0.1e-6
40
41real, parameter :: r3ccnmax = 125.e-12  ! ie rccnmax  = 500 microns
42real, parameter :: rccnmax  = 500.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
112DOUBLE PRECISION, intent(in)  :: qice,qccn,nccn
113real, intent(in)  :: coeff         ! this coeff is tauscaling if microphy = T (possibly ccn_factor^-1 otherwise)
114real, intent(out) :: rhocloudco2 ! rhocloud is needed for sedimentation and is also a good diagnostic variable
115double precision, intent(out) :: rice
116real nccn_true,qccn_true ! nombre et masse de CCN
117   
118nccn_true = max(nccn * coeff, 1.e-30)
119qccn_true = max(qccn * coeff, 1.e-30)
120
121
122  rhocloudco2 = (qice *rho_ice_co2 + qccn_true*rho_dust) / (qice + qccn_true)
123
124  rhocloudco2 = min(max(rhocloudco2,rho_ice_co2),rho_dust)
125
126  rice = (qice + qccn_true) * 0.75 / pi / rhocloudco2 / nccn_true
127
128  if (rice .le. r3iceco2min) then !r3icemin radius power 3 ?
129    rice = riceco2min
130  else if (rice .ge. r3iceco2max) then !r3icemin radius power 3 ?
131    rice = riceco2max
132  else
133    rice = rice**(1./3.) ! here rice is always positive
134  endif
135
136
137end subroutine updaterice_microco2
138
139
140
141!============================================================================
142!============================================================================
143!============================================================================
144! Update ice radius from a typical profile if microphys == false
145subroutine updaterice_typ(qice,tau,pzlay,rice)
146use tracer_mod, only: rho_ice
147USE comcstfi_h
148implicit none
149
150real, intent(in)  :: qice
151real, intent(in)  :: tau   ! tau for dust
152real, intent(in)  :: pzlay ! altitude at the middle of the layers
153real, intent(out) :: rice
154real rccn,nccn ! radius  and number of ice crystals
155
156!   Typical CCN profile following Montmessin et al. 2004
157!   (N0=2e6 m-3 has been converted into N0=1.3e8 kg-1, otherwise the equation for rice is not homogeneous...)
158     nccn  = 1.3e+8*max(tau,0.001)/0.1*exp(-pzlay/10000.)
159!   The previously used profile was not correct:
160!   Nccn=( epaisseur/masse ) * 2.e+6/0.1*max(tau,0.001)*exp(-pzlay/10000.)
161     
162if (nccn .le. 1) then
163
164  rice = ricemin
165
166else
167
168! Typical dust radius profile:
169  rccn = max(rdust0*exp(-pzlay/18000.),1.e-9)
170  rice  = qice * 0.75 / pi / rho_ice / nccn + rccn*rccn*rccn
171       
172  if (rice .le. r3icemin) then
173    rice = ricemin
174  else if (rice .ge. r3icemax) then
175    rice = ricemax
176  else
177    rice = rice**(1./3.) ! here rice is always positive
178  endif
179
180endif
181 
182end subroutine updaterice_typ
183!============================================================================
184!============================================================================
185!============================================================================
186
187
188
189
190
191
192!============================================================================
193!============================================================================
194!============================================================================
195! This subroutine computes the geometric mean radius(or number median radius)
196! For a lognormal distribution :
197! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
198! To be used with doubleq == true. otherwise, rdust is constant !!!
199subroutine updaterdust(qdust,ndust,rdust,tauscaling)
200use tracer_mod, only: r3n_q
201USE comcstfi_h
202implicit none
203
204real, intent(in) :: qdust,ndust ! needed if doubleq
205real, intent(in), optional :: tauscaling ! useful for realistic thresholds
206real, intent(out) :: rdust
207real coeff
208
209if (present(tauscaling)) then
210  coeff = tauscaling ! thresholds on realistic values
211else
212  coeff = 1. ! thresholds on virtual values
213endif
214
215if ((ndust .le. ndust_threshold/coeff) .or. (qdust .le. qdust_threshold/coeff)) then
216
217   rdust = rdustmin
218
219else
220
221   rdust  = r3n_q * qdust / ndust
222
223   if (rdust .le. r3dustmin) then
224       rdust = rdustmin
225   else if (rdust .ge. r3dustmax) then
226       rdust = rdustmax
227   else
228      rdust = rdust**(1./3.)
229   endif
230
231endif
232     
233
234end subroutine updaterdust
235!============================================================================
236!============================================================================
237!============================================================================
238
239
240
241
242
243
244!============================================================================
245!============================================================================
246!============================================================================
247! This subroutine computes the mass mean radius,
248! used for heterogenous nucleation on CCNs in microphysics.
249! For a lognormal distribution :
250! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
251subroutine updaterccn(qccn,nccn,rccn,tauscaling)
252use tracer_mod, only: rho_dust
253USE comcstfi_h
254implicit none
255
256real, intent(in) :: qccn,nccn ! needed if doubleq
257real, intent(in), optional :: tauscaling ! useful for realistic thresholds
258
259real, intent(out) :: rccn
260
261real coeff
262
263
264if (present(tauscaling)) then
265  coeff = tauscaling ! threshold on realistic values
266else
267  coeff = 1. ! threshold on virtual values
268endif
269
270if ((nccn .le. nccn_threshold/coeff) .or. (qccn .le. qccn_threshold/coeff)) then
271
272   rccn = rccnmin
273
274else
275
276   rccn  =  qccn * 0.75 / pi / rho_dust / nccn
277
278   if (rccn .le. r3ccnmin) then
279       rccn = rccnmin
280   else if (rccn .ge. r3ccnmax) then
281       rccn = rccnmax
282   else
283       rccn = rccn**(1./3.)
284   endif
285
286endif
287     
288end subroutine updaterccn
289!============================================================================
290!============================================================================
291!============================================================================
292
293
294
295end module updaterad
296
Note: See TracBrowser for help on using the repository browser.