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

Last change on this file since 2514 was 2494, checked in by cmathe, 4 years ago

Mars GCM:
co2_ice as scatterer in radiative transfert. Need co2clouds and

activeco2ice .eqv. True. Files involved:

  • aeropacity_mod.F
  • callradite_mod.F
  • physiq_mod.F
  • updatereffrad_mod.F
  • suaer.F90
  • determine co2_ice density from temperature. Used in riceco2 computation.

Files involved:

  • co2cloud.F90
  • improvedco2clouds_mod.F90
  • updaterad.F90
  • updatereffrad_mod.F
  • co2condens_mod4micro.F: variable initialization
  • initracer.F: add nuiceco2_ref = 0.2
  • phyredem.F: remove co2_ice from qsurf since co2_ice => co2ice
  • watercloud_mod.F: tiny typo

CM

File size: 10.7 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, temperature, coeff, rice, rhocloudco2)
103use tracer_mod, only: rho_dust
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
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 ! 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
121call density_co2_ice(dble(temperature), rho_ice_co2T)
122
123rhocloudco2 = (qice * rho_ice_co2T + qccn_true*rho_dust) / (qice + qccn_true)
124
125rhocloudco2 = min(max(rhocloudco2,rho_ice_co2T), rho_dust)
126
127rice = (qice + qccn_true) * 0.75 / pi / rhocloudco2 / nccn_true
128
129rhocloudco2 = (qice * rho_ice_co2T + qccn_true*rho_dust) / (qice + qccn_true)
130
131if (rice .le. r3iceco2min) then !r3icemin radius power 3 ?
132  rice = riceco2min
133else if (rice .ge. r3iceco2max) then !r3icemin radius power 3 ?
134  rice = riceco2max
135else
136  rice = rice**(1./3.) ! here rice is always positive
137endif
138
139
140end subroutine updaterice_microco2
141
142
143
144!============================================================================
145!============================================================================
146!============================================================================
147! Update ice radius from a typical profile if microphys == false
148subroutine updaterice_typ(qice,tau,pzlay,rice)
149use tracer_mod, only: rho_ice
150USE comcstfi_h
151implicit none
152
153real, intent(in)  :: qice
154real, intent(in)  :: tau   ! tau for dust
155real, intent(in)  :: pzlay ! altitude at the middle of the layers
156real, intent(out) :: rice
157real rccn,nccn ! radius  and number of ice crystals
158
159!   Typical CCN profile following Montmessin et al. 2004
160!   (N0=2e6 m-3 has been converted into N0=1.3e8 kg-1, otherwise the equation for rice is not homogeneous...)
161     nccn  = 1.3e+8*max(tau,0.001)/0.1*exp(-pzlay/10000.)
162!   The previously used profile was not correct:
163!   Nccn=( epaisseur/masse ) * 2.e+6/0.1*max(tau,0.001)*exp(-pzlay/10000.)
164     
165if (nccn .le. 1) then
166
167  rice = ricemin
168
169else
170
171! Typical dust radius profile:
172  rccn = max(rdust0*exp(-pzlay/18000.),1.e-9)
173  rice  = qice * 0.75 / pi / rho_ice / nccn + rccn*rccn*rccn
174       
175  if (rice .le. r3icemin) then
176    rice = ricemin
177  else if (rice .ge. r3icemax) then
178    rice = ricemax
179  else
180    rice = rice**(1./3.) ! here rice is always positive
181  endif
182
183endif
184 
185end subroutine updaterice_typ
186!============================================================================
187!============================================================================
188!============================================================================
189
190
191
192
193
194
195!============================================================================
196!============================================================================
197!============================================================================
198! This subroutine computes the geometric mean radius(or number median radius)
199! For a lognormal distribution :
200! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
201! To be used with doubleq == true. otherwise, rdust is constant !!!
202subroutine updaterdust(qdust,ndust,rdust,tauscaling)
203use tracer_mod, only: r3n_q
204USE comcstfi_h
205implicit none
206
207real, intent(in) :: qdust,ndust ! needed if doubleq
208real, intent(in), optional :: tauscaling ! useful for realistic thresholds
209real, intent(out) :: rdust
210real coeff
211
212if (present(tauscaling)) then
213  coeff = tauscaling ! thresholds on realistic values
214else
215  coeff = 1. ! thresholds on virtual values
216endif
217
218if ((ndust .le. ndust_threshold/coeff) .or. (qdust .le. qdust_threshold/coeff)) then
219
220   rdust = rdustmin
221
222else
223
224   rdust  = r3n_q * qdust / ndust
225
226   if (rdust .le. r3dustmin) then
227       rdust = rdustmin
228   else if (rdust .ge. r3dustmax) then
229       rdust = rdustmax
230   else
231      rdust = rdust**(1./3.)
232   endif
233
234endif
235     
236
237end subroutine updaterdust
238!============================================================================
239!============================================================================
240!============================================================================
241
242
243
244
245
246
247!============================================================================
248!============================================================================
249!============================================================================
250! This subroutine computes the mass mean radius,
251! used for heterogenous nucleation on CCNs in microphysics.
252! For a lognormal distribution :
253! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
254subroutine updaterccn(qccn,nccn,rccn,tauscaling)
255use tracer_mod, only: rho_dust
256USE comcstfi_h
257implicit none
258
259real, intent(in) :: qccn,nccn ! needed if doubleq
260real, intent(in), optional :: tauscaling ! useful for realistic thresholds
261
262real, intent(out) :: rccn
263
264real coeff
265
266
267if (present(tauscaling)) then
268  coeff = tauscaling ! threshold on realistic values
269else
270  coeff = 1. ! threshold on virtual values
271endif
272
273if ((nccn .le. nccn_threshold/coeff) .or. (qccn .le. qccn_threshold/coeff)) then
274
275   rccn = rccnmin
276
277else
278
279   rccn  =  qccn * 0.75 / pi / rho_dust / nccn
280
281   if (rccn .le. r3ccnmin) then
282       rccn = rccnmin
283   else if (rccn .ge. r3ccnmax) then
284       rccn = rccnmax
285   else
286       rccn = rccn**(1./3.)
287   endif
288
289endif
290     
291end subroutine updaterccn
292!============================================================================
293!============================================================================
294!============================================================================
295
296
297
298
299
300!============================================================================
301!============================================================================
302!============================================================================
303! This subroutine computes the mass mean radius,
304! used for heterogenous nucleation on CCNs in microphysics of CO2.
305! For a lognormal distribution :
306! geometric mean radius = mass mean radius x exp(-1.5 sigma0^2)
307subroutine updaterccnCO2(qccn,nccn,rccn,tauscaling)
308use tracer_mod, only: rho_dust
309USE comcstfi_h
310implicit none
311
312real, intent(in) :: qccn,nccn ! needed if doubleq
313real, intent(in), optional :: tauscaling ! useful for realistic thresholds
314
315real, intent(out) :: rccn
316
317real coeff
318
319
320if (present(tauscaling)) then
321  coeff = tauscaling ! threshold on realistic values
322else
323  coeff = 1. ! threshold on virtual values
324endif
325
326if ((nccn .le. nccn_threshold/coeff) .or. (qccn .le. qccn_threshold/coeff)) then
327
328   rccn = rccnCO2min
329
330else
331
332   rccn  =  qccn * 0.75 / pi / rho_dust / nccn
333   rccn = rccn**(1./3.)
334
335endif
336
337rccn=min(5.E-4,rccn)
338     
339end subroutine updaterccnCO2
340!============================================================================
341!============================================================================
342!============================================================================
343
344
345
346end module updaterad
347
Note: See TracBrowser for help on using the repository browser.