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

Last change on this file since 2156 was 1996, checked in by emillour, 6 years ago

Mars physics:

  • Turn watersat into a module.

CO2 cloud updates:

  • compute co2 condensation tendencies in the co2 cloud scheme and pass them on to vdifc (for tests; they might not be needed) and adapt newcondens.

DB+EM

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