1 | !*********************************************************************** |
---|
2 | |
---|
3 | module chimiedata_h |
---|
4 | |
---|
5 | !*********************************************************************** |
---|
6 | ! |
---|
7 | ! subject: |
---|
8 | ! -------- |
---|
9 | ! |
---|
10 | ! data for photochemistry |
---|
11 | ! |
---|
12 | ! update 06/03/2021 set to module - add global variable (Yassin Jaziri) |
---|
13 | ! |
---|
14 | ! |
---|
15 | !*********************************************************************** |
---|
16 | |
---|
17 | implicit none |
---|
18 | |
---|
19 | character*30, save, allocatable :: chemnoms(:) ! name of chemical tracers |
---|
20 | |
---|
21 | !-------------------------------------------- |
---|
22 | ! dimensions of photolysis lookup table |
---|
23 | !-------------------------------------------- |
---|
24 | |
---|
25 | integer, save :: nd ! species |
---|
26 | integer, save :: nz ! altitude |
---|
27 | integer, save :: nozo ! ozone |
---|
28 | integer, save :: nsza ! solar zenith angle |
---|
29 | integer, save :: ntemp ! temperature |
---|
30 | integer, save :: ntau ! dust |
---|
31 | integer, save :: nch4 ! ch4 |
---|
32 | |
---|
33 | !-------------------------------------------- |
---|
34 | |
---|
35 | real, parameter :: kb = 1.3806e-23 |
---|
36 | |
---|
37 | ! photolysis |
---|
38 | real, save, allocatable :: jphot(:,:,:,:,:,:,:) |
---|
39 | real, save, allocatable :: ephot(:,:,:,:,:,:,:) |
---|
40 | real, save, allocatable :: colairtab(:) |
---|
41 | real, save, allocatable :: szatab(:) |
---|
42 | real, save, allocatable :: table_ozo(:) |
---|
43 | real, save, allocatable :: tautab(:) |
---|
44 | real, save, allocatable :: table_ch4(:) |
---|
45 | |
---|
46 | ! deposition |
---|
47 | real, save, allocatable :: SF_value(:) ! SF_mode=1 (vmr) SF_mode=2 (cm.s-1) |
---|
48 | real, save, allocatable :: prod_rate(:) ! production flux in molecules/m2/s |
---|
49 | real, save, allocatable :: surface_flux(:,:) ! Surface flux from deposition molecules/m2/s |
---|
50 | real, save, allocatable :: surface_flux2(:,:) ! Surface flux mean molecules/m2/s |
---|
51 | real, save, allocatable :: escape(:,:) ! escape rate in molecules/m2/s |
---|
52 | integer, save, allocatable :: SF_mode(:) ! SF_mode=1 (fixed mixing ratio) / 2 (fixed sedimentation velocity in cm/s and/or flux in molecules/m^3) |
---|
53 | |
---|
54 | !-------------------------------------------- |
---|
55 | |
---|
56 | real, save, allocatable :: albedo_snow_chim(:) ! albedo snow on chemistry wavelenght grid |
---|
57 | real, save, allocatable :: albedo_co2_ice_chim(:) ! albedo co2 ice on chemistry wavelenght grid |
---|
58 | |
---|
59 | !-------------------------------------------- |
---|
60 | ! For hard coding reactions |
---|
61 | !-------------------------------------------- |
---|
62 | |
---|
63 | integer, parameter :: nphot_hard_coding = 0 |
---|
64 | integer, parameter :: n4_hard_coding = 0 |
---|
65 | integer, parameter :: n3_hard_coding = 0 |
---|
66 | |
---|
67 | contains |
---|
68 | |
---|
69 | subroutine indice_HC(nb_phot,nb_reaction_4,nb_reaction_3) |
---|
70 | |
---|
71 | use types_asis |
---|
72 | |
---|
73 | implicit none |
---|
74 | |
---|
75 | integer, intent(inout) :: nb_phot |
---|
76 | integer, intent(inout) :: nb_reaction_4 |
---|
77 | integer, intent(inout) :: nb_reaction_3 |
---|
78 | |
---|
79 | !!!!!!!!!!! Hard coding reaction !!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
80 | !=========================================================== |
---|
81 | ! d022 : HO2 + NO + M -> HNO3 + M |
---|
82 | !=========================================================== |
---|
83 | !nb_reaction_4 = nb_reaction_4 + 1 |
---|
84 | !indice_4(nb_reaction_4) = z4spec(1.0, indexchim('ho2'), 1.0, indexchim('no'), 1.0, indexchim('hno3'), 0.0, 1) |
---|
85 | |
---|
86 | !=========================================================== |
---|
87 | ! e001 : CO + OH -> CO2 + H |
---|
88 | !=========================================================== |
---|
89 | !nb_reaction_4 = nb_reaction_4 + 1 |
---|
90 | !indice_4(nb_reaction_4) = z4spec(1.0, indexchim('co'), 1.0, indexchim('oh'), 1.0, indexchim('co2'), 1.0, indexchim('h')) |
---|
91 | |
---|
92 | !=========================================================== |
---|
93 | ! f028 : CH + CH4 -> C2H4 + H |
---|
94 | !=========================================================== |
---|
95 | !nb_reaction_4 = nb_reaction_4 + 1 |
---|
96 | !indice_4(nb_reaction_4) = z4spec(1.0, indexchim('ch'), 1.0, indexchim('ch4'), 1.0, indexchim('c2h4'), 1.0, indexchim('h')) |
---|
97 | |
---|
98 | !=========================================================== |
---|
99 | ! r001 : HNO3 + rain -> H2O |
---|
100 | !=========================================================== |
---|
101 | !nb_phot = nb_phot + 1 |
---|
102 | !indice_phot(nb_phot) = z3spec(1.0, indexchim('hno3'), 1.0, indexchim('h2o_vap'), 0.0, 1) |
---|
103 | |
---|
104 | !=========================================================== |
---|
105 | ! e001 : CO + OH -> CO2 + H |
---|
106 | !=========================================================== |
---|
107 | !nb_reaction_4 = nb_reaction_4 + 1 |
---|
108 | !indice_4(nb_reaction_4) = z4spec(1.0, indexchim('co'), 1.0, indexchim('oh'), 1.0, indexchim('co2'), 1.0, indexchim('h')) |
---|
109 | |
---|
110 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
111 | ! photodissociation of NO : NO + hv -> N + O |
---|
112 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
113 | !nb_phot = nb_phot + 1 |
---|
114 | !indice_phot(nb_phot) = z3spec(1.0, indexchim('no'), 1.0, indexchim('n'), 1.0, indexchim('o')) |
---|
115 | |
---|
116 | !!!!!!!!!!! END Hard coding reaction !!!!!!!!!!!!!!!!!!!!!!! |
---|
117 | end subroutine indice_HC |
---|
118 | |
---|
119 | subroutine reactionrates_HC(nlayer,nesp,dens,t,press,zmmean,sza,c,v_phot,v_4,v_3,nb_phot,nb_reaction_4,nb_reaction_3) |
---|
120 | |
---|
121 | use comcstfi_mod, only: g, avocado |
---|
122 | use types_asis |
---|
123 | |
---|
124 | implicit none |
---|
125 | |
---|
126 | integer, intent(in) :: nlayer ! number of atmospheric layers |
---|
127 | integer, intent(in) :: nesp ! number of chemical species |
---|
128 | integer, intent(inout) :: nb_phot |
---|
129 | integer, intent(inout) :: nb_reaction_4 |
---|
130 | integer, intent(inout) :: nb_reaction_3 |
---|
131 | real, intent(in), dimension(nlayer) :: dens ! total number density (molecule.cm-3) |
---|
132 | real, intent(in), dimension(nlayer) :: press ! pressure (hPa) |
---|
133 | real, intent(in), dimension(nlayer) :: t ! temperature (K) |
---|
134 | real, intent(in), dimension(nlayer) :: zmmean ! mean molar mass (g/mole) |
---|
135 | real, intent(in) :: sza ! solar zenith angle (deg) |
---|
136 | real (kind = 8), intent(in), dimension(nlayer,nesp) :: c ! species number density (molecule.cm-3) |
---|
137 | real (kind = 8), intent(inout), dimension(nlayer, nb_phot_max) :: v_phot |
---|
138 | real (kind = 8), intent(inout), dimension(nlayer,nb_reaction_4_max) :: v_4 |
---|
139 | real (kind = 8), intent(inout), dimension(nlayer,nb_reaction_3_max) :: v_3 |
---|
140 | |
---|
141 | ! local |
---|
142 | |
---|
143 | real, dimension(nlayer) :: colo3 ! ozone columns (molecule.cm-2) |
---|
144 | integer :: ilev |
---|
145 | real :: k1a0, k1b0, k1ainf, k1a, k1b, fc, fx, x, y |
---|
146 | real :: ak0, ak1, rate, xpo, deq, rate1, rate2, xpo1, xpo2 |
---|
147 | real :: rain_h2o, rain_rate |
---|
148 | |
---|
149 | !!!!!!!!!!! Hard coding reaction !!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
150 | !---------------------------------------------------------------------- |
---|
151 | ! nitrogen reactions |
---|
152 | !---------------------------------------------------------------------- |
---|
153 | |
---|
154 | !--- d022: ho2 + no + m -> hno3 + m |
---|
155 | ! Butkovskaya et al. 2007 |
---|
156 | |
---|
157 | !d022(:) = 6.4e-14*exp(1644/T(:)) |
---|
158 | !nb_reaction_4 = nb_reaction_4 + 1 |
---|
159 | !v_4(:,nb_reaction_4) = 3.3e-12*exp(270./t(:))*((530/T(:)+6.4*1e-4*press(:)/1.33322-1.73))*1e-2 |
---|
160 | |
---|
161 | !---------------------------------------------------------------------- |
---|
162 | ! carbon reactions |
---|
163 | !---------------------------------------------------------------------- |
---|
164 | |
---|
165 | !--- e001: oh + co -> co2 + h |
---|
166 | |
---|
167 | !nb_reaction_4 = nb_reaction_4 + 1 |
---|
168 | |
---|
169 | ! jpl 2003 |
---|
170 | |
---|
171 | !e001(:) = 1.5e-13*(1 + 0.6*press(:)/1013.) |
---|
172 | |
---|
173 | ! mccabe et al., grl, 28, 3135, 2001 |
---|
174 | |
---|
175 | !e001(:) = 1.57e-13 + 3.54e-33*dens(:) |
---|
176 | |
---|
177 | ! jpl 2006 |
---|
178 | |
---|
179 | !ak0 = 1.5e-13*(t(:)/300.)**(0.6) |
---|
180 | !ak1 = 2.1e-9*(t(:)/300.)**(6.1) |
---|
181 | !rate1 = ak0/(1. + ak0/(ak1/dens(:))) |
---|
182 | !xpo1 = 1./(1. + alog10(ak0/(ak1/dens(:)))**2) |
---|
183 | |
---|
184 | !ak0 = 5.9e-33*(t(:)/300.)**(-1.4) |
---|
185 | !ak1 = 1.1e-12*(t(:)/300.)**(1.3) |
---|
186 | !rate2 = (ak0*dens(:))/(1. + ak0*dens(:)/ak1) |
---|
187 | !xpo2 = 1./(1. + alog10((ak0*dens(:))/ak1)**2) |
---|
188 | |
---|
189 | !e001(:) = rate1*0.6**xpo1 + rate2*0.6**xpo2 |
---|
190 | |
---|
191 | ! jpl 2015 |
---|
192 | |
---|
193 | !do ilev = 1,nlayer |
---|
194 | !!oh + co -> h + co2 |
---|
195 | ! rate1 = 1.5e-13*(t(ilev)/300.)**(0.0) |
---|
196 | !!oh + co + m -> hoco |
---|
197 | ! ak0 = 5.9e-33*(t(ilev)/300.)**(-1.0) |
---|
198 | ! ak1 = 1.1e-12*(t(ilev)/300.)**(1.3) |
---|
199 | ! rate2 = (ak0*dens(ilev))/(1. + ak0*dens(ilev)/ak1) |
---|
200 | ! xpo2 = 1./(1. + alog10((ak0*dens(ilev))/ak1)**2) |
---|
201 | ! |
---|
202 | ! v_4(ilev,nb_reaction_4) = rate1 + rate2*0.6**xpo2 |
---|
203 | !end do |
---|
204 | |
---|
205 | ! joshi et al., 2006 |
---|
206 | |
---|
207 | !do ilev = 1,nlayer |
---|
208 | ! k1a0 = 1.34*2.5*dens(ilev) & |
---|
209 | ! *1/(1/(3.62e-26*t(ilev)**(-2.739)*exp(-20./t(ilev))) & |
---|
210 | ! + 1/(6.48e-33*t(ilev)**(0.14)*exp(-57./t(ilev)))) ! typo in paper corrected |
---|
211 | ! k1b0 = 1.17e-19*t(ilev)**(2.053)*exp(139./t(ilev)) & |
---|
212 | ! + 9.56e-12*t(ilev)**(-0.664)*exp(-167./t(ilev)) |
---|
213 | ! k1ainf = 1.52e-17*t(ilev)**(1.858)*exp(28.8/t(ilev)) & |
---|
214 | ! + 4.78e-8*t(ilev)**(-1.851)*exp(-318./t(ilev)) |
---|
215 | ! x = k1a0/(k1ainf - k1b0) |
---|
216 | ! y = k1b0/(k1ainf - k1b0) |
---|
217 | ! fc = 0.628*exp(-1223./t(ilev)) + (1. - 0.628)*exp(-39./t(ilev)) & |
---|
218 | ! + exp(-t(ilev)/255.) |
---|
219 | ! fx = fc**(1./(1. + (alog(x))**2)) ! typo in paper corrected |
---|
220 | ! k1a = k1a0*((1. + y)/(1. + x))*fx |
---|
221 | ! k1b = k1b0*(1./(1.+x))*fx |
---|
222 | ! |
---|
223 | ! v_4(ilev,nb_reaction_4) = k1a + k1b |
---|
224 | !end do |
---|
225 | |
---|
226 | !--- f028: ch + ch4 + m -> c2h4 + h + m |
---|
227 | ! Romani 1993 |
---|
228 | |
---|
229 | !nb_reaction_4 = nb_reaction_4 + 1 |
---|
230 | !v_4(:,nb_reaction_4) = min(2.5e-11*exp(200./t(:)),1.7e-10) |
---|
231 | |
---|
232 | !---------------------------------------------------------------------- |
---|
233 | ! washout r001 : HNO3 + rain -> H2O |
---|
234 | !---------------------------------------------------------------------- |
---|
235 | |
---|
236 | !nb_phot = nb_phot + 1 |
---|
237 | ! |
---|
238 | !rain_h2o = 100.e-6 |
---|
239 | !!rain_rate = 1.e-6 ! 10 days |
---|
240 | !rain_rate = 1.e-8 |
---|
241 | ! |
---|
242 | !do ilev = 1,nlayer |
---|
243 | ! if (c(ilev,indexchim('h2o_vap'))/dens(ilev) >= rain_h2o) then |
---|
244 | ! v_phot(ilev,nb_phot) = rain_rate |
---|
245 | ! else |
---|
246 | ! v_phot(ilev,nb_phot) = 0. |
---|
247 | ! end if |
---|
248 | !end do |
---|
249 | |
---|
250 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
251 | ! photodissociation of NO |
---|
252 | !ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc |
---|
253 | |
---|
254 | !nb_phot = nb_phot + 1 |
---|
255 | ! |
---|
256 | !colo3(nlayer) = 0. |
---|
257 | !! ozone columns for other levels (molecule.cm-2) |
---|
258 | !do ilev = nlayer-1,1,-1 |
---|
259 | ! colo3(ilev) = colo3(ilev+1) + (c(ilev+1,indexchim('o3')) + c(ilev,indexchim('o3')))*0.5*avocado*1e-4*((press(ilev) - press(ilev+1))*100.)/(1.e-3*zmmean(ilev)*g*dens(ilev)) |
---|
260 | !end do |
---|
261 | !call jno(nlayer, c(nlayer:1:-1,indexchim('no')), c(nlayer:1:-1,indexchim('o2')), colo3(nlayer:1:-1), dens(nlayer:1:-1), press(nlayer:1:-1), sza, v_phot(nlayer:1:-1,nb_phot)) |
---|
262 | |
---|
263 | !!!!!!!!!!! END Hard coding reaction !!!!!!!!!!!!!!!!!!!!!!! |
---|
264 | end subroutine reactionrates_HC |
---|
265 | |
---|
266 | subroutine jno(nivbas, c_no, c_o2, o3t, hnm, pm, sza, tjno) |
---|
267 | |
---|
268 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
269 | ! ! |
---|
270 | ! parametrisation de la photodissociation de no ! |
---|
271 | ! d'apres minschwaner and siskind, a new calculation ! |
---|
272 | ! of nitric oxide photolysis in the stratosphere, ! |
---|
273 | ! mesosphere, and lower thermosphere, j. geophys. res., ! |
---|
274 | ! 98, 20401-20412, 1993 ! |
---|
275 | ! ! |
---|
276 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
277 | |
---|
278 | implicit none |
---|
279 | |
---|
280 | ! input |
---|
281 | |
---|
282 | integer, intent(in) :: nivbas |
---|
283 | real (kind = 8), dimension(nivbas), intent(in) :: c_no, c_o2 |
---|
284 | real, dimension(nivbas), intent(in) :: hnm, o3t, pm |
---|
285 | real, intent(in) :: sza |
---|
286 | |
---|
287 | ! output |
---|
288 | |
---|
289 | real, dimension(nivbas), intent(out) :: tjno |
---|
290 | |
---|
291 | ! local |
---|
292 | |
---|
293 | real, parameter :: factor = 3.141592/180. |
---|
294 | real, dimension(nivbas) :: colo2, colno, colo3 |
---|
295 | real, dimension(nivbas) :: r_o2, r_no |
---|
296 | real, dimension(nivbas) :: p, t_o3_1, t_o3_2, t_o3_3 |
---|
297 | real, dimension(nivbas) :: bd1, bd2, bd3 |
---|
298 | real, dimension(6,3) :: sigma_o2, sigma_no_1, sigma_no_2 |
---|
299 | real, dimension(6,3) :: w_no_1, w_no_2 |
---|
300 | real, dimension(3) :: delta_lambda, i0, sigma_o3 |
---|
301 | real :: a, d, kq, n2 |
---|
302 | real :: sec, zcmt, dp |
---|
303 | integer :: iniv |
---|
304 | |
---|
305 | ! sections efficaces |
---|
306 | |
---|
307 | sigma_o2(:,1) = (/1.12e-23, 2.45e-23, 7.19e-23, & |
---|
308 | 3.04e-22, 1.75e-21, 1.11e-20/) |
---|
309 | sigma_o2(:,2) = (/1.35e-22, 2.99e-22, 7.33e-22, & |
---|
310 | 3.07e-21, 1.69e-20, 1.66e-19/) |
---|
311 | sigma_o2(:,3) = (/2.97e-22, 5.83e-22, 2.05e-21, & |
---|
312 | 8.19e-21, 4.80e-20, 2.66e-19/) |
---|
313 | |
---|
314 | sigma_no_1(:,1) = (/0.00e+00, 1.32e-18, 6.35e-19, & |
---|
315 | 7.09e-19, 2.18e-19, 4.67e-19/) |
---|
316 | sigma_no_1(:,2) = (/0.00e+00, 0.00e+00, 3.05e-21, & |
---|
317 | 5.76e-19, 2.29e-18, 2.21e-18/) |
---|
318 | sigma_no_1(:,3) = (/1.80e-18, 1.50e-18, 5.01e-19, & |
---|
319 | 7.20e-20, 6.72e-20, 1.49e-21/) |
---|
320 | |
---|
321 | sigma_no_2(:,1) = (/0.00e+00, 4.41e-17, 4.45e-17, & |
---|
322 | 4.50e-17, 2.94e-17, 4.35e-17/) |
---|
323 | sigma_no_2(:,2) = (/0.00e+00, 0.00e+00, 3.20e-21, & |
---|
324 | 5.71e-17, 9.09e-17, 6.00e-17/) |
---|
325 | sigma_no_2(:,3) = (/1.40e-16, 1.52e-16, 7.00e-17, & |
---|
326 | 2.83e-17, 2.73e-17, 6.57e-18/) |
---|
327 | |
---|
328 | ! facteurs de ponderation |
---|
329 | |
---|
330 | w_no_1(:,1) = (/0.00e+00, 5.12e-02, 1.36e-01, & |
---|
331 | 1.65e-01, 1.41e-01, 4.50e-02/) |
---|
332 | w_no_1(:,2) = (/0.00e+00, 0.00e+00, 1.93e-03, & |
---|
333 | 9.73e-02, 9.75e-02, 3.48e-02/) |
---|
334 | w_no_1(:,3) = (/4.50e-02, 1.80e-01, 2.25e-01, & |
---|
335 | 2.25e-01, 1.80e-01, 4.50e-02/) |
---|
336 | |
---|
337 | w_no_2(:,1) = (/0.00e+00, 5.68e-03, 1.52e-02, & |
---|
338 | 1.83e-02, 1.57e-02, 5.00e-03/) |
---|
339 | w_no_2(:,2) = (/0.00e+00, 0.00e+00, 2.14e-04, & |
---|
340 | 1.08e-02, 1.08e-02, 3.86e-03/) |
---|
341 | w_no_2(:,3) = (/5.00e-03, 2.00e-02, 2.50e-02, & |
---|
342 | 2.50e-02, 2.00e-02, 5.00e-03/) |
---|
343 | |
---|
344 | ! largeurs spectrales pour les trois bandes considerees (nm) |
---|
345 | |
---|
346 | delta_lambda = (/2.3, 1.5, 1.5/) |
---|
347 | |
---|
348 | ! flux pour les trois bandes considerees (s-1 cm-2 nm-1) |
---|
349 | |
---|
350 | i0 = (/3.98e+11, 2.21e+11, 2.30e+11/) |
---|
351 | |
---|
352 | ! sections efficaces de l'ozone pour les trois bandes |
---|
353 | ! considerees (cm2) d'apres wmo 1985 |
---|
354 | |
---|
355 | sigma_o3 = (/4.6e-19, 6.7e-19, 7.1e-19/) |
---|
356 | |
---|
357 | ! zcmt: gas constant/boltzmann constant*g |
---|
358 | |
---|
359 | zcmt = 287.0/(1.38e-23*9.81) |
---|
360 | |
---|
361 | ! d: taux de predissociation spontanee (s-1) |
---|
362 | |
---|
363 | d = 1.65e+09 |
---|
364 | |
---|
365 | ! a: taux d'emission spontanee (s-1) |
---|
366 | |
---|
367 | a = 5.1e+07 |
---|
368 | |
---|
369 | ! kq: quenching rate constant (cm3 s-1) |
---|
370 | |
---|
371 | kq = 1.5e-09 |
---|
372 | |
---|
373 | ! rapports de melange |
---|
374 | |
---|
375 | r_no(:) = c_no(:)/hnm(:) |
---|
376 | r_o2(:) = c_o2(:)/hnm(:) |
---|
377 | |
---|
378 | !============================================================= |
---|
379 | ! calcul des colonnes de o2 et no |
---|
380 | !============================================================= |
---|
381 | |
---|
382 | ! premier niveau du modele (mol/cm2) |
---|
383 | |
---|
384 | colo2(1) = 6.8e+05*c_o2(1) |
---|
385 | colno(1) = 6.8e+05*c_no(1) |
---|
386 | colo3(1) = o3t(1) |
---|
387 | |
---|
388 | ! cas general |
---|
389 | |
---|
390 | do iniv = 2,nivbas |
---|
391 | dp = (pm(iniv) - pm(iniv - 1))*100. |
---|
392 | dp = max(dp, 0.) |
---|
393 | colo2(iniv) = colo2(iniv - 1) & |
---|
394 | + zcmt*(r_o2(iniv-1) + r_o2(iniv))*0.5*dp*1.e-4 |
---|
395 | colno(iniv) = colno(iniv - 1) & |
---|
396 | + zcmt*(r_no(iniv-1) + r_no(iniv))*0.5*dp*1.e-4 |
---|
397 | colo3(iniv) = o3t(iniv) |
---|
398 | end do |
---|
399 | |
---|
400 | !============================================================= |
---|
401 | ! boucle sur les niveaux |
---|
402 | !============================================================= |
---|
403 | |
---|
404 | do iniv = 1,nivbas |
---|
405 | if (sza <= 89.0) then |
---|
406 | |
---|
407 | sec = 1./cos(sza*factor) |
---|
408 | |
---|
409 | colo2(iniv) = colo2(iniv)*sec |
---|
410 | colno(iniv) = colno(iniv)*sec |
---|
411 | colo3(iniv) = colo3(iniv)*sec |
---|
412 | |
---|
413 | ! facteurs de transmission de l'ozone |
---|
414 | |
---|
415 | t_o3_1(iniv) = exp(-sigma_o3(1)*colo3(iniv)) |
---|
416 | t_o3_2(iniv) = exp(-sigma_o3(2)*colo3(iniv)) |
---|
417 | t_o3_3(iniv) = exp(-sigma_o3(3)*colo3(iniv)) |
---|
418 | |
---|
419 | ! calcul de la probabilite de predissociation |
---|
420 | |
---|
421 | n2 = hnm(iniv)*0.78 |
---|
422 | p(iniv) = d/(a + d + kq*n2) |
---|
423 | |
---|
424 | end if |
---|
425 | end do |
---|
426 | |
---|
427 | ! calcul proprement dit, pour les 3 bandes |
---|
428 | |
---|
429 | do iniv = 1,nivbas |
---|
430 | if (sza <= 89.0) then |
---|
431 | |
---|
432 | bd1(iniv) = delta_lambda(1)*i0(1)*t_o3_1(iniv)*p(iniv)*( & |
---|
433 | exp(-sigma_o2(1,1)*colo2(iniv)) & |
---|
434 | *(w_no_1(1,1)*sigma_no_1(1,1)*exp(-sigma_no_1(1,1)*colno(iniv)) & |
---|
435 | + w_no_2(1,1)*sigma_no_2(1,1)*exp(-sigma_no_2(1,1)*colno(iniv))) & |
---|
436 | + exp(-sigma_o2(2,1)*colo2(iniv)) & |
---|
437 | *(w_no_1(2,1)*sigma_no_1(2,1)*exp(-sigma_no_1(2,1)*colno(iniv)) & |
---|
438 | + w_no_2(2,1)*sigma_no_2(2,1)*exp(-sigma_no_2(2,1)*colno(iniv))) & |
---|
439 | + exp(-sigma_o2(3,1)*colo2(iniv)) & |
---|
440 | *(w_no_1(3,1)*sigma_no_1(3,1)*exp(-sigma_no_1(3,1)*colno(iniv)) & |
---|
441 | + w_no_2(3,1)*sigma_no_2(3,1)*exp(-sigma_no_2(3,1)*colno(iniv))) & |
---|
442 | + exp(-sigma_o2(4,1)*colo2(iniv)) & |
---|
443 | *(w_no_1(4,1)*sigma_no_1(4,1)*exp(-sigma_no_1(4,1)*colno(iniv)) & |
---|
444 | + w_no_2(4,1)*sigma_no_2(4,1)*exp(-sigma_no_2(4,1)*colno(iniv))) & |
---|
445 | + exp(-sigma_o2(5,1)*colo2(iniv)) & |
---|
446 | *(w_no_1(5,1)*sigma_no_1(5,1)*exp(-sigma_no_1(5,1)*colno(iniv)) & |
---|
447 | + w_no_2(5,1)*sigma_no_2(5,1)*exp(-sigma_no_2(5,1)*colno(iniv))) & |
---|
448 | + exp(-sigma_o2(6,1)*colo2(iniv)) & |
---|
449 | *(w_no_1(6,1)*sigma_no_1(6,1)*exp(-sigma_no_1(6,1)*colno(iniv)) & |
---|
450 | + w_no_2(6,1)*sigma_no_2(6,1)*exp(-sigma_no_2(6,1)*colno(iniv))) & |
---|
451 | ) |
---|
452 | |
---|
453 | bd2(iniv) = delta_lambda(2)*i0(2)*t_o3_2(iniv)*p(iniv)*( & |
---|
454 | exp(-sigma_o2(1,2)*colo2(iniv)) & |
---|
455 | *(w_no_1(1,2)*sigma_no_1(1,2)*exp(-sigma_no_1(1,2)*colno(iniv)) & |
---|
456 | + w_no_2(1,2)*sigma_no_2(1,2)*exp(-sigma_no_2(1,2)*colno(iniv))) & |
---|
457 | + exp(-sigma_o2(2,2)*colo2(iniv)) & |
---|
458 | *(w_no_1(2,2)*sigma_no_1(2,2)*exp(-sigma_no_1(2,2)*colno(iniv)) & |
---|
459 | + w_no_2(2,2)*sigma_no_2(2,2)*exp(-sigma_no_2(2,2)*colno(iniv))) & |
---|
460 | + exp(-sigma_o2(3,2)*colo2(iniv)) & |
---|
461 | *(w_no_1(3,2)*sigma_no_1(3,2)*exp(-sigma_no_1(3,2)*colno(iniv)) & |
---|
462 | + w_no_2(3,2)*sigma_no_2(3,2)*exp(-sigma_no_2(3,2)*colno(iniv))) & |
---|
463 | + exp(-sigma_o2(4,2)*colo2(iniv)) & |
---|
464 | *(w_no_1(4,2)*sigma_no_1(4,2)*exp(-sigma_no_1(4,2)*colno(iniv)) & |
---|
465 | + w_no_2(4,2)*sigma_no_2(4,2)*exp(-sigma_no_2(4,2)*colno(iniv))) & |
---|
466 | + exp(-sigma_o2(5,2)*colo2(iniv)) & |
---|
467 | *(w_no_1(5,2)*sigma_no_1(5,2)*exp(-sigma_no_1(5,2)*colno(iniv)) & |
---|
468 | + w_no_2(5,2)*sigma_no_2(5,2)*exp(-sigma_no_2(5,2)*colno(iniv))) & |
---|
469 | + exp(-sigma_o2(6,2)*colo2(iniv)) & |
---|
470 | *(w_no_1(6,2)*sigma_no_1(6,2)*exp(-sigma_no_1(6,2)*colno(iniv)) & |
---|
471 | + w_no_2(6,2)*sigma_no_2(6,2)*exp(-sigma_no_2(6,2)*colno(iniv))) & |
---|
472 | ) |
---|
473 | |
---|
474 | bd3(iniv) = delta_lambda(3)*i0(3)*t_o3_3(iniv)*p(iniv)*( & |
---|
475 | exp(-sigma_o2(1,3)*colo2(iniv)) & |
---|
476 | *(w_no_1(1,3)*sigma_no_1(1,3)*exp(-sigma_no_1(1,3)*colno(iniv)) & |
---|
477 | + w_no_2(1,3)*sigma_no_2(1,3)*exp(-sigma_no_2(1,3)*colno(iniv))) & |
---|
478 | + exp(-sigma_o2(2,3)*colo2(iniv)) & |
---|
479 | *(w_no_1(2,3)*sigma_no_1(2,3)*exp(-sigma_no_1(2,3)*colno(iniv)) & |
---|
480 | + w_no_2(2,3)*sigma_no_2(2,3)*exp(-sigma_no_2(2,3)*colno(iniv))) & |
---|
481 | + exp(-sigma_o2(3,3)*colo2(iniv)) & |
---|
482 | *(w_no_1(3,3)*sigma_no_1(3,3)*exp(-sigma_no_1(3,3)*colno(iniv)) & |
---|
483 | + w_no_2(3,3)*sigma_no_2(3,3)*exp(-sigma_no_2(3,3)*colno(iniv))) & |
---|
484 | + exp(-sigma_o2(4,3)*colo2(iniv)) & |
---|
485 | *(w_no_1(4,3)*sigma_no_1(4,3)*exp(-sigma_no_1(4,3)*colno(iniv)) & |
---|
486 | + w_no_2(4,3)*sigma_no_2(4,3)*exp(-sigma_no_2(4,3)*colno(iniv))) & |
---|
487 | + exp(-sigma_o2(5,3)*colo2(iniv)) & |
---|
488 | *(w_no_1(5,3)*sigma_no_1(5,3)*exp(-sigma_no_1(5,3)*colno(iniv)) & |
---|
489 | + w_no_2(5,3)*sigma_no_2(5,3)*exp(-sigma_no_2(5,3)*colno(iniv))) & |
---|
490 | + exp(-sigma_o2(6,3)*colo2(iniv)) & |
---|
491 | *(w_no_1(6,3)*sigma_no_1(6,3)*exp(-sigma_no_1(6,3)*colno(iniv)) & |
---|
492 | + w_no_2(6,3)*sigma_no_2(6,3)*exp(-sigma_no_2(6,3)*colno(iniv))) & |
---|
493 | ) |
---|
494 | |
---|
495 | tjno(iniv) = bd1(iniv) + bd2(iniv) + bd3(iniv) |
---|
496 | |
---|
497 | else ! night |
---|
498 | tjno(iniv) = 0. |
---|
499 | end if |
---|
500 | end do |
---|
501 | |
---|
502 | end subroutine jno |
---|
503 | |
---|
504 | function indexchim(traceurname) |
---|
505 | |
---|
506 | use tracer_h, only: noms, nqtot, is_chim |
---|
507 | |
---|
508 | implicit none |
---|
509 | |
---|
510 | !======================================================================= |
---|
511 | ! Get index of tracer in chemistry with its name using traceur tab |
---|
512 | ! |
---|
513 | ! version: April 2019 - Yassin Jaziri (update September 2020) |
---|
514 | !======================================================================= |
---|
515 | |
---|
516 | ! input/output |
---|
517 | |
---|
518 | integer :: indexchim ! index of tracer in chemistry |
---|
519 | character(len=*) :: traceurname ! name of traceur to find |
---|
520 | |
---|
521 | ! local variables |
---|
522 | |
---|
523 | integer :: iq, iesp |
---|
524 | |
---|
525 | |
---|
526 | iesp = 0 |
---|
527 | indexchim = 0 |
---|
528 | do iq = 1,nqtot |
---|
529 | if (is_chim(iq)==1) then |
---|
530 | iesp = iesp + 1 |
---|
531 | if (trim(traceurname)==trim(noms(iq))) then |
---|
532 | indexchim = iesp |
---|
533 | exit |
---|
534 | end if |
---|
535 | end if |
---|
536 | end do |
---|
537 | |
---|
538 | if (indexchim==0) then |
---|
539 | print*, 'Error indexchim: wrong traceur name ', trim(traceurname) |
---|
540 | print*, 'Check if the specie ', trim(traceurname),' is in traceur.def' |
---|
541 | print*, 'Check if the specie ', trim(traceurname),' is a chemical species' |
---|
542 | print*, '(option is_chim = 1)' |
---|
543 | call abort |
---|
544 | end if |
---|
545 | |
---|
546 | return |
---|
547 | end function indexchim |
---|
548 | |
---|
549 | end module chimiedata_h |
---|