source: trunk/LMDZ.VENUS/libf/phyvenus/chemparam_mod.F90 @ 3722

Last change on this file since 3722 was 3689, checked in by flefevre, 9 months ago

Implementation of SO dimer photochemistry. From the work of Joanna Egan et al. (GRL, e2024GL113090, 2015). Many thanks Jo!
Three extra species: osso_cis, osso_trans, s2o2_cyc
One less species: s2o2

Extra OSSO UV cross_section file: OSSO_cross_sections_formatted.txt (to be placed in INPUT/cross_sections)
download link : https://owncloud.latmos.ipsl.fr/index.php/s/7kmP3ANi9o37PoZ

  • Property svn:executable set to *
File size: 65.9 KB
Line 
1!============================================================================
2
3module chemparam_mod
4
5!============================================================================
6
7! 1) cloud microphysical parameters for the simplified scheme (cl_scheme = 1)
8! 2) indexes and molecular mass of chemical species
9
10implicit none
11
12!----------------------------------------------------------------------------
13!     chemical tracers
14!----------------------------------------------------------------------------
15
16integer, save :: i_co2, i_co, i_h2, i_h2o, i_o1d,        &
17                 i_o, i_o2, i_o2dg, i_o3, i_h,           &
18                 i_oh, i_ho2, i_h2o2, i_cl, i_clo,       &
19                 i_cl2, i_hcl, i_hocl, i_clco, i_clco3,  &
20                 i_cocl2, i_s, i_so, i_so2, i_so3,       &
21                 i_osso_cis, i_osso_trans, i_s2o2_cyc,   &
22                 i_ocs, i_hso3, i_h2so4, i_s2,           &
23                 i_clso2, i_oscl, i_n2, i_he, i_n, i_no, &
24                 i_no2, i_n2d,                           &
25                 i_co2plus, i_coplus, i_oplus, i_o2plus, &
26                 i_n2plus, i_hplus, i_h2oplus, i_nplus,  &
27                 i_ohplus, i_cplus, i_noplus, i_h3oplus, &
28                 i_hcoplus, i_hco2plus, i_elec
29
30integer, save :: i_h2oliq, i_h2so4liq
31
32integer, save :: i_m0_aer, i_m3_aer,                       &
33                 i_m0_mode1drop, i_m0_mode1ccn,            &
34                 i_m3_mode1sa, i_m3_mode1w, i_m3_mode1ccn, &
35                 i_m0_mode2drop, i_m0_mode2ccn,            &
36                 i_m3_mode2sa, i_m3_mode2w, i_m3_mode2ccn
37
38integer, save :: nmicro  ! number of species in the liquid phase
39
40real, dimension(:), save, allocatable :: m_tr           ! molecular mass of tracers
41real, dimension(:), save, allocatable :: type_tr        ! type of tracer
42
43real, dimension(:,:), save, allocatable :: no_emission
44real, dimension(:,:), save, allocatable :: o2_emission
45
46!----------------------------------------------------------------------------
47!     cloud parameters
48!----------------------------------------------------------------------------
49
50integer, save :: cloudmin, cloudmax
51
52!     qrad : ratio radius shell model of mode 3 (cimino, icarus, 1982)
53!     if qrad = 0, fully liquid, if qrad = 1 fully solid
54
55real, save :: qrad
56
57!     median radius and standard deviation in each mode
58
59real, save, dimension(:,:,:), allocatable :: r_median, stddev
60
61!     k_mass : defines how the condensed phase is distributed in each mode.
62!              sum of k_mass = 1
63
64real, save, dimension(:,:,:), allocatable :: k_mass
65
66real, save, dimension(:,:,:), allocatable :: nbrtot
67real, save, dimension(:,:), allocatable :: wh2so4
68real, save, dimension(:,:), allocatable :: rho_droplet
69
70contains
71
72!============================================================================
73
74subroutine cloud_ini(nbr_lon, nbr_lev, nbr_mode)
75
76!============================================================================
77
78!     sets cloud microphysical parameters for each mode:
79!     radius, standard deviation, mass distribution
80
81integer :: nbr_lon, nbr_lev, nbr_mode
82integer :: i_lev, ilon
83
84allocate(nbrtot(nbr_lon,nbr_lev,nbr_mode))
85allocate(r_median(nbr_lon,nbr_lev,nbr_mode))
86allocate(k_mass(nbr_lon,nbr_lev,nbr_mode))
87allocate(stddev(nbr_lon,nbr_lev,nbr_mode))
88allocate(wh2so4(nbr_lon,nbr_lev))
89allocate(rho_droplet(nbr_lon,nbr_lev))
90
91!     initialisation
92
93r_median(:,:,:)  = 0.     ! median radius
94stddev(:,:,:)    = 0.     ! geometric std deviation
95k_mass(:,:,:)    = 0.     ! coeff mass multimodal
96nbrtot(:,:,:)    = 0.
97wh2so4(:,:)      = 0.
98rho_droplet(:,:) = 0.
99
100!     minimum and maximum levels for the clouds
101
102cloudmin = 20   ! 20: 38 km
103cloudmax = 50   ! 50: 95 km
104
105print*,'================================'
106print*,'start initialisation cloud layer'
107print*,'================================'
108
109!       ===============================================
110!       knollenberg & hunten, 1980 and james et al 1997
111!       ===============================================
112!       initialisation unimodale
113!       ===============================================
114
115!     lower haze: mode 1
116!      do i_lev=cloudmin,20
117!      r_median(:,i_lev,1)=0.2e-6
118!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
119!      stddev(:,i_lev,1)=1.56
120!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
121!      k_mass(:,i_lev,1)=1.0
122!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
123!      end do
124
125!     lower cloud: mode 3
126!      do i_lev=21,23
127!      r_median(:,i_lev,1)=3.65e-6
128!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
129!      stddev(:,i_lev,1)=1.28
130!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
131!      k_mass(:,i_lev,1)=1.0
132!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
133!      end do
134
135!     middle cloud: mode 2 prime
136!      do i_lev=24,28
137!      r_median(:,i_lev,1)=1.4e-6
138!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
139!      stddev(:,i_lev,1)=1.23
140!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
141!      k_mass(:,i_lev,1)=1.0
142!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
143!      end do
144
145!     upper cloud: mode 2
146!      do i_lev=29,35
147!      r_median(:,i_lev,1)=1.0e-6
148!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
149!      stddev(:,i_lev,1)=1.29
150!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
151!      k_mass(:,i_lev,1)=1.0
152!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
153!      end do
154
155!     upper haze: mode 1
156!      do i_lev=36, cloudmax
157!      r_median(:,i_lev,1)=0.2e-6
158!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
159!      stddev(:,i_lev,1)=2.16
160!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
161!      k_mass(:,i_lev,1)=1.0
162!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
163!      end do
164
165!       ===============================================
166!       initialisation trimodale
167!       ===============================================
168
169!     lower haze: mode 1
170!      do i_lev=cloudmin,20
171!      r_median(:,i_lev,1)=0.3e-6
172!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
173!      stddev(:,i_lev,1)=1.56
174!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
175!      k_mass(:,i_lev,1)=1.0
176!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
177!      end do
178
179!     lower haze: mode 2
180   !   do i_lev=cloudmin,20
181   !   r_median(:,i_lev,2)=1.4e-6
182   !   print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
183   !   stddev(:,i_lev,2)=1.23
184   !   print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
185   !   k_mass(:,i_lev,2)=0.0
186   !   print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
187   !   end do
188
189!     lower haze: mode 3
190!      do i_lev=cloudmin,20
191!      r_median(:,i_lev,3)=3.65e-6
192!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
193!      stddev(:,i_lev,3)=1.28
194!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
195!      k_mass(:,i_lev,3)=0.
196!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
197!      end do
198
199!     lower cloud: mode 1
200!      do i_lev=21,23
201!      r_median(:,i_lev,1)=0.3e-6
202!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
203!      stddev(:,i_lev,1)=1.56
204!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
205!      k_mass(:,i_lev,1)=0.1
206!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
207!      end do
208
209!     lower cloud: mode 2 prime
210!      do i_lev=21,23
211!      r_median(:,i_lev,2)=1.4e-6
212!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
213!      stddev(:,i_lev,2)=1.23
214!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
215!      k_mass(:,i_lev,2)=0.4
216!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
217!      end do
218
219!     lower cloud: mode 3
220!      do i_lev=21,23
221!      r_median(:,i_lev,3)=3.65e-6
222!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
223!      stddev(:,i_lev,3)=1.28
224!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
225!      k_mass(:,i_lev,3)=0.5
226!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
227!      end do
228
229!     middle cloud: mode 1
230!      do i_lev=24,28
231!      r_median(:,i_lev,1)=0.3e-6
232!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
233!      stddev(:,i_lev,1)=1.56
234!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
235!      k_mass(:,i_lev,1)=0.0
236!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
237!      end do
238
239!     middle cloud: mode 2 prime
240!      do i_lev=24,28
241!      r_median(:,i_lev,2)=1.4e-6
242!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
243!      stddev(:,i_lev,2)=1.23
244!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
245!      k_mass(:,i_lev,2)=0.8
246!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
247!      end do
248
249!     middle cloud: mode 3
250!      do i_lev=24,28
251!      r_median(:,i_lev,3)=3.65e-6
252!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
253!      stddev(:,i_lev,3)=1.28
254!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
255!      k_mass(:,i_lev,3)=0.2
256!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
257!      end do
258
259
260!     upper cloud: mode 1
261!      do i_lev=29,35
262!      r_median(:,i_lev,1)=0.3e-6
263!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
264!      stddev(:,i_lev,1)=1.56
265!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
266!      k_mass(:,i_lev,1)=0.15
267!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
268!      end do
269
270!     upper cloud: mode 2
271!      do i_lev=29,35
272!      r_median(:,i_lev,2)=1.0e-6
273!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
274!      stddev(:,i_lev,2)=1.29
275!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
276!      k_mass(:,i_lev,2)=0.85
277!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
278!      end do
279
280!     upper cloud: mode 3
281!      do i_lev=29,35
282!      r_median(:,i_lev,3)=3.65e-6
283!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
284!      stddev(:,i_lev,3)=1.28
285!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
286!      k_mass(:,i_lev,3)=0.0
287!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
288!      end do
289
290!     upper haze: mode 1
291!      do i_lev=36, cloudmax
292!      r_median(:,i_lev,1)=0.3e-6
293!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
294!      stddev(:,i_lev,1)=1.56
295!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
296!      k_mass(:,i_lev,1)=1.0
297!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
298!      end do
299
300!     upper haze: mode 2
301!      do i_lev=36, cloudmax
302!      r_median(:,i_lev,2)=1.e-6
303!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
304!      stddev(:,i_lev,2)=1.29
305!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
306!      k_mass(:,i_lev,2)=0.0
307!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
308!      end do
309
310!     upper haze: mode 3
311!      do i_lev=36, cloudmax
312!      r_median(:,i_lev,3)=3.65e-6
313!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
314!      stddev(:,i_lev,3)=2.16
315!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
316!      k_mass(:,i_lev,3)=0.0
317!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
318!      end do
319!=============================================================
320
321!       ===============================================
322!       initialisation trimodale knollenberg
323!       ===============================================
324
325!     lower haze: mode 1
326      ! do i_lev=cloudmin,22
327      ! r_median(:,i_lev,1)=0.1e-6
328      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
329      ! stddev(:,i_lev,1)=1.57
330      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
331      ! k_mass(:,i_lev,1)=1.0
332      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
333      ! end do
334
335!     lower haze: mode 2
336      ! do i_lev=cloudmin,22
337      ! r_median(:,i_lev,2)=1.4e-6
338      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
339      ! stddev(:,i_lev,2)=1.23
340      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
341      ! k_mass(:,i_lev,2)=0.0
342      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
343      ! end do
344
345!     lower haze: mode 3
346      ! do i_lev=cloudmin,22
347      ! r_median(:,i_lev,3)=3.65e-6
348      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
349      ! stddev(:,i_lev,3)=1.28
350      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
351      ! k_mass(:,i_lev,3)=0.0
352      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
353      ! end do
354
355!     pre cloud: mode 1
356      ! do i_lev=23,23
357      ! r_median(:,i_lev,1)=0.15e-6
358      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
359      ! stddev(:,i_lev,1)=1.8
360      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
361      ! k_mass(:,i_lev,1)=0.04
362      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
363      ! end do
364
365!     pre cloud: mode 2
366      ! do i_lev=23,23
367      ! r_median(:,i_lev,2)=1.0e-6
368      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
369      ! stddev(:,i_lev,2)=1.29
370      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
371      ! k_mass(:,i_lev,2)=0.96
372      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
373      ! end do
374
375!     pre cloud: mode 3
376      ! do i_lev=23,23
377      ! r_median(:,i_lev,3)=3.65e-6
378      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
379      ! stddev(:,i_lev,3)=1.28
380      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
381      ! k_mass(:,i_lev,3)=0.0
382      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
383      ! end do
384
385!      lower cloud: mode 1
386      ! do i_lev=24,24
387      ! r_median(:,i_lev,1)=0.2e-6
388      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
389      ! stddev(:,i_lev,1)=1.8
390      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
391      ! k_mass(:,i_lev,1)=0.014
392      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
393      ! end do
394
395!     lower cloud: mode 2
396      ! do i_lev=24,24
397      ! r_median(:,i_lev,2)=1.0e-6
398      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
399      ! stddev(:,i_lev,2)=1.29
400      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
401      ! k_mass(:,i_lev,2)=0.02
402      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
403      ! end do
404
405!     lower cloud: mode 3
406      ! do i_lev=24,24
407      ! r_median(:,i_lev,3)=3.65e-6
408      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
409      ! stddev(:,i_lev,3)=1.28
410      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
411      ! k_mass(:,i_lev,3)=0.966
412      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
413      ! end do
414
415!     middle cloud: mode 1
416      ! do i_lev=25,28
417      ! r_median(:,i_lev,1)=0.15e-6
418      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
419      ! stddev(:,i_lev,1)=1.9
420      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
421      ! k_mass(:,i_lev,1)=0.0084
422      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
423      ! end do
424
425!     middle cloud: mode 2 prime
426      ! do i_lev=25,28
427      ! r_median(:,i_lev,2)=1.4e-6
428      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
429      ! stddev(:,i_lev,2)=1.23
430      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
431      ! k_mass(:,i_lev,2)=0.21
432      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
433      ! end do
434
435!     middle cloud: mode 3
436      ! do i_lev=25,28
437      ! r_median(:,i_lev,3)=3.65e-6
438      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
439      ! stddev(:,i_lev,3)=1.28
440      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
441      ! k_mass(:,i_lev,3)=0.7816
442      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
443      ! end do
444
445!      option: upper haze remplacee par extension upper cloud
446!         => 35 remplace par cloudmax et upper haze commentee
447!       ===============================================
448
449!     upper cloud: mode 1
450      ! do i_lev=29,35 !cloudmax
451      ! r_median(:,i_lev,1)=0.2e-6
452      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
453      ! stddev(:,i_lev,1)=2.16
454      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
455      ! k_mass(:,i_lev,1)=0.72
456      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
457      ! end do
458
459!     upper cloud: mode 2
460      ! do i_lev=29,35 !cloudmax
461      ! r_median(:,i_lev,2)=1.0e-6
462      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
463      ! stddev(:,i_lev,2)=1.29
464      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
465      ! k_mass(:,i_lev,2)=0.28
466      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
467      ! end do
468
469!     upper cloud: mode 3
470      ! do i_lev=29,35 !cloudmax
471      ! r_median(:,i_lev,3)=3.65e-6
472      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
473      ! stddev(:,i_lev,3)=1.28
474      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
475      ! k_mass(:,i_lev,3)=0.0
476      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
477      ! end do
478
479!     upper haze: mode 1
480      ! do i_lev=36, cloudmax
481      ! r_median(:,i_lev,1)=0.2e-6
482      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
483      ! stddev(:,i_lev,1)=2.16
484      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
485      ! k_mass(:,i_lev,1)=1.0
486      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
487      ! end do
488
489!     upper haze: mode 2
490      ! do i_lev=36, cloudmax
491      ! r_median(:,i_lev,2)=1.e-6
492      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
493      ! stddev(:,i_lev,2)=1.29
494      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
495      ! k_mass(:,i_lev,2)=0.0
496      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
497      ! end do
498
499!     upper haze: mode 3
500      ! do i_lev=36, cloudmax
501      ! r_median(:,i_lev,3)=3.65e-6
502      ! print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
503      ! stddev(:,i_lev,3)=2.16
504      ! print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
505      ! k_mass(:,i_lev,3)=0.0
506      ! print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
507      ! end do
508
509!=============================================================
510
511!       ===============================================================
512!       initialisation trimodale "knollenberg" sans mode3, mode2 etendu
513!       ===============================================================
514
515!     lower haze: mode 1
516!      do i_lev=cloudmin,22
517!      r_median(:,i_lev,1)=0.1e-6
518!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
519!      stddev(:,i_lev,1)=1.57
520!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
521!      k_mass(:,i_lev,1)=1.0
522!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
523!      end do
524
525!     lower haze: mode 2
526!      do i_lev=cloudmin,22
527!      r_median(:,i_lev,2)=1.4e-6
528!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
529!      stddev(:,i_lev,2)=1.23
530!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
531!      k_mass(:,i_lev,2)=0.0
532!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
533!      end do
534
535!     lower haze: mode 3
536!      do i_lev=cloudmin,22
537!      r_median(:,i_lev,3)=3.65e-6
538!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
539!      stddev(:,i_lev,3)=1.28
540!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
541!      k_mass(:,i_lev,3)=0.0
542!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
543!      end do
544
545!     pre cloud: mode 1
546!      do i_lev=23,23
547!      r_median(:,i_lev,1)=0.15e-6
548!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
549!      stddev(:,i_lev,1)=1.8
550!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
551!      k_mass(:,i_lev,1)=0.04
552!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
553!      end do
554
555!     pre cloud: mode 2
556!      do i_lev=23,23
557!      r_median(:,i_lev,2)=1.0e-6
558!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
559!      stddev(:,i_lev,2)=1.29
560!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
561!      k_mass(:,i_lev,2)=0.96
562!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
563!      end do
564
565!     pre cloud: mode 3
566!      do i_lev=23,23
567!      r_median(:,i_lev,3)=3.65e-6
568!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
569!      stddev(:,i_lev,3)=1.28
570!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
571!      k_mass(:,i_lev,3)=0.0
572!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
573!      end do
574
575!      lower cloud: mode 1
576!      do i_lev=24,24
577!      r_median(:,i_lev,1)=0.2e-6
578!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
579!      stddev(:,i_lev,1)=1.8
580!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
581!      k_mass(:,i_lev,1)=0.014
582!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
583!      end do
584
585!     lower cloud: mode 2
586!      do i_lev=24,24
587!      r_median(:,i_lev,2)=1.0e-6
588!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
589!      stddev(:,i_lev,2)=1.6
590!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
591!      k_mass(:,i_lev,2)=0.986
592!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
593!      end do
594
595!     lower cloud: mode 3
596!      do i_lev=24,24
597!      r_median(:,i_lev,3)=3.65e-6
598!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
599!      stddev(:,i_lev,3)=1.28
600!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
601!      k_mass(:,i_lev,3)=0.
602!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
603!      end do
604
605!     middle cloud: mode 1
606!      do i_lev=25,28
607!      r_median(:,i_lev,1)=0.15e-6
608!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
609!      stddev(:,i_lev,1)=1.9
610!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
611!      k_mass(:,i_lev,1)=0.0084
612!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
613!      end do
614
615!     middle cloud: mode 2 prime
616!      do i_lev=25,28
617!      r_median(:,i_lev,2)=1.4e-6
618!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
619!      stddev(:,i_lev,2)=1.6
620!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
621!      k_mass(:,i_lev,2)=0.9916
622!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
623!      end do
624
625!     middle cloud: mode 3
626!      do i_lev=25,28
627!      r_median(:,i_lev,3)=3.65e-6
628!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
629!      stddev(:,i_lev,3)=1.28
630!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
631!      k_mass(:,i_lev,3)=0.0
632!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
633!      end do
634
635
636!     upper cloud: mode 1
637!      do i_lev=29,35
638!      r_median(:,i_lev,1)=0.2e-6
639!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
640!      stddev(:,i_lev,1)=2.16
641!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
642!      k_mass(:,i_lev,1)=0.72
643!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
644!      end do
645
646!     upper cloud: mode 2
647!      do i_lev=29,35
648!      r_median(:,i_lev,2)=1.0e-6
649!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
650!      stddev(:,i_lev,2)=1.29
651!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
652!      k_mass(:,i_lev,2)=0.28
653!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
654!      end do
655
656!     upper cloud: mode 3
657!      do i_lev=29,35
658!      r_median(:,i_lev,3)=3.65e-6
659!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
660!      stddev(:,i_lev,3)=1.28
661!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
662!      k_mass(:,i_lev,3)=0.0
663!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
664!      end do
665
666!     upper haze: mode 1
667!      do i_lev=36, cloudmax
668!      r_median(:,i_lev,1)=0.2e-6
669!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
670!      stddev(:,i_lev,1)=2.16
671!      print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
672!      k_mass(:,i_lev,1)=1.0
673!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
674!      end do
675
676!     upper haze: mode 2
677!      do i_lev=36, cloudmax
678!      r_median(:,i_lev,2)=1.e-6
679!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
680!      stddev(:,i_lev,2)=1.29
681!      print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
682!      k_mass(:,i_lev,2)=0.0
683!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
684!      end do
685
686!     upper haze: mode 3
687!      do i_lev=36, cloudmax
688!      r_median(:,i_lev,3)=3.65e-6
689!      print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
690!      stddev(:,i_lev,3)=2.16
691!      print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
692!      k_mass(:,i_lev,3)=0.0
693!      print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
694!      end do
695
696   ! ========================================================
697   ! initialisation bimodale k&h 1980 with mode 3 fully solid
698   ! ========================================================
699   !    ! mode 3 fully solid
700   !    qrad=1
701   !    ! normally nb_mode=2 in physiq.def !!!
702   !    do ilon=1,nbr_lon
703   ! !     mode 1
704   !       do i_lev=cloudmin,20
705   !          r_median(ilon,i_lev,1)=0.125e-6
706   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
707   !          stddev(ilon,i_lev,1)=1.57
708   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
709   !          k_mass(ilon,i_lev,1)=1.0
710   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
711   !       end do
712
713   !       r_median(ilon,21,1)=0.125e-6
714   !       print*,'level',21,'r r_median',r_median(1,21,1)
715   !       stddev(ilon,21,1)=1.57
716   !       print*,'level',21,'dev std',stddev(1,21,1)
717   !       k_mass(ilon,21,1)=0.02
718   !       print*,'level',21,'coeff mass: k_mass',k_mass(1,21,1)
719
720   !       r_median(ilon,22,1)=0.2e-6
721   !       print*,'level',22,'r r_median',r_median(1,22,1)
722   !       stddev(ilon,22,1)=1.8
723   !       print*,'level',22,'dev std',stddev(1,22,1)
724   !       k_mass(ilon,22,1)=0.02
725   !       print*,'level',22,'coeff mass: k_mass',k_mass(1,22,1)
726
727   !       r_median(ilon,23,1)=0.15e-6
728   !       print*,'level',23,'r r_median',r_median(1,23,1)
729   !       stddev(ilon,23,1)=1.8
730   !       print*,'level',23,'dev std',stddev(1,23,1)
731   !       k_mass(ilon,23,1)=0.02
732   !       print*,'level',23,'coeff mass: k_mass',k_mass(1,23,1)
733
734   !       do i_lev=24,25
735   !          r_median(ilon,i_lev,1)=0.15e-6
736   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
737   !          stddev(ilon,i_lev,1)=1.9
738   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
739   !          k_mass(ilon,i_lev,1)=0.02
740   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
741   !       end do
742
743   !       r_median(ilon,26,1)=0.175e-6
744   !       print*,'level',26,'r r_median',r_median(1,26,1)
745   !       stddev(ilon,26,1)=2.16
746   !       print*,'level',26,'dev std',stddev(1,26,1)
747   !       k_mass(ilon,26,1)=0.175
748   !       print*,'level',26,'coeff mass: k_mass',k_mass(1,26,1)
749
750   !       do i_lev=27,33
751   !          r_median(ilon,i_lev,1)=0.175e-6
752   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
753   !          stddev(ilon,i_lev,1)=2.16
754   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
755   !          k_mass(ilon,i_lev,1)=0.25
756   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
757   !       end do
758
759   !       do i_lev=34,cloudmax
760   !          r_median(ilon,i_lev,1)=0.175e-6
761   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
762   !          stddev(ilon,i_lev,1)=2.16
763   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
764   !          k_mass(ilon,i_lev,1)=1.0
765   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
766   !       end do
767
768   !    !     mode 2
769   !       do i_lev=cloudmin,20
770   !          r_median(ilon,i_lev,2)=1.4e-6
771   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
772   !          stddev(ilon,i_lev,2)=1.35
773   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
774   !          k_mass(ilon,i_lev,2)=0.0
775   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
776   !       end do
777
778   !       do i_lev=21,22
779   !          r_median(ilon,i_lev,2)=1.4e-6
780   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
781   !          stddev(ilon,i_lev,2)=1.35
782   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
783   !          k_mass(ilon,i_lev,2)=0.98
784   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
785   !       end do
786
787   !       r_median(ilon,23,2)=1.35e-6
788   !       print*,'level',23,'r r_median',r_median(1,23,2)
789   !       stddev(ilon,23,2)=1.25
790   !       print*,'level',23,'dev std',stddev(1,23,2)
791   !       k_mass(ilon,23,2)=0.98
792   !       print*,'level',23,'coeff mass: k_mass',k_mass(1,23,2)
793
794
795   !       r_median(ilon,24,2)=1.375e-6
796   !       print*,'level',24,'r r_median',r_median(1,24,2)
797   !       stddev(ilon,24,2)=1.2
798   !       print*,'level',24,'dev std',stddev(1,24,2)
799   !       k_mass(ilon,24,2)=0.98
800   !       print*,'level',24,'coeff mass: k_mass',k_mass(1,24,2)
801
802
803   !       r_median(ilon,25,2)=1.4e-6
804   !       print*,'level',25,'r r_median',r_median(1,25,2)
805   !       stddev(ilon,25,2)=1.16
806   !       print*,'level',25,'dev std',stddev(1,25,2)
807   !       k_mass(ilon,25,2)=0.98
808   !       print*,'level',25,'coeff mass: k_mass',k_mass(1,25,2)
809
810
811   !       r_median(ilon,26,2)=1.15e-6
812   !       print*,'level',26,'r r_median',r_median(1,26,2)
813   !       stddev(ilon,26,2)=1.34
814   !       print*,'level',26,'dev std',stddev(1,26,2)
815   !       k_mass(ilon,26,2)=0.825
816   !       print*,'level',26,'coeff mass: k_mass',k_mass(1,26,2)
817
818
819   !       r_median(ilon,27,2)=1.14e-6
820   !       print*,'level',27,'r r_median',r_median(1,27,2)
821   !       stddev(ilon,27,2)=1.33
822   !       print*,'level',27,'dev std',stddev(1,27,2)
823   !       k_mass(ilon,27,2)=0.75
824   !       print*,'level',27,'coeff mass: k_mass',k_mass(1,27,2)
825
826
827   !       r_median(ilon,28,2)=1.35e-6
828   !       print*,'level',28,'r r_median',r_median(1,28,2)
829   !       stddev(ilon,28,2)=1.32
830   !       print*,'level',28,'dev std',stddev(1,28,2)
831   !       k_mass(ilon,28,2)=0.75
832   !       print*,'level',28,'coeff mass: k_mass',k_mass(1,28,2)
833
834
835   !       r_median(ilon,29,2)=1.125e-6
836   !       print*,'level',29,'r r_median',r_median(1,29,2)
837   !       stddev(ilon,29,2)=1.31
838   !       print*,'level',29,'dev std',stddev(1,29,2)
839   !       k_mass(ilon,29,2)=0.75
840   !       print*,'level',29,'coeff mass: k_mass',k_mass(1,29,2)
841
842
843   !       r_median(ilon,30,2)=1.118e-6
844   !       print*,'level',30,'r r_median',r_median(1,30,2)
845   !       stddev(ilon,30,2)=1.30
846   !       print*,'level',30,'dev std',stddev(1,30,2)
847   !       k_mass(ilon,30,2)=0.75
848   !       print*,'level',30,'coeff mass: k_mass',k_mass(1,30,2)
849
850
851   !       r_median(ilon,31,2)=1.11e-6
852   !       print*,'level',31,'r r_median',r_median(1,31,2)
853   !       stddev(ilon,31,2)=1.29
854   !       print*,'level',31,'dev std',stddev(1,31,2)
855   !       k_mass(ilon,31,2)=0.75
856   !       print*,'level',31,'coeff mass: k_mass',k_mass(1,31,2)
857
858   !       do i_lev=32,33
859   !          r_median(ilon,i_lev,2)=1.1e-6
860   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
861   !          stddev(ilon,i_lev,2)=1.28
862   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
863   !          k_mass(ilon,i_lev,2)=0.75
864   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
865   !       end do
866
867   !       do i_lev=34,cloudmax
868   !          r_median(ilon,i_lev,2)=1.1e-6
869   !          print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
870   !          stddev(ilon,i_lev,2)=1.28
871   !          print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
872   !          k_mass(ilon,i_lev,2)=0.0
873   !          print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
874   !       end do
875   !    end do
876
877! ======================================================
878! knollenberg & hunten (1980) with mode 3 97% solid
879! ======================================================
880
881! mode 3 97% solid => qrad=0.97
882
883qrad = 0.97
884
885do ilon = 1,nbr_lon
886
887   ! mode 1
888
889   do i_lev = cloudmin,21
890      r_median(ilon,i_lev,1) = 0.125e-6
891      stddev(ilon,i_lev,1)   = 1.57
892      k_mass(ilon,i_lev,1)   = 1.0
893   end do
894
895   r_median(ilon,22,1)       = 0.125e-6
896   stddev(ilon,22,1)         = 1.57
897   k_mass(ilon,22,1)         = 1.0
898
899   r_median(ilon,23,1)       = 0.2e-6
900   stddev(ilon,23,1)         = 1.8
901   k_mass(ilon,23,1)         = 0.01
902
903   r_median(ilon,24,1)       = 0.15e-6
904   stddev(ilon,24,1)         = 1.8
905   k_mass(ilon,24,1)         = 0.01
906
907   do i_lev = 25,26
908      r_median(ilon,i_lev,1) = 0.15e-6
909      stddev(ilon,i_lev,1)   = 1.9
910      k_mass(ilon,i_lev,1)   = 0.01
911   end do
912
913   r_median(ilon,27,1)       = 0.175e-6
914   stddev(ilon,27,1)         = 2.16
915   k_mass(ilon,27,1)         = 0.175
916
917   do i_lev = 28,34
918      r_median(ilon,i_lev,1) = 0.175e-6
919      stddev(ilon,i_lev,1)   = 2.16
920      k_mass(ilon,i_lev,1)   = 0.25
921   end do
922
923   do i_lev = 35,cloudmax
924      r_median(ilon,i_lev,1) = 0.175e-6
925      stddev(ilon,i_lev,1)   = 2.16
926      k_mass(ilon,i_lev,1)   = 0.25
927   end do
928
929   !  mode 2
930
931   do i_lev = cloudmin,21
932      r_median(ilon,i_lev,2) = 1.4e-6
933      stddev(ilon,i_lev,2)   = 1.35
934      k_mass(ilon,i_lev,2)   = 0.0
935   end do
936
937   r_median(ilon,22,2)       = 1.4e-6
938   stddev(ilon,22,2)         = 1.25
939   k_mass(ilon,22,2)         = 0.0
940
941   r_median(ilon,23,2)       = 1.4e-6
942   stddev(ilon,23,2)         = 1.25
943   k_mass(ilon,23,2)         = 0.35
944
945   r_median(ilon,24,2)       = 1.35e-6
946   stddev(ilon,24,2)         = 1.25
947   k_mass(ilon,24,2)         = 0.35
948
949   r_median(ilon,25,2)       = 1.375e-6
950   stddev(ilon,25,2)         = 1.2
951   k_mass(ilon,25,2)         = 0.35
952
953   r_median(ilon,26,2)       = 1.4e-6
954   stddev(ilon,26,2)         = 1.16
955   k_mass(ilon,26,2)         = 0.35
956
957   r_median(ilon,27,2)       = 1.15e-6
958   stddev(ilon,27,2)         = 1.34
959   k_mass(ilon,27,2)         = 0.825
960
961   r_median(ilon,28,2)       = 1.14e-6
962   stddev(ilon,28,2)         = 1.33
963   k_mass(ilon,28,2)         = 0.75
964
965   r_median(ilon,29,2)       = 1.35e-6
966   stddev(ilon,29,2)         = 1.32
967   k_mass(ilon,29,2)         = 0.75
968
969   r_median(ilon,30,2)       = 1.125e-6
970   stddev(ilon,30,2)         = 1.31
971   k_mass(ilon,30,2)         = 0.75
972
973   r_median(ilon,31,2)       = 1.118e-6
974   stddev(ilon,31,2)         = 1.31
975   k_mass(ilon,31,2)         = 0.75
976
977   r_median(ilon,32,2)       = 1.11e-6
978   stddev(ilon,32,2)         = 1.29
979   k_mass(ilon,32,2)         = 0.75
980
981   do i_lev = 33,34
982      r_median(ilon,i_lev,2) = 1.1e-6
983      stddev(ilon,i_lev,2)   = 1.28
984      k_mass(ilon,i_lev,2)   = 0.75
985   end do
986
987   do i_lev = 35,cloudmax
988      r_median(ilon,i_lev,2) = 1.1e-6
989      stddev(ilon,i_lev,2)   = 1.28
990      k_mass(ilon,i_lev,2)   = 0.75
991   end do
992
993   ! mode 3
994
995   do i_lev = cloudmin,22
996      r_median(ilon,i_lev,3) = 3.65e-6
997      stddev(ilon,i_lev,3)   = 1.28
998      k_mass(ilon,i_lev,3)   = 0.0
999   end do
1000
1001   do i_lev = 23,26
1002      r_median(ilon,i_lev,3) = 3.65e-6
1003      stddev(ilon,i_lev,3)   = 1.28
1004      k_mass(ilon,i_lev,3)   = 0.64
1005   end do
1006
1007   do i_lev = 27,cloudmax
1008      r_median(ilon,i_lev,3) = 3.65e-6
1009      stddev(ilon,i_lev,3)   = 1.28
1010      k_mass(ilon,i_lev,3)   = 0.0
1011   end do
1012end do
1013
1014      ! ! ==================================================================
1015      ! ! initialisation bimodale k&h 1980 with mode 3 0% solid fully liquid
1016      ! ! ==================================================================
1017      ! ! mode 3 0% solid, fully liquid
1018      ! qrad=0.0
1019      ! ! normally nb_mode=3 in physiq.def !!!
1020      ! do ilon=1,nbr_lon
1021      !    ! mode 1
1022      !    do i_lev=cloudmin,20
1023      !       r_median(ilon,i_lev,1)=0.125e-6
1024      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1025      !       stddev(ilon,i_lev,1)=1.57
1026      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1027      !       k_mass(ilon,i_lev,1)=1.0
1028      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1029      !    end do
1030
1031      !    r_median(ilon,21,1)=0.125e-6
1032      !    print*,'level',21,'r r_median',r_median(1,21,1)
1033      !    stddev(ilon,21,1)=1.57
1034      !    print*,'level',21,'dev std',stddev(1,21,1)
1035      !    k_mass(ilon,21,1)=1.0
1036      !    print*,'level',21,'coeff mass: k_mass',k_mass(1,21,1)
1037
1038      !    r_median(ilon,22,1)=0.2e-6
1039      !    print*,'level',22,'r r_median',r_median(1,22,1)
1040      !    stddev(ilon,22,1)=1.8
1041      !    print*,'level',22,'dev std',stddev(1,22,1)
1042      !    k_mass(ilon,22,1)=0.01
1043      !    print*,'level',22,'coeff mass: k_mass',k_mass(1,22,1)
1044
1045      !    r_median(ilon,23,1)=0.15e-6
1046      !    print*,'level',23,'r r_median',r_median(1,23,1)
1047      !    stddev(ilon,23,1)=1.8
1048      !    print*,'level',23,'dev std',stddev(1,23,1)
1049      !    k_mass(ilon,23,1)=0.01
1050      !    print*,'level',23,'coeff mass: k_mass',k_mass(1,23,1)
1051
1052      !    do i_lev=24,25
1053      !       r_median(ilon,i_lev,1)=0.15e-6
1054      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1055      !       stddev(ilon,i_lev,1)=1.9
1056      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1057      !       k_mass(ilon,i_lev,1)=0.01
1058      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1059      !    end do
1060
1061      !    r_median(ilon,26,1)=0.175e-6
1062      !    print*,'level',26,'r r_median',r_median(1,26,1)
1063      !    stddev(ilon,26,1)=2.16
1064      !    print*,'level',26,'dev std',stddev(1,26,1)
1065      !    k_mass(ilon,26,1)=0.175
1066      !    print*,'level',26,'coeff mass: k_mass',k_mass(1,26,1)
1067
1068      !    do i_lev=27,33
1069      !       r_median(ilon,i_lev,1)=0.175e-6
1070      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1071      !       stddev(ilon,i_lev,1)=2.16
1072      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1073      !       k_mass(ilon,i_lev,1)=0.25
1074      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1075      !    end do
1076
1077      !    do i_lev=34,cloudmax
1078      !       r_median(ilon,i_lev,1)=0.175e-6
1079      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1080      !       stddev(ilon,i_lev,1)=2.16
1081      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1082      !       k_mass(ilon,i_lev,1)=0.25
1083      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1084      !    end do
1085
1086      !    !  mode 2
1087      !    do i_lev=cloudmin,20
1088      !       r_median(ilon,i_lev,2)=1.4e-6
1089      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
1090      !       stddev(ilon,i_lev,2)=1.35
1091      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
1092      !       k_mass(ilon,i_lev,2)=0.0
1093      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
1094      !    end do
1095
1096      !    r_median(ilon,21,2)=1.4e-6
1097      !    print*,'level',21,'r r_median',r_median(1,21,2)
1098      !    stddev(ilon,21,2)=1.25
1099      !    print*,'level',21,'dev std',stddev(1,21,2)
1100      !    k_mass(ilon,21,2)=0.0
1101      !    print*,'level',21,'coeff mass: k_mass',k_mass(1,21,2)
1102
1103      !    r_median(ilon,22,2)=1.4e-6
1104      !    print*,'level',22,'r r_median',r_median(1,22,2)
1105      !    stddev(ilon,22,2)=1.25
1106      !    print*,'level',22,'dev std',stddev(1,22,2)
1107      !    k_mass(ilon,22,2)=0.04
1108      !    print*,'level',22,'coeff mass: k_mass',k_mass(1,22,2)
1109
1110      !    r_median(ilon,23,2)=1.35e-6
1111      !    print*,'level',23,'r r_median',r_median(1,23,2)
1112      !    stddev(ilon,23,2)=1.25
1113      !    print*,'level',23,'dev std',stddev(1,23,2)
1114      !    k_mass(ilon,23,2)=0.04
1115      !    print*,'level',23,'coeff mass: k_mass',k_mass(1,23,2)
1116
1117
1118      !    r_median(ilon,24,2)=1.375e-6
1119      !    print*,'level',24,'r r_median',r_median(1,24,2)
1120      !    stddev(ilon,24,2)=1.2
1121      !    print*,'level',24,'dev std',stddev(1,24,2)
1122      !    k_mass(ilon,24,2)=0.04
1123      !    print*,'level',24,'coeff mass: k_mass',k_mass(1,24,2)
1124
1125
1126      !    r_median(ilon,25,2)=1.4e-6
1127      !    print*,'level',25,'r r_median',r_median(1,25,2)
1128      !    stddev(ilon,25,2)=1.16
1129      !    print*,'level',25,'dev std',stddev(1,25,2)
1130      !    k_mass(ilon,25,2)=0.04
1131      !    print*,'level',25,'coeff mass: k_mass',k_mass(1,25,2)
1132
1133
1134      !    r_median(ilon,26,2)=1.15e-6
1135      !    print*,'level',26,'r r_median',r_median(1,26,2)
1136      !    stddev(ilon,26,2)=1.34
1137      !    print*,'level',26,'dev std',stddev(1,26,2)
1138      !    k_mass(ilon,26,2)=0.825
1139      !    print*,'level',26,'coeff mass: k_mass',k_mass(1,26,2)
1140
1141
1142      !    r_median(ilon,27,2)=1.14e-6
1143      !    print*,'level',27,'r r_median',r_median(1,27,2)
1144      !    stddev(ilon,27,2)=1.33
1145      !    print*,'level',27,'dev std',stddev(1,27,2)
1146      !    k_mass(ilon,27,2)=0.75
1147      !    print*,'level',27,'coeff mass: k_mass',k_mass(1,27,2)
1148
1149
1150      !    r_median(ilon,28,2)=1.35e-6
1151      !    print*,'level',28,'r r_median',r_median(1,28,2)
1152      !    stddev(ilon,28,2)=1.32
1153      !    print*,'level',28,'dev std',stddev(1,28,2)
1154      !    k_mass(ilon,28,2)=0.75
1155      !    print*,'level',28,'coeff mass: k_mass',k_mass(1,28,2)
1156
1157
1158      !    r_median(ilon,29,2)=1.125e-6
1159      !    print*,'level',29,'r r_median',r_median(1,29,2)
1160      !    stddev(ilon,29,2)=1.31
1161      !    print*,'level',29,'dev std',stddev(1,29,2)
1162      !    k_mass(ilon,29,2)=0.75
1163      !    print*,'level',29,'coeff mass: k_mass',k_mass(1,29,2)
1164
1165
1166      !    r_median(ilon,30,2)=1.118e-6
1167      !    print*,'level',30,'r r_median',r_median(1,30,2)
1168      !    stddev(ilon,30,2)=1.30
1169      !    print*,'level',30,'dev std',stddev(1,30,2)
1170      !    k_mass(ilon,30,2)=0.75
1171      !    print*,'level',30,'coeff mass: k_mass',k_mass(1,30,2)
1172
1173
1174      !    r_median(ilon,31,2)=1.11e-6
1175      !    print*,'level',31,'r r_median',r_median(1,31,2)
1176      !    stddev(ilon,31,2)=1.29
1177      !    print*,'level',31,'dev std',stddev(1,31,2)
1178      !    k_mass(ilon,31,2)=0.75
1179      !    print*,'level',31,'coeff mass: k_mass',k_mass(1,31,2)
1180
1181      !    do i_lev=32,33
1182      !       r_median(ilon,i_lev,2)=1.1e-6
1183      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
1184      !       stddev(ilon,i_lev,2)=1.28
1185      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
1186      !       k_mass(ilon,i_lev,2)=0.75
1187      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
1188      !    end do
1189
1190      !    ! if k_mass > 0 it means we have a bimodal upper haze.
1191      !    do i_lev=34,cloudmax
1192      !       r_median(ilon,i_lev,2)=1.1e-6
1193      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
1194      !       stddev(ilon,i_lev,2)=1.28
1195      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
1196      !       k_mass(ilon,i_lev,2)=0.75
1197      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
1198      !    end do
1199
1200      !    ! mode 3
1201      !    do i_lev=cloudmin,21
1202      !       r_median(ilon,i_lev,3)=3.65e-6
1203      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
1204      !       stddev(ilon,i_lev,3)=1.28
1205      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
1206      !       k_mass(ilon,i_lev,3)=0.0
1207      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
1208      !    end do
1209
1210      !    do i_lev=22,25
1211      !       r_median(ilon,i_lev,3)=3.65e-6
1212      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
1213      !       stddev(ilon,i_lev,3)=1.28
1214      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
1215      !       k_mass(ilon,i_lev,3)=0.95
1216      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
1217      !    end do
1218
1219      !    do i_lev=26,cloudmax
1220      !       r_median(ilon,i_lev,3)=3.65e-6
1221      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
1222      !       stddev(ilon,i_lev,3)=1.28
1223      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
1224      !       k_mass(ilon,i_lev,3)=0.0
1225      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
1226      !    end do
1227      ! end do
1228
1229      ! print*,'r_median size',size(r_median)
1230
1231      ! print*,'r_median(1,28,3)',r_median(1,28,3)
1232      ! print*,'r_median(1,28,3) ne devrait pas exister puisque les dim sont (285,78,2)'
1233      ! print*,'r_median(1,28,3) a la valeur de k_mass(1,28,1)', k_mass(1,28,1)
1234      ! print*,'r_median(1,28,1) just to see',r_median(1,28,1)
1235      ! print*,'r_median size after adding some shenanigans',size(r_median)
1236
1237      ! ! ==================================================================
1238      ! ! initialisation stupid with mode 3 100 um
1239      ! ! ==================================================================
1240      ! ! mode 3 0% solid, fully liquid
1241      ! qrad=0.0
1242      ! ! normally nb_mode=3 in physiq.def !!!
1243      ! do ilon=1,nbr_lon
1244      !    ! mode 1
1245      !    do i_lev=cloudmin,20
1246      !       r_median(ilon,i_lev,1)=0.125e-6
1247      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1248      !       stddev(ilon,i_lev,1)=1.57
1249      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1250      !       k_mass(ilon,i_lev,1)=1.0
1251      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1252      !    end do
1253
1254      !    r_median(ilon,21,1)=0.125e-6
1255      !    print*,'level',21,'r r_median',r_median(1,21,1)
1256      !    stddev(ilon,21,1)=1.57
1257      !    print*,'level',21,'dev std',stddev(1,21,1)
1258      !    k_mass(ilon,21,1)=1.0
1259      !    print*,'level',21,'coeff mass: k_mass',k_mass(1,21,1)
1260
1261      !    r_median(ilon,22,1)=0.2e-6
1262      !    print*,'level',22,'r r_median',r_median(1,22,1)
1263      !    stddev(ilon,22,1)=1.8
1264      !    print*,'level',22,'dev std',stddev(1,22,1)
1265      !    k_mass(ilon,22,1)=0.0
1266      !    print*,'level',22,'coeff mass: k_mass',k_mass(1,22,1)
1267
1268      !    r_median(ilon,23,1)=0.15e-6
1269      !    print*,'level',23,'r r_median',r_median(1,23,1)
1270      !    stddev(ilon,23,1)=1.8
1271      !    print*,'level',23,'dev std',stddev(1,23,1)
1272      !    k_mass(ilon,23,1)=0.0
1273      !    print*,'level',23,'coeff mass: k_mass',k_mass(1,23,1)
1274
1275      !    do i_lev=24,25
1276      !       r_median(ilon,i_lev,1)=0.15e-6
1277      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1278      !       stddev(ilon,i_lev,1)=1.9
1279      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1280      !       k_mass(ilon,i_lev,1)=0.0
1281      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1282      !    end do
1283
1284      !    r_median(ilon,26,1)=0.175e-6
1285      !    print*,'level',26,'r r_median',r_median(1,26,1)
1286      !    stddev(ilon,26,1)=2.16
1287      !    print*,'level',26,'dev std',stddev(1,26,1)
1288      !    k_mass(ilon,26,1)=0.0
1289      !    print*,'level',26,'coeff mass: k_mass',k_mass(1,26,1)
1290
1291      !    do i_lev=27,33
1292      !       r_median(ilon,i_lev,1)=0.175e-6
1293      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1294      !       stddev(ilon,i_lev,1)=2.16
1295      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1296      !       k_mass(ilon,i_lev,1)=0.0
1297      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1298      !    end do
1299
1300      !    do i_lev=34,cloudmax
1301      !       r_median(ilon,i_lev,1)=0.175e-6
1302      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,1)
1303      !       stddev(ilon,i_lev,1)=2.16
1304      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,1)
1305      !       k_mass(ilon,i_lev,1)=0.0
1306      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,1)
1307      !    end do
1308
1309      !    !  mode 2
1310      !    do i_lev=cloudmin,20
1311      !       r_median(ilon,i_lev,2)=1.4e-6
1312      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
1313      !       stddev(ilon,i_lev,2)=1.35
1314      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
1315      !       k_mass(ilon,i_lev,2)=0.0
1316      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
1317      !    end do
1318
1319      !    r_median(ilon,21,2)=1.4e-6
1320      !    print*,'level',21,'r r_median',r_median(1,21,2)
1321      !    stddev(ilon,21,2)=1.25
1322      !    print*,'level',21,'dev std',stddev(1,21,2)
1323      !    k_mass(ilon,21,2)=0.0
1324      !    print*,'level',21,'coeff mass: k_mass',k_mass(1,21,2)
1325
1326      !    r_median(ilon,22,2)=1.4e-6
1327      !    print*,'level',22,'r r_median',r_median(1,22,2)
1328      !    stddev(ilon,22,2)=1.25
1329      !    print*,'level',22,'dev std',stddev(1,22,2)
1330      !    k_mass(ilon,22,2)=0.0
1331      !    print*,'level',22,'coeff mass: k_mass',k_mass(1,22,2)
1332
1333      !    r_median(ilon,23,2)=1.35e-6
1334      !    print*,'level',23,'r r_median',r_median(1,23,2)
1335      !    stddev(ilon,23,2)=1.25
1336      !    print*,'level',23,'dev std',stddev(1,23,2)
1337      !    k_mass(ilon,23,2)=0.0
1338      !    print*,'level',23,'coeff mass: k_mass',k_mass(1,23,2)
1339
1340
1341      !    r_median(ilon,24,2)=1.375e-6
1342      !    print*,'level',24,'r r_median',r_median(1,24,2)
1343      !    stddev(ilon,24,2)=1.2
1344      !    print*,'level',24,'dev std',stddev(1,24,2)
1345      !    k_mass(ilon,24,2)=0.0
1346      !    print*,'level',24,'coeff mass: k_mass',k_mass(1,24,2)
1347
1348
1349      !    r_median(ilon,25,2)=1.4e-6
1350      !    print*,'level',25,'r r_median',r_median(1,25,2)
1351      !    stddev(ilon,25,2)=1.16
1352      !    print*,'level',25,'dev std',stddev(1,25,2)
1353      !    k_mass(ilon,25,2)=0.0
1354      !    print*,'level',25,'coeff mass: k_mass',k_mass(1,25,2)
1355
1356
1357      !    r_median(ilon,26,2)=1.15e-6
1358      !    print*,'level',26,'r r_median',r_median(1,26,2)
1359      !    stddev(ilon,26,2)=1.34
1360      !    print*,'level',26,'dev std',stddev(1,26,2)
1361      !    k_mass(ilon,26,2)=0.0
1362      !    print*,'level',26,'coeff mass: k_mass',k_mass(1,26,2)
1363
1364
1365      !    r_median(ilon,27,2)=1.14e-6
1366      !    print*,'level',27,'r r_median',r_median(1,27,2)
1367      !    stddev(ilon,27,2)=1.33
1368      !    print*,'level',27,'dev std',stddev(1,27,2)
1369      !    k_mass(ilon,27,2)=0.0
1370      !    print*,'level',27,'coeff mass: k_mass',k_mass(1,27,2)
1371
1372
1373      !    r_median(ilon,28,2)=1.35e-6
1374      !    print*,'level',28,'r r_median',r_median(1,28,2)
1375      !    stddev(ilon,28,2)=1.32
1376      !    print*,'level',28,'dev std',stddev(1,28,2)
1377      !    k_mass(ilon,28,2)=0.0
1378      !    print*,'level',28,'coeff mass: k_mass',k_mass(1,28,2)
1379
1380
1381      !    r_median(ilon,29,2)=1.125e-6
1382      !    print*,'level',29,'r r_median',r_median(1,29,2)
1383      !    stddev(ilon,29,2)=1.31
1384      !    print*,'level',29,'dev std',stddev(1,29,2)
1385      !    k_mass(ilon,29,2)=0.0
1386      !    print*,'level',29,'coeff mass: k_mass',k_mass(1,29,2)
1387
1388
1389      !    r_median(ilon,30,2)=1.118e-6
1390      !    print*,'level',30,'r r_median',r_median(1,30,2)
1391      !    stddev(ilon,30,2)=1.30
1392      !    print*,'level',30,'dev std',stddev(1,30,2)
1393      !    k_mass(ilon,30,2)=0.0
1394      !    print*,'level',30,'coeff mass: k_mass',k_mass(1,30,2)
1395
1396
1397      !    r_median(ilon,31,2)=1.11e-6
1398      !    print*,'level',31,'r r_median',r_median(1,31,2)
1399      !    stddev(ilon,31,2)=1.29
1400      !    print*,'level',31,'dev std',stddev(1,31,2)
1401      !    k_mass(ilon,31,2)=0.0
1402      !    print*,'level',31,'coeff mass: k_mass',k_mass(1,31,2)
1403
1404      !    do i_lev=32,33
1405      !       r_median(ilon,i_lev,2)=1.1e-6
1406      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
1407      !       stddev(ilon,i_lev,2)=1.28
1408      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
1409      !       k_mass(ilon,i_lev,2)=0.0
1410      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
1411      !    end do
1412
1413      !    ! if k_mass > 0 it means we have a bimodal upper haze.
1414      !    do i_lev=34,cloudmax
1415      !       r_median(ilon,i_lev,2)=1.1e-6
1416      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,2)
1417      !       stddev(ilon,i_lev,2)=1.28
1418      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,2)
1419      !       k_mass(ilon,i_lev,2)=0.0
1420      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,2)
1421      !    end do
1422
1423      !    ! mode 3
1424      !    do i_lev=cloudmin,21
1425      !       r_median(ilon,i_lev,3)=100.0e-6
1426      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
1427      !       stddev(ilon,i_lev,3)=1.28
1428      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
1429      !       k_mass(ilon,i_lev,3)=1.0
1430      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
1431      !    end do
1432
1433      !    do i_lev=22,25
1434      !       r_median(ilon,i_lev,3)=100.0e-6
1435      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
1436      !       stddev(ilon,i_lev,3)=1.28
1437      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
1438      !       k_mass(ilon,i_lev,3)=1.0
1439      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
1440      !    end do
1441
1442      !    do i_lev=26,cloudmax
1443      !       r_median(ilon,i_lev,3)=100.0e-6
1444      !       print*,'level',i_lev,'r r_median',r_median(1,i_lev,3)
1445      !       stddev(ilon,i_lev,3)=1.28
1446      !       print*,'level',i_lev,'dev std',stddev(1,i_lev,3)
1447      !       k_mass(ilon,i_lev,3)=1.0
1448      !       print*,'level',i_lev,'coeff mass: k_mass',k_mass(1,i_lev,3)
1449      !    end do
1450      ! end do
1451
1452      ! print*,'r_median size',size(r_median)
1453
1454      ! print*,'r_median(1,28,3)',r_median(1,28,3)
1455      ! print*,'r_median(1,28,3) ne devrait pas exister puisque les dim sont (285,78,2)'
1456      ! print*,'r_median(1,28,3) a la valeur de k_mass(1,28,1)', k_mass(1,28,1)
1457      ! print*,'r_median(1,28,1) just to see',r_median(1,28,1)
1458      ! print*,'r_median size after adding some shenanigans',size(r_median)
1459
1460! check if sum of kmass = 1 at each level
1461
1462print*, 'chemparam_mod: start checking k_mass=1'
1463
1464do i_lev = cloudmin,cloudmax
1465   if (nbr_mode == 3) then
1466      if ((k_mass(1,i_lev,1) + k_mass(1,i_lev,2) + k_mass(1,i_lev,3)) /= 1.0) then
1467         print*, 'kmass total is not 1.0'
1468         print*, k_mass(1,i_lev,1) + k_mass(1,i_lev,2) + k_mass(1,i_lev,3)
1469         print*, 'at level: ',i_lev
1470         print*, 'check chemparam_mod cloud structure definition'
1471         stop
1472      end if
1473   else
1474      if ((k_mass(1,i_lev,1) + k_mass(1,i_lev,2)) /= 1.0) then
1475         print*, 'kmass total is not 1.0'
1476         print*, k_mass(1,i_lev,1) + k_mass(1,i_lev,2)
1477         print*, 'at level: ',i_lev
1478         print*, 'check chemparam_mod cloud structure definition'
1479      end if
1480   end if
1481end do
1482
1483print*, 'chemparam_mod: end checking k_mass=1'
1484print*, 'k_mass is fine'
1485print*
1486print*,'==============================='
1487print*,'end initialisation cloud layer'
1488print*,'==============================='
1489
1490end subroutine cloud_ini
1491
1492!============================================================================
1493
1494  subroutine chemparam_ini
1495
1496!============================================================================
1497
1498      use infotrac_phy, only: nqtot, tname
1499
1500      implicit none
1501
1502      integer :: i
1503
1504      allocate(m_tr(nqtot))    ! molecular mass of tracers
1505      allocate(type_tr(nqtot)) ! type of chemical tracers 1: neutral, 2: ion, 3: liquid
1506
1507      ! index for chemical tracers
1508
1509      ! neutrals
1510
1511      i_co2      = 0
1512      i_co       = 0
1513      i_h2       = 0
1514      i_h2o      = 0
1515      i_o1d      = 0
1516      i_o        = 0
1517      i_o2       = 0
1518      i_o2dg     = 0
1519      i_o3       = 0
1520      i_h        = 0
1521      i_oh       = 0
1522      i_ho2      = 0
1523      i_h2o2     = 0
1524      i_cl       = 0
1525      i_clo      = 0
1526      i_cl2      = 0
1527      i_hcl      = 0
1528      i_hocl     = 0
1529      i_clco     = 0
1530      i_clco3    = 0
1531      i_cocl2    = 0
1532      i_s        = 0
1533      i_so       = 0
1534      i_so2      = 0
1535      i_so3      = 0
1536      i_osso_cis = 0
1537      i_osso_trans = 0
1538      i_s2o2_cyc = 0
1539      i_ocs      = 0
1540      i_hso3     = 0
1541      i_h2so4    = 0
1542      i_s2       = 0
1543      i_clso2    = 0
1544      i_oscl     = 0
1545      i_n2       = 0
1546      i_he       = 0
1547      i_n2d      = 0
1548      i_n        = 0
1549      i_no       = 0
1550      i_no2      = 0
1551
1552      ! ions
1553
1554      i_co2plus  = 0
1555      i_coplus   = 0
1556      i_oplus    = 0
1557      i_o2plus   = 0
1558      i_n2plus   = 0
1559      i_hplus    = 0
1560      i_h2oplus  = 0
1561      i_nplus    = 0
1562      i_ohplus   = 0
1563      i_cplus    = 0
1564      i_noplus   = 0
1565      i_h3oplus  = 0
1566      i_hcoplus  = 0
1567      i_hco2plus = 0
1568      i_elec     = 0
1569
1570      ! liquid
1571
1572      i_h2oliq   = 0
1573      i_h2so4liq = 0
1574
1575      do i = 1,nqtot
1576         print*,'tname(i)',tname(i)
1577         select case(tname(i))
1578
1579            ! neutrals
1580
1581            case('co2')
1582               i_co2 = i
1583               print*,'co2',i_co2
1584               m_tr(i_co2) = 44.0095
1585               type_tr(i_co2) = 1
1586            case('co')
1587               i_co = i
1588               print*,'co',i_co
1589               m_tr(i_co) = 28.0101
1590               type_tr(i_co) = 1
1591            case('h2')
1592               i_h2 = i
1593               print*,'h2',i_h2
1594               m_tr(i_h2) = 2.01588
1595               type_tr(i_h2) = 1
1596            case('h2o')
1597               i_h2o = i
1598               print*,'h2o',i_h2o
1599               m_tr(i_h2o) = 18.0153
1600               type_tr(i_h2o) = 1
1601            case('o1d')
1602               i_o1d = i
1603               print*,'o1d',i_o1d
1604               m_tr(i_o1d) = 15.994
1605               type_tr(i_o1d) = 1
1606            case('o')
1607               i_o = i
1608               print*,'o',i_o
1609               m_tr(i_o) = 15.994
1610               type_tr(i_o) = 1
1611            case('o2')
1612               i_o2 = i
1613               print*,'o2',i_o2
1614               m_tr(i_o2) = 31.9988
1615               type_tr(i_o2) = 1
1616            case('o2dg')
1617               i_o2dg = i
1618               print*,'o2dg',i_o2dg
1619               m_tr(i_o2dg) = 31.9988
1620               type_tr(i_o2dg) = 1
1621            case('o3')
1622               i_o3 = i
1623               print*,'o3',i_o3
1624               m_tr(i_o3) = 47.9982
1625               type_tr(i_o3) = 1
1626            case('h')
1627               i_h = i
1628               print*,'h',i_h
1629               m_tr(i_h) = 1.00794
1630               type_tr(i_h) = 1
1631            case('oh')
1632               i_oh = i
1633               print*,'oh',i_oh
1634               m_tr(i_oh) = 17.0073
1635               type_tr(i_oh) = 1
1636            case('ho2')
1637               i_ho2 = i
1638               print*,'ho2',i_ho2
1639               m_tr(i_ho2) = 33.0067
1640               type_tr(i_ho2) = 1
1641            case('h2o2')
1642               i_h2o2 = i
1643               print*,'h2o2',i_h2o2
1644               m_tr(i_h2o2) = 34.0147
1645               type_tr(i_h2o2) = 1
1646            case('cl')
1647               i_cl = i
1648               print*,'cl',i_cl
1649               m_tr(i_cl) = 35.453
1650               type_tr(i_cl) = 1
1651            case('clo')
1652               i_clo = i
1653               print*,'clo',i_clo
1654               m_tr(i_clo) = 51.452
1655               type_tr(i_clo) = 1
1656            case('cl2')
1657               i_cl2 = i
1658               print*,'cl2',i_cl2
1659               m_tr(i_cl2) = 70.906
1660               type_tr(i_cl2) = 1
1661            case('hcl')
1662               i_hcl = i
1663               print*,'hcl',i_hcl
1664               m_tr(i_hcl) = 36.461
1665               type_tr(i_hcl) = 1
1666            case('hocl')
1667               i_hocl = i
1668               print*,'hocl',i_hocl
1669               m_tr(i_hocl) = 52.46
1670               type_tr(i_hocl) = 1
1671            case('clco')
1672               i_clco = i
1673               print*,'clco',i_clco
1674               m_tr(i_clco) = 63.463
1675               type_tr(i_clco) = 1
1676            case('clco3')
1677               i_clco3 = i
1678               print*,'clco3',i_clco3
1679               m_tr(i_clco3) = 95.462
1680               type_tr(i_clco3) = 1
1681            case('cocl2')
1682               i_cocl2 = i
1683               print*,'cocl2',i_cocl2
1684               m_tr(i_cocl2) = 98.916
1685               type_tr(i_cocl2) = 1
1686            case('s')
1687               i_s = i
1688               print*,'s',i_s
1689               m_tr(i_s) = 32.065
1690               type_tr(i_s) = 1
1691            case('so')
1692               i_so = i
1693               print*,'so',i_so
1694               m_tr(i_so) = 48.0644
1695               type_tr(i_so) = 1
1696            case('so2')
1697               i_so2 = i
1698               print*,'so2',i_so2
1699               m_tr(i_so2) = 64.064
1700               type_tr(i_so2) = 1
1701            case('so3')
1702               i_so3 = i
1703               print*,'so3',i_so3
1704               m_tr(i_so3) = 80.063
1705               type_tr(i_so3) = 1
1706            case('osso_cis')
1707               i_osso_cis = i
1708               print*,'osso_cis',i_osso_cis
1709               m_tr(i_osso_cis)= 96.1288
1710               type_tr(i_osso_cis) = 1
1711            case('osso_trans')
1712               i_osso_trans = i
1713               print*,'osso_trans',i_osso_trans
1714               m_tr(i_osso_trans)= 96.1288
1715               type_tr(i_osso_trans) = 1
1716            case('s2o2_cyc')
1717               i_s2o2_cyc = i
1718               print*,'s2o2_cyc',i_s2o2_cyc
1719               m_tr(i_s2o2_cyc)= 96.1288
1720               type_tr(i_s2o2_cyc) = 1
1721            case('ocs')
1722               i_ocs = i
1723               print*,'ocs',i_ocs
1724               m_tr(i_ocs) = 60.0751
1725               type_tr(i_ocs) = 1
1726            case('hso3')
1727               i_hso3 = i
1728               print*,'hso3',i_hso3
1729               m_tr(i_hso3) = 81.071
1730               type_tr(i_hso3) = 1
1731            case('h2so4')
1732               i_h2so4 = i
1733               print*,'h2so4',i_h2so4
1734               m_tr(i_h2so4) = 98.078
1735               type_tr(i_h2so4) = 1
1736            case('s2')
1737               i_s2 = i
1738               print*,'s2',i_s2
1739               m_tr(i_s2) = 64.13
1740               type_tr(i_s2) = 1
1741            case('clso2')
1742               i_clso2 = i
1743               print*,'clso2',i_clso2
1744               m_tr(i_clso2) = 99.517
1745               type_tr(i_clso2) = 1
1746            case('oscl')
1747               i_oscl = i
1748               print*,'oscl',i_oscl
1749               m_tr(i_oscl) = 83.517
1750               type_tr(i_oscl) = 1
1751            case('n2')
1752               i_n2 = i
1753               print*,'n2',i_n2
1754               m_tr(i_n2) = 28.013
1755               type_tr(i_n2) = 1
1756            case('he')
1757               i_he = i
1758               print*,'he',i_he
1759               m_tr(i_he) = 4.0026
1760               type_tr(i_he) = 1
1761            case('n2d')
1762               i_n2d = i
1763               print*,'n2d',i_n2d
1764               m_tr(i_n2d) = 14.0067
1765               type_tr(i_n2d) = 1
1766            case('n')
1767               i_n = i
1768               print*,'n',i_n
1769               m_tr(i_n) = 14.0067
1770               type_tr(i_n) = 1
1771            case('no')
1772               i_no = i
1773               print*,'no',i_no
1774               m_tr(i_no) = 30.0061
1775               type_tr(i_no) = 1
1776            case('no2')
1777               i_no2 = i
1778               print*,'no2',i_no2
1779               m_tr(i_no2) = 46.0055
1780               type_tr(i_no2) = 1
1781
1782            ! ions
1783
1784            case('co2plus')
1785               i_co2plus = i
1786               print*,'co2plus',i_co2plus
1787               m_tr(i_co2plus) = 44.0095
1788               type_tr(i_co2plus) = 2
1789            case('coplus')
1790               i_coplus = i
1791               print*,'coplus',i_coplus
1792               m_tr(i_coplus) = 28.0101
1793               type_tr(i_coplus) = 2
1794            case('oplus')
1795               i_oplus = i
1796               print*,'oplus',i_oplus
1797               m_tr(i_oplus) = 15.994
1798               type_tr(i_oplus) = 2
1799            case('o2plus')
1800               i_o2plus = i
1801               print*,'o2plus',i_o2plus
1802               m_tr(i_o2plus) = 31.9988
1803               type_tr(i_o2plus) = 2
1804            case('n2plus')
1805               i_n2plus = i
1806               print*,'n2plus',i_n2plus
1807               m_tr(i_n2plus) = 28.013
1808               type_tr(i_n2plus) = 2
1809            case('hplus')
1810               i_hplus = i
1811               print*,'hplus',i_hplus
1812               m_tr(i_hplus) = 1.00794
1813               type_tr(i_hplus) = 2
1814            case('h2oplus')
1815               i_h2oplus = i
1816               print*,'h2oplus',i_h2oplus
1817               m_tr(i_h2oplus) = 18.0153
1818               type_tr(i_h2oplus) = 2
1819            case('nplus')
1820               i_nplus = i
1821               print*,'nplus',i_nplus
1822               m_tr(i_nplus) = 14.0067
1823               type_tr(i_nplus) = 2
1824            case('ohplus')
1825               i_ohplus = i
1826               print*,'ohplus',i_ohplus
1827               m_tr(i_ohplus) = 17.0073
1828               type_tr(i_ohplus) = 2
1829            case('cplus')
1830               i_cplus = i
1831               print*,'cplus',i_cplus
1832               m_tr(i_cplus) = 12.011
1833               type_tr(i_cplus) = 2
1834            case('noplus')
1835               i_noplus = i
1836               print*,'noplus',i_noplus
1837               m_tr(i_noplus) = 30.0061
1838               type_tr(i_noplus) = 2
1839            case('h3oplus')
1840               i_h3oplus = i
1841               print*,'h3oplus',i_h3oplus
1842               m_tr(i_h3oplus) = 19.0232
1843               type_tr(i_h3oplus) = 2
1844            case('hcoplus')
1845               i_hcoplus = i
1846               print*,'hcoplus',i_hcoplus
1847               m_tr(i_hcoplus) = 29.0180
1848               type_tr(i_hcoplus) = 2
1849            case('hco2plus')
1850               i_hco2plus = i
1851               print*,'hco2plus',i_hco2plus
1852               m_tr(i_hco2plus) = 45.
1853               type_tr(i_hco2plus) = 2
1854            case('elec')
1855               i_elec = i
1856               print*,'elec',i_elec
1857               m_tr(i_elec) = 1./1822.89
1858               type_tr(i_elec) = 2
1859
1860            ! liquid tracers (cl_scheme = 1)
1861
1862            case('h2oliq')
1863               i_h2oliq = i
1864               print*,'h2oliq',i_h2oliq
1865               m_tr(i_h2oliq) = 18.0153
1866               type_tr(i_h2oliq) = 3
1867            case('h2so4liq')
1868               i_h2so4liq = i
1869               print*,'h2so4liq',i_h2so4liq
1870               m_tr(i_h2so4liq) = 98.078
1871               type_tr(i_h2so4liq) = 3
1872
1873            ! liquid tracers (cl_scheme = 2)
1874
1875            case('m0_aer')
1876               i_m0_aer = i
1877               print*,'m0_aer',i_m0_aer
1878               type_tr(i_m0_aer) = 10
1879            case('m3_aer')
1880               i_m3_aer = i
1881               print*,'m3_aer',i_m3_aer
1882               type_tr(i_m3_aer) = 10
1883            case('m0_m1drop')
1884               i_m0_mode1drop = i
1885               print*,'m0_m1drop',i_m0_mode1drop
1886               type_tr(i_m0_mode1drop) = 10
1887            case('m0_m1ccn')
1888               i_m0_mode1ccn = i
1889               print*,'m0_m1ccn',i_m0_mode1ccn
1890               type_tr(i_m0_mode1ccn) = 10
1891            case('m3_m1sa')
1892               i_m3_mode1sa = i
1893               print*,'m3_m1sa',i_m3_mode1sa
1894               type_tr(i_m3_mode1sa) = 10
1895            case('m3_m1w')
1896               i_m3_mode1w = i
1897               print*,'m3_m1w',i_m3_mode1w
1898               type_tr(i_m3_mode1w) = 10
1899            case('m3_m1ccn')
1900               i_m3_mode1ccn = i
1901               print*,'m3_m1ccn',i_m3_mode1ccn
1902               type_tr(i_m3_mode1ccn) = 10
1903            case('m0_m2drop')
1904               i_m0_mode2drop = i
1905               print*,'m0_m2drop',i_m0_mode2drop
1906               type_tr(i_m0_mode2drop) = 10
1907            case('m0_m2ccn')
1908               i_m0_mode2ccn = i
1909               print*,'m0_m2ccn',i_m0_mode2ccn
1910               type_tr(i_m0_mode2ccn) = 10
1911            case('m3_m2sa')
1912               i_m3_mode2sa = i
1913               print*,'m3_m2sa',i_m3_mode2sa
1914               type_tr(i_m3_mode2sa) = 10
1915            case('m3_m2w')
1916               i_m3_mode2w = i
1917               print*,'m3_m2w',i_m3_mode2w
1918               type_tr(i_m3_mode2w) = 10
1919            case('m3_m2ccn')
1920               i_m3_mode2ccn = i
1921               print*,'m3_m2ccn',i_m3_mode2ccn
1922               type_tr(i_m3_mode2ccn) = 10
1923         end select
1924      end do
1925
1926  end subroutine chemparam_ini
1927
1928!============================================================================
1929
1930  subroutine vapors4muphy_ini(nlon,nlev,trac)
1931
1932!============================================================================
1933
1934  use infotrac_phy, only: nqtot, tname
1935
1936  integer :: nlon, nlev
1937  real    :: trac(nlon,nlev,nqtot) ! traceur ( en vmr)
1938
1939!  integer :: i
1940!  real    :: trac1d(nlev,2) ! traceur lu ( en vmr)
1941
1942! lecture d'un fichier texte contenant les profils de trac1d(:1) = h2o et trac1d(:,2) = h2so4
1943!  do i=1,nlon
1944!     trac(i,:,i_h2o) = trac1d(:,1)
1945!     trac(i,:,i_h2so4) = trac1d(:,2)
1946!  enddo
1947
1948!  intitialisation profils altitude h2o et h2so4
1949!  profil h2o initial vap+liq == que vap
1950   trac(:,1:24,i_h2o) = 30.e-6 !
1951   trac(:,25:50,i_h2o) = 1.e-6 !
1952
1953   trac(:,:,i_h2so4) = 3.e-9 ! limite sup sandor 2012
1954   trac(:,23:50,i_h2so4) = 2.e-6 ! profil h2so4 initial => vap+liq
1955
1956  end subroutine vapors4muphy_ini
1957
1958end module chemparam_mod
Note: See TracBrowser for help on using the repository browser.