source: trunk/LMDZ.MARS/libf/phymars/nltecool.F @ 308

Last change on this file since 308 was 38, checked in by emillour, 14 years ago

Ajout du modè Martien (mon LMDZ.MARS.BETA, du 28/01/2011) dans le rértoire mars, pour pouvoir suivre plus facilement les modifs.
EM

File size: 10.9 KB
Line 
1c**************************************************************************
2c
3      subroutine nltecool(ngrid,nlayer,pplay,pt,dtnlte)
4c
5c This code was designed as a delivery for the "Martian Environment Models"
6c project ( ESA contract 11369/95/nl/jg CCN2 )
7c Computes non-LTE heating rates from CO2 emission at 15 um
8c in the Martian upper atmosphere.
9c Uses a simplified model consisting of two excited levels with two
10c emission bands, one of them stronger than the other, which correspond
11c to the behaviours of the 626 fundamental band and the isotopic fund.bands.
12c It uses a cool-to-space approximation with tabulated escape functions.
13c These escape functions have been precomputed for the strong and weak bands,
14c and are given as a function of pressure in separate files.
15c The output values are the heating rates (actually, cooling, since they
16c are always negative) for the two bands, i.e., the total cooling is the
17c sum of them.
18c Miguel A. Lopez Valverde
19c Instituto de Astrofisica de Andalucia (CSIC), Granada, Spain
20c
21c Version 1b.  See description above.  22-March-2000.
22c Adapted as a subroutine for use in GCM -- PLR/SRL 6/2000
23c Version 1c.  Inclusion of VMR in the tabulation of escape functions.
24c              Table contains now only 1 input file -- Miguel 11/Jul/2000
25c Version 1d  data contained in original input file "nlte_escape.dat"
26c now stored in include file "nltedata.h" Y.Wanherdrick 09/2000
27c     
28c***************************************************************************
29
30      implicit none
31
32#include "nltedata.h" ! (Equivalent to the reading of the "nlte_escape.dat" file)
33 
34c Input and output variables
35c
36      integer     ngrid                           ! no. of horiz. gridpoints
37      integer     nlayer                          ! no. of atmospheric layers
38      real        pplay(ngrid,nlayer)             ! input pressure grid
39      real        pt(ngrid,nlayer)                ! input temperatures
40      real        dtnlte(ngrid,nlayer)            ! output temp. tendencies
41
42c
43c Standard atmosphere variables
44c
45      real        nt                              ! number density [cm-3]
46      real        co2(nlayer)                     !  "   of CO2
47      real        o3p(nlayer)                     !  "   of atomic oxygen
48      real        n2co(nlayer)                    !  "  of N2 + CO
49      real        pyy(nlayer)                     ! auxiliary pressure grid
50
51c
52c Vectors and indexes for the tabulation of escape functions and VMR
53c
54c                 np                              ! # data points in tabulation
55c                 pnb(np)                         ! Pressure in tabulation
56c                 ef1(np)                         ! Esc.funct.#1, tabulated
57c                 ef2(np)                         ! Esc.funct.#2, tabulated
58c                 co2vmr(np)                      ! CO2 VMR tabulated
59c                 o3pvmr(np)                      ! CO2 VMR tabulated
60c                 n2covmr(np)                     ! N2+CO VMR tabulated
61      real        escf1(nlayer)                   ! Esc.funct.#1, interpolated
62      real        escf2(nlayer)                   ! Esc.funct.#2, interpolated
63
64
65c
66c Local Constants
67c
68      real       nu1, nu2                         ! freq. of energy levels
69      real       imr1, imr2                       ! isotopic abundances
70      real       hplanck, gamma, vlight           ! physical constants
71      real       ee
72      real       rfvt                             ! collisional rate
73      real       rfvto3p                          !     "
74      real       rfvv                             !     "
75
76c
77c Local variables for the main loop
78c
79      real       n1, n2, co2t                     ! ground populations
80      real       l1, p1, p12                      ! prod & losses
81      real       l2, p2, p21
82      real       tt                               ! dummy variable
83      real       c1, c2                           ! molecular constants
84      real       ae1, ae2                         ! einstein spontaneous emission
85      real       a1, a2, a12, a21
86      real       pl1, pl2
87      real       el1, el2
88      real       hr1, hr2                         ! heating rate due to each band
89      real       hr(nlayer)                       ! total heating rate
90
91c
92c Indexes
93c
94      integer    i
95      integer    j,ii
96
97c
98c Rate coefficients
99c
100      real       k19xca, k19xcb
101      real       k19cap1, k19cap2
102      real       k19cbp1, k19cbp2
103      real       d19c, d19cp1, d19cp2
104      real       k20xc, k20cp1, k20cp2
105      real       k21xc, k21cp2
106
107      logical    firstcall
108      data       firstcall/.true./
109      save       firstcall,ef1,ef2,co2vmr,n2covmr,o3pvmr,pnb
110
111c
112c Data
113c
114      data       nu1, nu2, hplanck, gamma, vlight, ee/
115     1     667.38, 662.3734, 6.6261e-27, 1.191e-5, 3.e10, 1.438769/
116
117c*************************************************************************
118c       PROGRAM  STARTS
119c*************************************************************************
120
121      imr1 = 0.987
122      imr2 = 0.00408 + 0.0112
123      rfvt = 0.1
124      rfvto3p = 1.0
125      rfvv = 0.1
126
127      if(firstcall) then
128
129         do i=1,np
130            pnb(i)=1.0e-4*exp(pnb(i)) ! p into Pa
131         end do
132
133         firstcall = .false.
134
135      endif
136
137c
138c MAIN LOOP, for each gridpoint and altitude:
139c
140      do j=1,ngrid  ! loop over grid points
141c
142c set up local pressure grid
143c
144         do ii=1,nlayer
145            pyy(ii)=pplay(j,ii)
146         enddo
147!
148! Interpolate escape functions and VMR to the desired grid
149!
150         call interp1(escf2,pyy,nlayer,ef2,pnb,np)
151         call interp1(escf1,pyy,nlayer,ef1,pnb,np)
152         call interp3(co2,o3p,n2co,pyy,nlayer,
153     &        co2vmr,o3pvmr,n2covmr,pnb,np)
154       
155         do i=1,nlayer  ! loop over layers
156C
157C test if p lies outside range (p > 3.5 Pa)
158C changed to 1 Pa since transition will always be higher than this
159C
160            if(pyy(i) .gt. 1.0 .or. pyy(i) .lt. 4.0e-6) then
161               hr(i)=0.0
162               dtnlte(j,i)=0.0
163            else
164c
165c           if(pt(j,i).lt.1.0)print*,pt(j,i)
166               nt = pyy(i)/(1.381e-17*pt(j,i)) ! nt in cm-3
167               co2(i)=co2(i)*nt                ! CO2 density in cm-3
168               o3p(i)=o3p(i)*nt                ! O3p density in cm-3
169               n2co(i)=n2co(i)*nt              ! N2+CO in cm-3
170c molecular populations
171               n1 = co2(i) * imr1
172               n2 = co2(i) * imr2
173               co2t = n1 + n2
174
175c intermediate collisional rates
176               tt = pt(j,i)*pt(j,i)
177
178               if (pt(j,i).le.175.0) then
179                  k19xca = 3.3e-15
180                  k19xcb = 7.6e-16
181               else
182                  k19xca = 4.2e-12 * exp( -2988.0/pt(j,i) + 303930.0/tt)
183                  k19xcb = 2.1e-12 * exp( -2659.0/pt(j,i) + 223052.0/tt)
184               endif
185               k19xca = k19xca * rfvt
186               k19xcb = k19xcb * rfvt
187               k19cap1 = k19xca * 2.0 * exp( -ee*nu1/pt(j,i) )
188               k19cap2 = k19xca * 2.0 * exp( -ee*nu2/pt(j,i) )
189               k19cbp1 = k19xcb * 2.0 * exp( -ee*nu1/pt(j,i) )
190               k19cbp2 = k19xcb * 2.0 * exp( -ee*nu2/pt(j,i) )
191               d19c = k19xca*co2t + k19xcb*n2co(i)
192               d19cp1 = k19cap1*co2t + k19cbp1*n2co(i)
193               d19cp2 = k19cap2*co2t + k19cbp2*n2co(i)
194                                !
195               k20xc = 3.e-12 * rfvto3p
196               k20cp1 = k20xc * 2.0 * exp( -ee/pt(j,i) * nu1 )
197               k20cp2 = k20xc * 2.0 * exp( -ee/pt(j,i) * nu2 )
198                                !
199               k21xc = 2.49e-11 * 0.5 * rfvv
200               k21cp2 = k21xc * exp( - ee/pt(j,i) * (nu2-nu1) )
201                                !
202               l1 = d19c + k20xc*o3p(i) + k21cp2*n2
203               p1 = ( d19cp1 + k20cp1*o3p(i) ) * n1
204               p12 = k21xc*n1
205                                !
206               l2 = d19c + k20xc*o3p(i) + k21xc*n1
207               p2 = ( d19cp2 + k20cp2*o3p(i) ) * n2
208               p21 = k21cp2*n2
209
210c radiative rates
211               ae1 = 1.3546 * 1.66 / 4.0 * escf1(i)
212               ae2 = ( 1.3452 + 1.1878 ) * 1.66 / 4.0 * escf2(i)
213               l1 = l1 + ae1
214               l2 = l2 + ae2
215
216c solving the system
217               c1 = gamma*nu1**3. * 0.5
218               c2 = gamma*nu2**3. * 0.5
219               a1 = c1 * p1 / (n1*l1)
220               a2 = c2 * p2 / (n2*l2)
221               a12 = (nu1/nu2)**3. * n2/n1 * p12/l1
222               a21 = (nu2/nu1)**3. * n1/n2 * p21/l2
223               el2 = (a2 + a21 * a1 ) / ( 1.0 - a21 * a12 )
224               el1 = a1 + a12 * el2
225               pl1 = el1 * n1 / c1
226               pl2 = el2 * n2 / c2
227
228c  heating rate
229               hr1 = - hplanck*vlight * nu1 * ae1 * pl1
230               hr2 = - hplanck*vlight * nu2 * ae2 * pl2
231               hr(i) = hr1 + hr2
232               dtnlte(j,i)=0.1*hr(i)*pt(j,i)/(4.4*pyy(i)) ! dtnlte in K s-1
233c              write(7,25)pxx(i),hr1,hr2,hr(i),qt
234c  25         format(' ',1p5e12.4)
235
236            endif
237
238         enddo  ! end loop over layers
239      enddo     ! end loop over grid points
240c     close(7)
241c
242        return
243        end
244
245c***********************************************************************
246
247      subroutine interp1(escout,p,nlayer,escin,pin,nl)
248C
249C subroutine to perform linear interpolation in pressure from 1D profile
250C escin(nl) sampled on pressure grid pin(nl) to profile
251C escout(nlayer) on pressure grid p(nlayer).
252C
253      real escout(nlayer),p(nlayer)
254      real escin(nl),pin(nl),wm,wp
255      integer nl,nlayer,n1,n,nm,np
256      do n1=1,nlayer
257         if(p(n1) .gt. 3.5 .or. p(n1) .lt. 4.0e-6) then
258            escout(n1) = 0.0
259         else
260            do n = 1,nl-1
261               if (p(n1).le.pin(n).and.p(n1).ge.pin(n+1)) then
262                  nm=n
263                  np=n+1
264                  wm=abs(pin(np)-p(n1))/(pin(nm)-pin(np))
265                  wp=1.0 - wm
266               endif
267            enddo
268            escout(n1) = escin(nm)*wm + escin(np)*wp
269         endif
270      enddo
271      return
272      end
273
274c***********************************************************************
275
276      subroutine interp3(esco1,esco2,esco3,p,nlayer,
277     1     esci1,esci2,esci3,pin,nl)
278C
279C subroutine to perform 3 simultaneous linear interpolations in pressure from
280C 1D profiles esci1-3(nl) sampled on pressure grid pin(nl) to 1D profiles
281C esco1-3(nlayer) on pressure grid p(ngrid,nlayer).
282C
283      real esco1(nlayer),esco2(nlayer),esco3(nlayer),p(nlayer)
284      real esci1(nl),    esci2(nl),    esci3(nl),    pin(nl),wm,wp
285      integer nl,nlayer,n1,n,nm,np
286      do n1=1,nlayer
287         if (p(n1).gt. 3.5 .or. p(n1) .lt. 4.0e-6) then
288            esco1(n1)=0.0
289            esco2(n1)=0.0
290            esco3(n1)=0.0
291         else
292            do n = 1,nl-1
293               if (p(n1).le.pin(n).and.p(n1).ge.pin(n+1)) then
294                  nm=n
295                  np=n+1
296                  wm=abs(pin(np)-p(n1))/(pin(nm)-pin(np))
297                  wp=1.0 - wm
298               endif
299            enddo
300            esco1(n1) = esci1(nm)*wm + esci1(np)*wp
301            esco2(n1) = esci2(nm)*wm + esci2(np)*wp
302            esco3(n1) = esci3(nm)*wm + esci3(np)*wp
303         endif
304      enddo
305      return
306      end
Note: See TracBrowser for help on using the repository browser.