source: trunk/LMDZ.GENERIC/libf/aeronostd/chimiedata_h.F90 @ 3529

Last change on this file since 3529 was 3309, checked in by yjaziri, 7 months ago

GENERIC PCM:

  • Cosmetic + clarifying some variables and comments in chemistry
  • add controle variable for actinique UV flux in photochemistry and output surface UV flux in diagspecUV.nc

YJ

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