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

Last change on this file since 1266 was 1266, checked in by aslmd, 11 years ago

LMDZ.MARS
IMPORTANT CHANGE

  • Remove all reference/use of nlayermx and dimphys.h
  • Made use of automatic arrays whenever arrays are needed with dimension nlayer
  • Remove lots of obsolete reference to dimensions.h
  • Converted iono.h and param_v4.h into corresponding modules

(with embedded subroutine to allocate arrays)
(no arrays allocated if thermosphere not used)

  • Deleted param.h and put contents into module param_v4_h
  • Adapted testphys1d, newstart, etc...
  • Made DATA arrays in param_read to be initialized by subroutine

fill_data_thermos in module param_v4_h

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