source: trunk/LMDZ.MARS/libf/aeronomars/jthermcalc_e107.F

Last change on this file was 3726, checked in by emillour, 3 months ago

Mars PCM:
Turn "callkeys.h" into module "callkeys_mod.F90"
EM

File size: 37.5 KB
RevLine 
[3464]1      module jthermcalc_e107_mod
2     
3      implicit none
4     
5      contains
6
[705]7c**********************************************************************
8
9      subroutine jthermcalc_e107
[1266]10     $     (ig,nlayer,chemthermod,rm,nesptherm,tx,iz,zenit,zday)
[705]11
12
13c     feb 2002        fgg           first version
14c     nov 2002        fgg           second version
15c
16c modified from paramhr.F
17c MAC July 2003
18c**********************************************************************
19
[1266]20      use param_v4_h, only: jfotsout,crscabsi2,
21     .    c1_16,c17_24,c25_29,c30_31,c32,c33,c34,c35,c36,
22     .    co2crsc195,co2crsc295,t0,
23     .    jabsifotsintpar,ninter,nz2,
24     .    nabs,e107,date_e107,e107_tab,
25     .    coefit0,coefit1,coefit2,coefit3,coefit4
[2429]26      use comsaison_h, only: dist_sol
[3464]27      use jthermcalc_util, only: column, interfast
[3726]28      use callkeys_mod, only: solvarmod, fixed_euv_value
[1266]29
[705]30      implicit none
31
32c     input and output variables
33
[2429]34      integer,intent(in) :: ig,nlayer
35      integer,intent(in) :: chemthermod
36      integer,intent(in) :: nesptherm          !Number of species considered
37      real,intent(in) :: rm(nlayer,nesptherm)  !Densities (cm-3)
38      real,intent(in) :: tx(nlayer)            !temperature
39      real,intent(in) :: zenit         !SZA
40      real,intent(in) :: iz(nlayer)    !Local altitude
41      real,intent(in) :: zday          !Martian day after Ls=0
[705]42
[2429]43! NB: no output variable! (computed jfotsout() is in module param_v4_h)
[705]44
45c    local parameters and variables
46
[1266]47      real       co2colx(nlayer)              !column density of CO2 (cm^-2)
48      real       o2colx(nlayer)               !column density of O2(cm^-2)
49      real       o3pcolx(nlayer)              !column density of O(3P)(cm^-2)
50      real       h2colx(nlayer)               !H2 column density (cm-2)
51      real       h2ocolx(nlayer)              !H2O column density (cm-2)
52      real       h2o2colx(nlayer)             !column density of H2O2(cm^-2)
53      real       o3colx(nlayer)               !O3 column density (cm-2)
54      real       n2colx(nlayer)               !N2 column density (cm-2)
55      real       ncolx(nlayer)                !N column density (cm-2)
56      real       nocolx(nlayer)               !NO column density (cm-2)
57      real       cocolx(nlayer)               !CO column density (cm-2)
58      real       hcolx(nlayer)                !H column density (cm-2)
59      real       no2colx(nlayer)              !NO2 column density (cm-2)
60      real       t2(nlayer)
61      real       coltemp(nlayer)
62      real       sigma(ninter,nlayer)
63      real       alfa(ninter,nlayer)
[705]64      real       realday
65     
66      integer    i,j,k,indexint                 !indexes
67      character  dn
68      integer    tinf,tsup
69
70
71
72c     variables used in interpolation
73
74      real*8      auxcoltab(nz2)
75      real*8      auxjco2(nz2)
76      real*8      auxjo2(nz2)
77      real*8      auxjo3p(nz2)
78      real*8      auxjh2o(nz2)
79      real*8      auxjh2(nz2)
80      real*8      auxjh2o2(nz2)
81      real*8      auxjo3(nz2)
82      real*8      auxjn2(nz2)
83      real*8      auxjn(nz2)
84      real*8      auxjno(nz2)
85      real*8      auxjco(nz2)
86      real*8      auxjh(nz2)
87      real*8      auxjno2(nz2)
[1266]88      real*8      wp(nlayer),wm(nlayer)
89      real*8      auxcolinp(nlayer)
90      integer     auxind(nlayer)
[705]91      integer     auxi
92      integer     ind
[1266]93      real*8      cortemp(nlayer)
[705]94
95      real*8      limdown                      !limits for interpolation
96      real*8      limup                        !        ""
97
98      !!!ATTENTION. Here i_co2 has to have the same value than in chemthermos.F90
99      !!!If the value is changed there, if has to be changed also here !!!!
100      integer,parameter :: i_co2=1
101
102
103c*************************PROGRAM STARTS*******************************
104     
105      if(zenit.gt.140.) then
106         dn='n'
107         else
108         dn='d'
109      end if
110      if(dn.eq.'n') then
111        return
112      endif
113     
114      !Initializing the photoabsorption coefficients
115      jfotsout(:,:,:)=0.
116
117      !Auxiliar temperature to take into account the temperature dependence
118      !of CO2 cross section
[1266]119      do i=1,nlayer
[705]120         t2(i)=tx(i)
121         if(t2(i).lt.195.0) t2(i)=195.0
122         if(t2(i).gt.295.0) t2(i)=295.0
123      end do
124
125      !Calculation of column amounts
[1266]126      call column(ig,nlayer,chemthermod,rm,nesptherm,tx,iz,zenit,
[705]127     $     co2colx,o2colx,o3pcolx,h2colx,h2ocolx,
128     $     h2o2colx,o3colx,n2colx,ncolx,nocolx,cocolx,hcolx,no2colx)
129
130      !Auxiliar column to include the temperature dependence
131      !of CO2 cross section
[1266]132      coltemp(nlayer)=co2colx(nlayer)*abs(t2(nlayer)-t0(nlayer))
133      do i=nlayer-1,1,-1
[705]134        coltemp(i)=!coltemp(i+1)+     PQ SE ELIMINA? REVISAR
135     $         ( rm(i,i_co2) + rm(i+1,i_co2) ) * 0.5
136     $         * 1e5 * (iz(i+1)-iz(i)) * abs(t2(i)-t0(i))
137      end do
138     
139      !Calculation of CO2 cross section at temperature t0(i)
[1266]140      do i=1,nlayer
[705]141         do indexint=24,32
142           sigma(indexint,i)=co2crsc195(indexint-23)
143           alfa(indexint,i)=((co2crsc295(indexint-23)
144     $          /sigma(indexint,i))-1.)/(295.-t0(i))
145        end do
146      end do
147
[1684]148      if (solvarmod==0) then
149        e107=fixed_euv_value
150      else
151        !E10.7 for the day: linear interpolation to tabulated values
152        realday=mod(zday,669.)
153        if(realday.lt.date_e107(1)) then
[705]154         e107=e107_tab(1)
[1684]155        else if(realday.ge.date_e107(669)) then
[762]156         e107=e107_tab(669)   
[1684]157        else if(realday.ge.date_e107(1).and.
[762]158     $        realday.lt.date_e107(669)) then
159         do i=1,668
[705]160            if(realday.ge.date_e107(i).and.
161     $           realday.lt.date_e107(i+1)) then
162               tinf=i
163               tsup=i+1
[762]164               e107=e107_tab(tinf)+(realday-date_e107(tinf))*
[705]165     $              (e107_tab(tsup)-e107_tab(tinf))
166            endif
167         enddo
[1684]168        endif
169      endif ! of if (solvarmod==0)
[762]170
[705]171      !Photoabsorption coefficients at TOA as a function of E10.7
172      do j=1,nabs
173         do indexint=1,ninter
[1266]174            jfotsout(indexint,j,nlayer)=coefit0(indexint,j)+
[705]175     $           coefit1(indexint,j)*e107+coefit2(indexint,j)*e107**2+
176     $           coefit3(indexint,j)*e107**3+coefit4(indexint,j)*e107**4
177         enddo
178      enddo
179! Interpolation to the tabulated photoabsorption coefficients for each species
180! in each spectral interval
181
182
183c     auxcolinp-> Actual atmospheric column
184c     auxj*-> Tabulated photoabsorption coefficients
185c     auxcoltab-> Tabulated atmospheric columns
186
187ccccccccccccccccccccccccccccccc
188c     0.1,5.0 (int 1)
189c
190c     Absorption by:
191c     CO2, O2, O, H2, N
192ccccccccccccccccccccccccccccccc
193
194c     Input atmospheric column
195      indexint=1
[1266]196      do i=1,nlayer
197         auxcolinp(nlayer-i+1) = co2colx(i)*crscabsi2(1,indexint) +
[705]198     $        o2colx(i)*crscabsi2(2,indexint) +
199     $        o3pcolx(i)*crscabsi2(3,indexint) +
200     $        h2colx(i)*crscabsi2(5,indexint) +
201     $        ncolx(i)*crscabsi2(9,indexint)
202      end do
203      limdown=1.e-20
204      limup=1.e26
205
206
207c     Interpolations
208
209      do i=1,nz2
210         auxi = nz2-i+1
211         !CO2 tabulated coefficient
212         auxjco2(i) = jabsifotsintpar(auxi,1,indexint)
213         !O2 tabulated coefficient
214         auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
215         !O3p tabulated coefficient
216         auxjo3p(i) = jabsifotsintpar(auxi,3,indexint)
217         !H2 tabulated coefficient
218         auxjh2(i) = jabsifotsintpar(auxi,5,indexint)
219         !Tabulated column
220         auxcoltab(i) = c1_16(auxi,indexint)
221      enddo
222      !Only if chemthermod.ge.2
223      !N tabulated coefficient
224      if(chemthermod.ge.2) then
225         do i=1,nz2
226            auxjn(i) = jabsifotsintpar(nz2-i+1,9,indexint)
227         enddo
228      endif
229
230      call interfast
[1266]231     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]232     
[1266]233      do i=1,nlayer
[705]234         ind=auxind(i)
[1266]235         auxi=nlayer-i+1
[705]236         !CO2 interpolated coefficient
[1266]237         jfotsout(indexint,1,auxi) = jfotsout(indexint,1,nlayer) *
[705]238     $        (wm(i)*auxjco2(ind+1) + wp(i)*auxjco2(ind))
239         !O2 interpolated coefficient
[1266]240         jfotsout(indexint,2,auxi) = jfotsout(indexint,2,nlayer) *
[705]241     $        (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind))
242         !O3p interpolated coefficient
[1266]243         jfotsout(indexint,3,auxi) = jfotsout(indexint,3,nlayer) *
[705]244     $        (wm(i)*auxjo3p(ind+1) + wp(i)*auxjo3p(ind))
245         !H2 interpolated coefficient
[2808]246         jfotsout(indexint,5,auxi) = jfotsout(indexint,5,nlayer) *
[705]247     $        (wm(i)*auxjh2(ind+1) + wp(i)*auxjh2(ind))
248      enddo
249      !Only if chemthermod.ge.2
250      !N interpolated coefficient
251      if(chemthermod.ge.2) then
[1266]252         do i=1,nlayer
[705]253            ind=auxind(i)
[1266]254            jfotsout(indexint,9,nlayer-i+1) = 
255     $           jfotsout(indexint,9,nlayer) *
[705]256     $           (wm(i)*auxjn(ind+1) + wp(i)*auxjn(ind))
257         enddo
258      endif
259         
260
261c     End interval 1
262
263
264ccccccccccccccccccccccccccccccc
265c     5-80.5nm (int 2-15)
266c
267c     Absorption by:
268c     CO2, O2, O, H2, N2, N,
269c     NO, CO, H, NO2
270ccccccccccccccccccccccccccccccc
271
272c     Input atmospheric column
273      do indexint=2,15
[1266]274         do i=1,nlayer
275            auxcolinp(nlayer-i+1) = co2colx(i)*crscabsi2(1,indexint)+
[705]276     $           o2colx(i)*crscabsi2(2,indexint)+
277     $           o3pcolx(i)*crscabsi2(3,indexint)+
278     $           h2colx(i)*crscabsi2(5,indexint)+
279     $           n2colx(i)*crscabsi2(8,indexint)+
280     $           ncolx(i)*crscabsi2(9,indexint)+
281     $           nocolx(i)*crscabsi2(10,indexint)+
282     $           cocolx(i)*crscabsi2(11,indexint)+
283     $           hcolx(i)*crscabsi2(12,indexint)+
284     $           no2colx(i)*crscabsi2(13,indexint)
285         end do
286
287c     Interpolations
288
289         do i=1,nz2
290            auxi = nz2-i+1
291            !O2 tabulated coefficient
292            auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
293            !O3p tabulated coefficient
294            auxjo3p(i) = jabsifotsintpar(auxi,3,indexint)
295            !CO2 tabulated coefficient
296            auxjco2(i) = jabsifotsintpar(auxi,1,indexint)
297            !H2 tabulated coefficient
298            auxjh2(i) = jabsifotsintpar(auxi,5,indexint)
299            !N2 tabulated coefficient
300            auxjn2(i) = jabsifotsintpar(auxi,8,indexint)
301            !CO tabulated coefficient
302            auxjco(i) = jabsifotsintpar(auxi,11,indexint)
303            !H tabulated coefficient
304            auxjh(i) = jabsifotsintpar(auxi,12,indexint)
305            !tabulated column
306            auxcoltab(i) = c1_16(auxi,indexint)
307         enddo
308         !Only if chemthermod.ge.2
309         if(chemthermod.ge.2) then
310            do i=1,nz2
311               auxi = nz2-i+1
312               !N tabulated coefficient
313               auxjn(i) = jabsifotsintpar(auxi,9,indexint)
314               !NO tabulated coefficient
315               auxjno(i) = jabsifotsintpar(auxi,10,indexint)
316               !NO2 tabulated coefficient
317               auxjno2(i) = jabsifotsintpar(auxi,13,indexint)
318            enddo
319         endif
320
[1266]321          call interfast(wm,wp,auxind,auxcolinp,nlayer,
[705]322     $        auxcoltab,nz2,limdown,limup)
[2615]323
[1266]324          do i=1,nlayer
[705]325             ind=auxind(i)
[1266]326             auxi = nlayer-i+1
[705]327             !O2 interpolated coefficient
328             jfotsout(indexint,2,auxi) =
[1266]329     $            jfotsout(indexint,2,nlayer) *
[705]330     $            (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind))
331             !O3p interpolated coefficient
332             jfotsout(indexint,3,auxi) =
[1266]333     $            jfotsout(indexint,3,nlayer) *
[705]334     $            (wm(i)*auxjo3p(ind+1) + wp(i)*auxjo3p(ind))
335             !CO2 interpolated coefficient
336             jfotsout(indexint,1,auxi) =
[1266]337     $            jfotsout(indexint,1,nlayer) *
[705]338     $            (wm(i)*auxjco2(ind+1) + wp(i)*auxjco2(ind))
339             !H2 interpolated coefficient
340             jfotsout(indexint,5,auxi) =
[1266]341     $            jfotsout(indexint,5,nlayer) *
[705]342     $            (wm(i)*auxjh2(ind+1) + wp(i)*auxjh2(ind))
343             !N2 interpolated coefficient
344             jfotsout(indexint,8,auxi) =
[1266]345     $            jfotsout(indexint,8,nlayer) *
[705]346     $            (wm(i)*auxjn2(ind+1) + wp(i)*auxjn2(ind))
347             !CO interpolated coefficient
348             jfotsout(indexint,11,auxi) =
[1266]349     $            jfotsout(indexint,11,nlayer) *
[705]350     $            (wm(i)*auxjco(ind+1) + wp(i)*auxjco(ind))
351             !H interpolated coefficient
352             jfotsout(indexint,12,auxi) =
[1266]353     $            jfotsout(indexint,12,nlayer) *
[705]354     $            (wm(i)*auxjh(ind+1) + wp(i)*auxjh(ind))
355          enddo
356          !Only if chemthermod.ge.2
357          if(chemthermod.ge.2) then
[1266]358             do i=1,nlayer
[705]359                ind=auxind(i)
[1266]360                auxi = nlayer-i+1
[705]361                !N interpolated coefficient
362                jfotsout(indexint,9,auxi) =
[1266]363     $               jfotsout(indexint,9,nlayer) *
[705]364     $               (wm(i)*auxjn(ind+1) + wp(i)*auxjn(ind))
365                !NO interpolated coefficient
366                jfotsout(indexint,10,auxi)=
[1266]367     $               jfotsout(indexint,10,nlayer) *
[705]368     $               (wm(i)*auxjno(ind+1) + wp(i)*auxjno(ind))
369                !NO2 interpolated coefficient
370                jfotsout(indexint,13,auxi)=
[1266]371     $               jfotsout(indexint,13,nlayer) *
[705]372     $               (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind))
373             enddo
374          endif   
375      end do
376
377c     End intervals 2-15
378
379
380ccccccccccccccccccccccccccccccc
381c     80.6-90.8nm (int16)
382c
383c     Absorption by:
384c     CO2, O2, O, N2, N, NO,
385c     CO, H, NO2
386ccccccccccccccccccccccccccccccc
387
388c     Input atmospheric column
389      indexint=16
[1266]390      do i=1,nlayer
391         auxcolinp(nlayer-i+1) = co2colx(i)*crscabsi2(1,indexint)+
[705]392     $        o2colx(i)*crscabsi2(2,indexint)+
393     $        o3pcolx(i)*crscabsi2(3,indexint)+
394     $        n2colx(i)*crscabsi2(8,indexint)+
395     $        ncolx(i)*crscabsi2(9,indexint)+
396     $        nocolx(i)*crscabsi2(10,indexint)+
397     $        cocolx(i)*crscabsi2(11,indexint)+
398     $        hcolx(i)*crscabsi2(12,indexint)+
399     $        no2colx(i)*crscabsi2(13,indexint)
400      end do
401
402c     Interpolations
403
404      do i=1,nz2
405         auxi = nz2-i+1
406         !O2 tabulated coefficient
407         auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
408         !CO2 tabulated coefficient
409         auxjco2(i) = jabsifotsintpar(auxi,1,indexint)
410         !O3p tabulated coefficient
411         auxjo3p(i) = jabsifotsintpar(auxi,3,indexint)
412         !N2 tabulated coefficient
413         auxjn2(i) = jabsifotsintpar(auxi,8,indexint)
414         !CO tabulated coefficient
415         auxjco(i) = jabsifotsintpar(auxi,11,indexint)
416         !H tabulated coefficient
417         auxjh(i) = jabsifotsintpar(auxi,12,indexint)
418         !NO2 tabulated coefficient
419         auxjno2(i) = jabsifotsintpar(auxi,13,indexint)
420         !Tabulated column
421         auxcoltab(i) = c1_16(auxi,indexint)
422      enddo
423      !Only if chemthermod.ge.2
424      if(chemthermod.ge.2) then
425         do i=1,nz2
426            auxi = nz2-i+1
427            !N tabulated coefficient
428            auxjn(i) = jabsifotsintpar(auxi,9,indexint)
429            !NO tabulated coefficient
430            auxjno(i) = jabsifotsintpar(auxi,10,indexint)
431            !NO2 tabulated coefficient
432            auxjno2(i) = jabsifotsintpar(auxi,13,indexint)
433         enddo
434      endif
435
436      call interfast
[1266]437     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]438
[1266]439      do i=1,nlayer
[705]440         ind=auxind(i)
[1266]441         auxi = nlayer-i+1
[705]442         !O2 interpolated coefficient
[1266]443         jfotsout(indexint,2,auxi) = jfotsout(indexint,2,nlayer) *
[705]444     $        (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind))
445         !CO2 interpolated coefficient
[1266]446         jfotsout(indexint,1,auxi) = jfotsout(indexint,1,nlayer) *
[705]447     $        (wm(i)*auxjco2(ind+1) + wp(i)*auxjco2(ind))
448         !O3p interpolated coefficient
[1266]449         jfotsout(indexint,3,auxi) = jfotsout(indexint,3,nlayer) *
[705]450     $        (wm(i)*auxjo3p(ind+1) + wp(i)*auxjo3p(ind))
451         !N2 interpolated coefficient
[1266]452         jfotsout(indexint,8,auxi) = jfotsout(indexint,8,nlayer) *
[705]453     $        (wm(i)*auxjn2(ind+1) + wp(i)*auxjn2(ind))
454         !CO interpolated coefficient
455         jfotsout(indexint,11,auxi) =
[1266]456     $        jfotsout(indexint,11,nlayer) *
[705]457     $        (wm(i)*auxjco(ind+1) + wp(i)*auxjco(ind))
458         !H interpolated coefficient
459         jfotsout(indexint,12,auxi) =
[1266]460     $        jfotsout(indexint,12,nlayer) *
[705]461     $        (wm(i)*auxjh(ind+1) + wp(i)*auxjh(ind))
462      enddo
463      !Only if chemthermod.ge.2
464      if(chemthermod.ge.2) then
[1266]465         do i=1,nlayer
[705]466            ind=auxind(i)
[1266]467            auxi = nlayer-i+1
[705]468            !N interpolated coefficient
469            jfotsout(indexint,9,auxi) =
[1266]470     $           jfotsout(indexint,9,nlayer) *
[705]471     $           (wm(i)*auxjn(ind+1) + wp(i)*auxjn(ind))
472            !NO interpolated coefficient
473            jfotsout(indexint,10,auxi) =
[1266]474     $           jfotsout(indexint,10,nlayer) *
[705]475     $           (wm(i)*auxjno(ind+1) + wp(i)*auxjno(ind))
476            !NO2 interpolated coefficient
477            jfotsout(indexint,13,auxi) =
[1266]478     $           jfotsout(indexint,13,nlayer) *
[705]479     $           (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind))
480         enddo
481      endif
482c     End interval 16
483
484
485ccccccccccccccccccccccccccccccc
486c     90.9-119.5nm (int 17-24)
487c
488c     Absorption by:
489c     CO2, O2, N2, NO, CO, NO2
490ccccccccccccccccccccccccccccccc
491
492c     Input column
493
[1266]494      do i=1,nlayer
495         auxcolinp(nlayer-i+1) = co2colx(i) + o2colx(i) + n2colx(i) +
[705]496     $        nocolx(i) + cocolx(i) + no2colx(i)
497      end do
498
499      do indexint=17,24
500
501c     Interpolations
502
503         do i=1,nz2
504            auxi = nz2-i+1
505            !CO2 tabulated coefficient
506            auxjco2(i) = jabsifotsintpar(auxi,1,indexint)
507            !O2 tabulated coefficient
508            auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
509            !N2 tabulated coefficient
510            auxjn2(i) = jabsifotsintpar(auxi,8,indexint)
511            !CO tabulated coefficient
512            auxjco(i) = jabsifotsintpar(auxi,11,indexint)           
513            !Tabulated column
514            auxcoltab(i) = c17_24(auxi)
515         enddo
516         !Only if chemthermod.ge.2
517         if(chemthermod.ge.2) then
518            do i=1,nz2
519               auxi = nz2-i+1
520               !NO tabulated coefficient
521               auxjno(i) = jabsifotsintpar(auxi,10,indexint)
522               !NO2 tabulated coefficient
523               auxjno2(i) = jabsifotsintpar(auxi,13,indexint)
524            enddo
525         endif
526
527         call interfast
[1266]528     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]529
[705]530         !Correction to include T variation of CO2 cross section
531         if(indexint.eq.24) then
[1266]532            do i=1,nlayer
533               auxi = nlayer-i+1
[705]534               if(sigma(indexint,auxi)*
535     $              alfa(indexint,auxi)*coltemp(auxi)
536     $              .lt.60.) then
537                  cortemp(i)=exp(-sigma(indexint,auxi)*
538     $                alfa(indexint,auxi)*coltemp(auxi))
539               else
540                  cortemp(i)=0.
541               end if
542            enddo
543         else
[1266]544            do i=1,nlayer
[705]545               cortemp(i)=1.
546            enddo
547         end if
[1266]548         do i=1,nlayer           
[705]549            ind=auxind(i)
[1266]550            auxi = nlayer-i+1
[705]551            !O2 interpolated coefficient
552            jfotsout(indexint,2,auxi) =
[1266]553     $           jfotsout(indexint,2,nlayer) *
[705]554     $           (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind)) *
555     $           cortemp(i)
556            !CO2 interpolated coefficient
557            jfotsout(indexint,1,auxi) =
[1266]558     $           jfotsout(indexint,1,nlayer) *
[705]559     $           (wm(i)*auxjco2(ind+1) + wp(i)*auxjco2(ind))
560     $           * cortemp(i)
561            if(indexint.eq.24) jfotsout(indexint,1,auxi)=
562     $           jfotsout(indexint,1,auxi)*
563     $           (1+alfa(indexint,auxi)*
564     $           (t2(auxi)-t0(auxi)))
565            !N2 interpolated coefficient
566            jfotsout(indexint,8,auxi) =
[1266]567     $           jfotsout(indexint,8,nlayer) *
[705]568     $           (wm(i)*auxjn2(ind+1) + wp(i)*auxjn2(ind)) *
569     $           cortemp(i)           
570            !CO interpolated coefficient
571            jfotsout(indexint,11,auxi) =
[1266]572     $           jfotsout(indexint,11,nlayer) *
[705]573     $           (wm(i)*auxjco(ind+1) + wp(i)*auxjco(ind)) *
574     $           cortemp(i)           
575         enddo
576         !Only if chemthermod.ge.2
577         if(chemthermod.ge.2) then
[1266]578            do i=1,nlayer
[705]579               ind=auxind(i)
[1266]580               auxi = nlayer-i+1
[705]581               !NO interpolated coefficient
582               jfotsout(indexint,10,auxi)=
[1266]583     $              jfotsout(indexint,10,nlayer) *
[705]584     $              (wm(i)*auxjno(ind+1) + wp(i)*auxjno(ind)) *
585     $              cortemp(i)
586               !NO2 interpolated coefficient
587               jfotsout(indexint,13,auxi)=
[1266]588     $              jfotsout(indexint,13,nlayer) *
[705]589     $              (wm(i)*auxjno2(ind+1)+ wp(i)*auxjno2(ind)) *
590     $              cortemp(i)
591            enddo
592         endif               
593      end do
594c     End intervals 17-24
595
596
597ccccccccccccccccccccccccccccccc
598c     119.6-167.0nm (int 25-29)
599c
600c     Absorption by:
601c     CO2, O2, H2O, H2O2, NO,
602c     CO, NO2
603ccccccccccccccccccccccccccccccc
604
605c     Input atmospheric column
606
[1266]607      do i=1,nlayer
608         auxcolinp(nlayer-i+1) = co2colx(i) + o2colx(i) + h2ocolx(i) +
[705]609     $        h2o2colx(i) + nocolx(i) + cocolx(i) + no2colx(i)
610      end do
611
612      do indexint=25,29
613
614c     Interpolations
615
616         do i=1,nz2
617            auxi = nz2-i+1
618            !CO2 tabulated coefficient
619            auxjco2(i) = jabsifotsintpar(auxi,1,indexint)
620            !O2 tabulated coefficient
621            auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
622            !H2O tabulated coefficient
623            auxjh2o(i) = jabsifotsintpar(auxi,4,indexint)
624            !H2O2 tabulated coefficient
625            auxjh2o2(i) = jabsifotsintpar(auxi,6,indexint)           
626            !CO tabulated coefficient
627            auxjco(i) = jabsifotsintpar(auxi,11,indexint)           
628            !Tabulated column
629            auxcoltab(i) = c25_29(auxi)
630         enddo
631         !Only if chemthermod.ge.2
632         if(chemthermod.ge.2) then
633            do i=1,nz2
634               auxi = nz2-i+1
635               !NO tabulated coefficient
636               auxjno(i) = jabsifotsintpar(auxi,10,indexint)
637               !NO2 tabulated coefficient
638               auxjno2(i) = jabsifotsintpar(auxi,13,indexint)
639            enddo
640         endif
641         call interfast
[1266]642     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]643
[1266]644         do i=1,nlayer
[705]645            ind=auxind(i)
[1266]646            auxi = nlayer-i+1
[705]647            !Correction to include T variation of CO2 cross section
648            if(sigma(indexint,auxi)*alfa(indexint,auxi)*
649     $           coltemp(auxi).lt.60.) then
650               cortemp(i)=exp(-sigma(indexint,auxi)*
651     $              alfa(indexint,auxi)*coltemp(auxi))
652            else
653               cortemp(i)=0.
654            end if
655            !CO2 interpolated coefficient
656            jfotsout(indexint,1,auxi) =
[1266]657     $           jfotsout(indexint,1,nlayer) *
[705]658     $           (wm(i)*auxjco2(ind+1) + wp(i)*auxjco2(ind)) *
659     $           cortemp(i) *
660     $           (1+alfa(indexint,auxi)*
661     $           (t2(auxi)-t0(auxi)))
662            !O2 interpolated coefficient
663            jfotsout(indexint,2,auxi) =
[1266]664     $           jfotsout(indexint,2,nlayer) *
[705]665     $           (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind)) *
666     $           cortemp(i)
667            !H2O interpolated coefficient
668            jfotsout(indexint,4,auxi) =
[1266]669     $           jfotsout(indexint,4,nlayer) *
[705]670     $           (wm(i)*auxjh2o(ind+1) + wp(i)*auxjh2o(ind)) *
671     $           cortemp(i)
672            !H2O2 interpolated coefficient
673            jfotsout(indexint,6,auxi) =
[1266]674     $           jfotsout(indexint,6,nlayer) *
[705]675     $           (wm(i)*auxjh2o2(ind+1) + wp(i)*auxjh2o2(ind)) *
676     $           cortemp(i)           
677            !CO interpolated coefficient
678            jfotsout(indexint,11,auxi) =
[1266]679     $           jfotsout(indexint,11,nlayer) *
[705]680     $           (wm(i)*auxjco(ind+1) + wp(i)*auxjco(ind)) *
681     $           cortemp(i)
682         enddo
683         !Only if chemthermod.ge.2
684         if(chemthermod.ge.2) then
[1266]685            do i=1,nlayer
[705]686               ind=auxind(i)
[1266]687               auxi = nlayer-i+1
[705]688               !NO interpolated coefficient
689               jfotsout(indexint,10,auxi)=
[1266]690     $              jfotsout(indexint,10,nlayer) *
[705]691     $              (wm(i)*auxjno(ind+1) + wp(i)*auxjno(ind)) *
692     $              cortemp(i)
693               !NO2 interpolated coefficient
694               jfotsout(indexint,13,auxi)=
[1266]695     $              jfotsout(indexint,13,nlayer) *
[705]696     $              (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind)) *
697     $              cortemp(i)
698            enddo
699         endif
700
701      end do
702
703c     End intervals 25-29
704
705
706cccccccccccccccccccccccccccccccc
707c     167.1-202.5nm (int 30-31)
708c   
709c     Absorption by:
710c     CO2, O2, H2O, H2O2, NO,
711c     NO2
712cccccccccccccccccccccccccccccccc
713
714c     Input atmospheric column
715
[1266]716      do i=1,nlayer
717         auxcolinp(nlayer-i+1) = co2colx(i) + o2colx(i) + h2ocolx(i) +
[705]718     $        h2o2colx(i) + nocolx(i) + no2colx(i)
719      end do
720
721c     Interpolation
722
723      do indexint=30,31
724
725         do i=1,nz2
726            auxi = nz2-i+1
727            !CO2 tabulated coefficient
728            auxjco2(i) = jabsifotsintpar(auxi,1,indexint)
729            !O2 tabulated coefficient
730            auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
731            !H2O tabulated coefficient
732            auxjh2o(i) = jabsifotsintpar(auxi,4,indexint)
733            !H2O2 tabulated coefficient
734            auxjh2o2(i) = jabsifotsintpar(auxi,6,indexint)           
735            !Tabulated column
736            auxcoltab(i) = c30_31(auxi)
737         enddo
738         !Only if chemthermod.ge.2
739         if(chemthermod.ge.2) then
740            do i=1,nz2
741               auxi = nz2-i+1
742               !NO tabulated coefficient
743               auxjno(i) = jabsifotsintpar(auxi,10,indexint)
744               !NO2 tabulated coefficient
745               auxjno2(i) = jabsifotsintpar(auxi,13,indexint)
746            enddo
747         endif
748
749         call interfast
[1266]750     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]751
[1266]752         do i=1,nlayer
[705]753            ind=auxind(i)
[1266]754            auxi = nlayer-i+1
[705]755            !Correction to include T variation of CO2 cross section
756            if(sigma(indexint,auxi)*alfa(indexint,auxi)*
757     $           coltemp(auxi).lt.60.) then
758               cortemp(i)=exp(-sigma(indexint,auxi)*
759     $              alfa(indexint,auxi)*coltemp(auxi))
760            else
761               cortemp(i)=0.
762            end if
763            !CO2 interpolated coefficient
764            jfotsout(indexint,1,auxi) =
[1266]765     $           jfotsout(indexint,1,nlayer) *
[705]766     $           (wm(i)*auxjco2(ind+1) + wp(i)*auxjco2(ind)) *
767     $           cortemp(i) *
768     $           (1+alfa(indexint,auxi)*
769     $           (t2(auxi)-t0(auxi)))
770            !O2 interpolated coefficient
771            jfotsout(indexint,2,auxi) =
[1266]772     $           jfotsout(indexint,2,nlayer) *
[705]773     $           (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind)) *
774     $           cortemp(i)
775            !H2O interpolated coefficient
776            jfotsout(indexint,4,auxi) =
[1266]777     $           jfotsout(indexint,4,nlayer) *
[705]778     $           (wm(i)*auxjh2o(ind+1) + wp(i)*auxjh2o(ind)) *
779     $           cortemp(i)
780            !H2O2 interpolated coefficient
781            jfotsout(indexint,6,auxi) =
[1266]782     $           jfotsout(indexint,6,nlayer) *
[705]783     $           (wm(i)*auxjh2o2(ind+1) + wp(i)*auxjh2o2(ind)) *
784     $           cortemp(i)           
785         enddo
786         !Only if chemthermod.ge.2
787         if(chemthermod.ge.2) then
[1266]788            do i=1,nlayer
[705]789               ind=auxind(i)
[1266]790               auxi = nlayer-i+1
[705]791               !NO interpolated coefficient
792               jfotsout(indexint,10,auxi)=
[1266]793     $              jfotsout(indexint,10,nlayer) *
[705]794     $              (wm(i)*auxjno(ind+1) +wp(i)*auxjno(ind)) *
795     $              cortemp(i)
796               !NO2 interpolated coefficient
797               jfotsout(indexint,13,auxi)=
[2808]798     $              jfotsout(indexint,13,nlayer) *
[705]799     $              (wm(i)*auxjno2(ind+1)+wp(i)*auxjno2(ind)) *
800     $              cortemp(i)
801            enddo
802         endif
803
804      end do
805
806c     End intervals 30-31
807
808
809ccccccccccccccccccccccccccccccc
810c     202.6-210.0nm (int 32)
811c
812c     Absorption by:
813c     CO2, O2, H2O2, NO, NO2
814ccccccccccccccccccccccccccccccc
815
816c     Input atmospheric column
817
818      indexint=32
[1266]819      do i=1,nlayer
820         auxcolinp(nlayer-i+1) =co2colx(i) + o2colx(i) + h2o2colx(i) +
[705]821     $        nocolx(i) + no2colx(i)
822      end do
823
824c     Interpolation
825
826      do i=1,nz2
827         auxi = nz2-i+1
828         !CO2 tabulated coefficient
829         auxjco2(i) = jabsifotsintpar(auxi,1,indexint)
830         !O2 tabulated coefficient
831         auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
832         !H2O2 tabulated coefficient
833         auxjh2o2(i) = jabsifotsintpar(auxi,6,indexint)         
834         !Tabulated column
835         auxcoltab(i) = c32(auxi)
836      enddo
837      !Only if chemthermod.ge.2
838      if(chemthermod.ge.2) then
839         do i=1,nz2
840            auxi = nz2-i+1
841            !NO tabulated coefficient
842            auxjno(i) = jabsifotsintpar(auxi,10,indexint)
843            !NO2 tabulated coefficient
844            auxjno2(i) = jabsifotsintpar(auxi,13,indexint)
845         enddo
846      endif
847      call interfast
[1266]848     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]849
[1266]850      do i=1,nlayer
[705]851         ind=auxind(i)
[1266]852         auxi = nlayer-i+1
[705]853         !Correction to include T variation of CO2 cross section
[1266]854         if(sigma(indexint,nlayer-i+1)*alfa(indexint,auxi)*
[705]855     $        coltemp(auxi).lt.60.) then
856            cortemp(i)=exp(-sigma(indexint,auxi)*
857     $           alfa(indexint,auxi)*coltemp(auxi))
858         else
859            cortemp(i)=0.
860         end if
861         !CO2 interpolated coefficient
862         jfotsout(indexint,1,auxi) =
[1266]863     $        jfotsout(indexint,1,nlayer) *
[705]864     $        (wm(i)*auxjco2(ind+1)+wp(i)*auxjco2(ind)) *
865     $        cortemp(i) *
866     $        (1+alfa(indexint,auxi)*
867     $        (t2(auxi)-t0(auxi)))
868         !O2 interpolated coefficient
869         jfotsout(indexint,2,auxi) =
[1266]870     $        jfotsout(indexint,2,nlayer) *
[705]871     $        (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind)) *
872     $        cortemp(i)
873         !H2O2 interpolated coefficient
874         jfotsout(indexint,6,auxi) =
[1266]875     $        jfotsout(indexint,6,nlayer) *
[705]876     $        (wm(i)*auxjh2o2(ind+1) + wp(i)*auxjh2o2(ind)) *
877     $        cortemp(i)         
878      enddo
879      !Only if chemthermod.ge.2
880      if(chemthermod.ge.2) then
[1266]881         do i=1,nlayer
882            auxi = nlayer-i+1
[705]883            ind=auxind(i)
884            !NO interpolated coefficient
885            jfotsout(indexint,10,auxi) =
[1266]886     $           jfotsout(indexint,10,nlayer) *
[705]887     $           (wm(i)*auxjno(ind+1) + wp(i)*auxjno(ind)) *
888     $           cortemp(i)
889           !NO2 interpolated coefficient
890            jfotsout(indexint,13,auxi) =
[1266]891     $           jfotsout(indexint,13,nlayer) *
[705]892     $           (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind)) *
893     $           cortemp(i)
894         enddo
895      endif
896
897c     End of interval 32
898
899
900ccccccccccccccccccccccccccccccc
901c     210.1-231.0nm (int 33)
902c     
903c     Absorption by:
904c     O2, H2O2, NO2
905ccccccccccccccccccccccccccccccc
906
907c     Input atmospheric column
908
909      indexint=33
[1266]910      do i=1,nlayer
911         auxcolinp(nlayer-i+1) = o2colx(i) + h2o2colx(i) + no2colx(i)
[705]912      end do
913
914c     Interpolation
915
916      do i=1,nz2
917         auxi = nz2-i+1
918         !O2 tabulated coefficient
919         auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
920         !H2O2 tabulated coefficient
921         auxjh2o2(i) = jabsifotsintpar(auxi,6,indexint)
922         !Tabulated column
923         auxcoltab(i) = c33(auxi)
924      enddo
925      !Only if chemthermod.ge.2
926      if(chemthermod.ge.2) then
927         do i=1,nz2
928            !NO2 tabulated coefficient
929            auxjno2(i) = jabsifotsintpar(nz2-i+1,13,indexint)
930         enddo
931      endif
932      call interfast
[1266]933     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]934
[1266]935      do i=1,nlayer
[705]936         ind=auxind(i)
[1266]937         auxi = nlayer-i+1
[705]938         !O2 interpolated coefficient
[1266]939         jfotsout(indexint,2,auxi) = jfotsout(indexint,2,nlayer) *
[705]940     $        (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind))
941         !H2O2 interpolated coefficient
[1266]942         jfotsout(indexint,6,auxi) = jfotsout(indexint,6,nlayer) *
[705]943     $        (wm(i)*auxjh2o2(ind+1) + wp(i)*auxjh2o2(ind))
944      enddo
945      !Only if chemthermod.ge.2
946      if(chemthermod.ge.2) then
[1266]947         do i=1,nlayer
[705]948            ind=auxind(i)
949            !NO2 interpolated coefficient
[1266]950            jfotsout(indexint,13,nlayer-i+1) =
951     $           jfotsout(indexint,13,nlayer) *
[705]952     $           (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind))
953         enddo
954      endif
955
956c     End of interval 33
957
958
959ccccccccccccccccccccccccccccccc
960c     231.1-240.0nm (int 34)
961c
962c     Absorption by:
963c     O2, H2O2, O3, NO2
964ccccccccccccccccccccccccccccccc
965
966c     Input atmospheric column
967
968      indexint=34
[1266]969      do i=1,nlayer
970         auxcolinp(nlayer-i+1) = h2o2colx(i) + o2colx(i) + o3colx(i) +
[705]971     $        no2colx(i)
972      end do
973
974c     Interpolation
975
976      do i=1,nz2
977         auxi = nz2-i+1
978         !O2 tabulated coefficient
979         auxjo2(i) = jabsifotsintpar(auxi,2,indexint)
980         !H2O2 tabulated coefficient
981         auxjh2o2(i) = jabsifotsintpar(auxi,6,indexint)
982         !O3 tabulated coefficient
983         auxjo3(i) = jabsifotsintpar(auxi,7,indexint)         
984         !Tabulated column
985         auxcoltab(i) = c34(nz2-i+1)
986      enddo
987      !Only if chemthermod.ge.2
988      if(chemthermod.ge.2) then
989         do i=1,nz2
990            !NO2 tabulated coefficient
991            auxjno2(i) = jabsifotsintpar(nz2-i+1,13,indexint)
992         enddo
993      endif
994      call interfast
[1266]995     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]996
[1266]997      do i=1,nlayer
[705]998         ind=auxind(i)
[1266]999         auxi = nlayer-i+1
[705]1000         !O2 interpolated coefficient
[1266]1001         jfotsout(indexint,2,auxi) = jfotsout(indexint,2,nlayer) *
[705]1002     $        (wm(i)*auxjo2(ind+1) + wp(i)*auxjo2(ind))
1003         !H2O2 interpolated coefficient
[1266]1004         jfotsout(indexint,6,auxi) = jfotsout(indexint,6,nlayer) *
[705]1005     $        (wm(i)*auxjh2o2(ind+1) + wp(i)*auxjh2o2(ind))
1006         !O3 interpolated coefficient
[1266]1007         jfotsout(indexint,7,auxi) = jfotsout(indexint,7,nlayer) *
[705]1008     $        (wm(i)*auxjo3(ind+1) + wp(i)*auxjo3(ind))
1009      enddo
1010      !Only if chemthermod.ge.2
1011      if(chemthermod.ge.2) then
[1266]1012         do i=1,nlayer
[705]1013            ind=auxind(i)
1014            !NO2 interpolated coefficient
[1266]1015            jfotsout(indexint,13,nlayer-i+1) =
1016     $           jfotsout(indexint,13,nlayer) *
[705]1017     $           (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind))
1018         enddo
1019      endif
1020
1021c     End of interval 34     
1022
1023
1024ccccccccccccccccccccccccccccccc
1025c     240.1-337.7nm (int 35)
1026c
1027c     Absorption by:
1028c     H2O2, O3, NO2
1029ccccccccccccccccccccccccccccccc
1030
1031c     Input atmospheric column
1032
1033      indexint=35
[1266]1034      do i=1,nlayer
1035         auxcolinp(nlayer-i+1) = h2o2colx(i) + o3colx(i) + no2colx(i)
[705]1036      end do
1037
1038c     Interpolation
1039
1040      do i=1,nz2
1041         auxi = nz2-i+1
1042         !H2O2 tabulated coefficient
1043         auxjh2o2(i) = jabsifotsintpar(auxi,6,indexint)
1044         !O3 tabulated coefficient
1045         auxjo3(i) = jabsifotsintpar(auxi,7,indexint)
1046         !Tabulated column
1047         auxcoltab(i) = c35(auxi)
1048      enddo
1049      !Only if chemthermod.ge.2
1050      if(chemthermod.ge.2) then
1051         do i=1,nz2
1052            !NO2 tabulated coefficient
1053            auxjno2(i) = jabsifotsintpar(nz2-i+1,13,indexint)
1054         enddo
1055      endif
1056      call interfast
[1266]1057     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]1058
[1266]1059      do i=1,nlayer
[705]1060         ind=auxind(i)
[1266]1061         auxi = nlayer-i+1
[705]1062         !H2O2 interpolated coefficient
[1266]1063         jfotsout(indexint,6,auxi) = jfotsout(indexint,6,nlayer) *
[705]1064     $        (wm(i)*auxjh2o2(ind+1) + wp(i)*auxjh2o2(ind))
1065         !O3 interpolated coefficient
[1266]1066         jfotsout(indexint,7,auxi) = jfotsout(indexint,7,nlayer) *
[705]1067     $        (wm(i)*auxjo3(ind+1) + wp(i)*auxjo3(ind))
1068      enddo
1069      if(chemthermod.ge.2) then
[1266]1070         do i=1,nlayer
[705]1071            ind=auxind(i)
1072            !NO2 interpolated coefficient
[1266]1073            jfotsout(indexint,13,nlayer-i+1) =
1074     $           jfotsout(indexint,13,nlayer) *
[705]1075     $           (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind))
1076         enddo
1077      endif
1078
1079c     End of interval 35
1080
1081ccccccccccccccccccccccccccccccc
1082c     337.8-800.0 nm (int 36)
1083c     
1084c     Absorption by:
1085c     O3, NO2
1086ccccccccccccccccccccccccccccccc
1087
1088c     Input atmospheric column
1089
1090      indexint=36
[1266]1091      do i=1,nlayer
1092         auxcolinp(nlayer-i+1) = o3colx(i) + no2colx(i)
[705]1093      end do
1094
1095c     Interpolation
1096
1097      do i=1,nz2
1098         auxi = nz2-i+1
1099         !O3 tabulated coefficient
1100         auxjo3(i) = jabsifotsintpar(auxi,7,indexint)         
1101         !Tabulated column
1102         auxcoltab(i) = c36(auxi)
1103      enddo
1104      !Only if chemthermod.ge.2
1105      if(chemthermod.ge.2) then
1106         do i=1,nz2
1107            !NO2 tabulated coefficient
1108            auxjno2(i) = jabsifotsintpar(nz2-i+1,13,indexint)
1109         enddo
1110      endif
1111      call interfast
[1266]1112     $     (wm,wp,auxind,auxcolinp,nlayer,auxcoltab,nz2,limdown,limup)
[2615]1113
[1266]1114      do i=1,nlayer
[705]1115         ind=auxind(i)
1116         !O3 interpolated coefficient
[1266]1117         jfotsout(indexint,7,nlayer-i+1) =
1118     $        jfotsout(indexint,7,nlayer) *
[705]1119     $        (wm(i)*auxjo3(ind+1) + wp(i)*auxjo3(ind))
1120      enddo
1121      !Only if chemthermod.ge.2
1122      if(chemthermod.ge.2) then
[1266]1123         do i=1,nlayer
[705]1124            ind=auxind(i)
1125            !NO2 interpolated coefficient
[1266]1126            jfotsout(indexint,13,nlayer-i+1) =
1127     $           jfotsout(indexint,13,nlayer) *
[705]1128     $           (wm(i)*auxjno2(ind+1) + wp(i)*auxjno2(ind))
1129         enddo
1130      endif
1131
1132c     End of interval 36
1133
1134c     End of interpolation to obtain photoabsorption rates
1135
[2429]1136c     Correction to the actual Sun-Mars distance
[705]1137
[2429]1138      jfotsout(:,:,:)=jfotsout(:,:,:)*(1.52/dist_sol)**2
[705]1139
[3464]1140      end subroutine jthermcalc_e107
1141     
1142      end module jthermcalc_e107_mod
[705]1143
1144
1145
1146
1147
1148
Note: See TracBrowser for help on using the repository browser.