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

Last change on this file since 1036 was 1036, checked in by emillour, 11 years ago

Mars GCM: (a first step towards using parallel dynamics)

  • IMPORTANT CHANGE: Implemented dynamic tracers. It is no longer necessary to compile the model with the '-t #' option, number of tracers is simply read from tracer.def file (as before). Adapted makegcm_* scripts (and co.) accordingly. Technical aspects of the switch to dynamic tracers are:
    • advtrac.h (in dyn3d) removed and replaced by module infotrac.F
    • tracer.h (in phymars) removed and replaced by module tracer_mod.F90 (which contains nqmx, the number of tracers, etc. and can be used anywhere in the physics).
  • Included some side cleanups: removed unused files (in dyn3d) anldoppler2.F, anl_mcdstats.F and anl_stats-diag.F, and all the unecessary dimensions.* files in grid/dimension.
  • Checked that changes are clean and that GCM yields identical results (in debug mode) to previous svn version.

EM

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