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

Last change on this file since 1242 was 1047, checked in by emillour, 11 years ago

Mars GCM:

  • IMPORTANT CHANGE: Removed all reference/use of ngridmx (dimphys.h) in routines (necessary prerequisite to using parallel dynamics); in most cases this just means adding 'ngrid' as routine argument, and making local saved variables allocatable (and allocated at first call). In the process, had to convert many *.h files to equivalent modules: yomaer.h => yomaer_h.F90 , surfdat.h => surfdat_h.F90 , comsaison.h => comsaison_h.F90 , yomlw.h => yomlw_h.F90 , comdiurn.h => comdiurn_h.F90 , dimradmars.h => dimradmars_mod.F90 , comgeomfi.h => comgeomfi_h.F90, comsoil.h => comsoil_h.F90 , slope.h => slope_mod.F90
  • Also updated EOF routines, everything is now in eofdump_mod.F90
  • Removed unused routine lectfux.F (in dyn3d)

EM

File size: 11.8 KB
Line 
1c**************************************************************************
2c
3      subroutine nltecool(ngrid,nlayer,nq,pplay,pt,pq,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
27
28c       jul 2011 fgg   Modified to allow variable O
29c     
30c***************************************************************************
31
32      use tracer_mod, only: igcm_co2, igcm_co, igcm_o, igcm_n2, mmol
33      use conc_mod, only: mmean
34      implicit none
35
36#include "nltedata.h" ! (Equivalent to the reading of the "nlte_escape.dat" file)
37#include "dimensions.h"
38#include "dimphys.h"
39#include "chimiedata.h"
40!#include "conc.h" !Added to have "dynamic composition" in the scheme
41!#include "tracer.h" !"
42#include "callkeys.h"
43
44c Input and output variables
45c
46      integer     ngrid                           ! no. of horiz. gridpoints
47      integer     nlayer                          ! no. of atmospheric layers
48      integer     nq                              ! no. of tracers
49      real        pplay(ngrid,nlayer)             ! input pressure grid
50      real        pt(ngrid,nlayer)                ! input temperatures
51      real        pq(ngrid,nlayer,nq)                ! input mmrs
52      real        dtnlte(ngrid,nlayer)            ! output temp. tendencies
53
54c
55c Standard atmosphere variables
56c
57      real        nt                              ! number density [cm-3]
58      real        co2(nlayer)                     !  "   of CO2
59      real        o3p(nlayer)                     !  "   of atomic oxygen
60      real        n2co(nlayer)                    !  "  of N2 + CO
61      real        pyy(nlayer)                     ! auxiliary pressure grid
62
63c
64c Vectors and indexes for the tabulation of escape functions and VMR
65c
66c                 np                              ! # data points in tabulation
67c                 pnb(np)                         ! Pressure in tabulation
68c                 ef1(np)                         ! Esc.funct.#1, tabulated
69c                 ef2(np)                         ! Esc.funct.#2, tabulated
70c                 co2vmr(np)                      ! CO2 VMR tabulated
71c                 o3pvmr(np)                      ! CO2 VMR tabulated
72c                 n2covmr(np)                     ! N2+CO VMR tabulated
73      real        escf1(nlayer)                   ! Esc.funct.#1, interpolated
74      real        escf2(nlayer)                   ! Esc.funct.#2, interpolated
75
76
77c
78c Local Constants
79c
80      real       nu1, nu2                         ! freq. of energy levels
81      real       imr1, imr2                       ! isotopic abundances
82      real       hplanck, gamma, vlight           ! physical constants
83      real       ee
84      real       rfvt                             ! collisional rate
85      real       rfvto3p                          !     "
86      real       rfvv                             !     "
87
88c
89c Local variables for the main loop
90c
91      real       n1, n2, co2t                     ! ground populations
92      real       l1, p1, p12                      ! prod & losses
93      real       l2, p2, p21
94      real       tt                               ! dummy variable
95      real       c1, c2                           ! molecular constants
96      real       ae1, ae2                         ! einstein spontaneous emission
97      real       a1, a2, a12, a21
98      real       pl1, pl2
99      real       el1, el2
100      real       hr1, hr2                         ! heating rate due to each band
101      real       hr(nlayer)                       ! total heating rate
102
103c
104c Indexes
105c
106      integer    i
107      integer    j,ii
108
109c
110c Rate coefficients
111c
112      real       k19xca, k19xcb
113      real       k19cap1, k19cap2
114      real       k19cbp1, k19cbp2
115      real       d19c, d19cp1, d19cp2
116      real       k20xc, k20cp1, k20cp2
117      real       k21xc, k21cp2
118
119      logical    firstcall
120      data       firstcall/.true./
121      save       firstcall,ef1,ef2,co2vmr,n2covmr,o3pvmr,pnb
122
123c
124c Data
125c
126      data       nu1, nu2, hplanck, gamma, vlight, ee/
127     1     667.38, 662.3734, 6.6261e-27, 1.191e-5, 3.e10, 1.438769/
128
129c*************************************************************************
130c       PROGRAM  STARTS
131c*************************************************************************
132
133      imr1 = 0.987
134      imr2 = 0.00408 + 0.0112
135      rfvt = 0.1
136      rfvto3p = 1.0
137      rfvv = 0.1
138
139      if(firstcall) then
140
141         do i=1,np
142            pnb(i)=1.0e-4*exp(pnb(i)) ! p into Pa
143         end do
144
145         firstcall = .false.
146
147      endif
148
149c
150c MAIN LOOP, for each gridpoint and altitude:
151c
152      do j=1,ngrid  ! loop over grid points
153c
154c set up local pressure grid
155c
156         do ii=1,nlayer
157            pyy(ii)=pplay(j,ii)
158         enddo
159!
160! Interpolate escape functions and VMR to the desired grid
161!
162         call interp1(escf2,pyy,nlayer,ef2,pnb,np)
163         call interp1(escf1,pyy,nlayer,ef1,pnb,np)
164         if(nltemodel.eq.0) then
165            call interp3(co2,o3p,n2co,pyy,nlayer,
166     &           co2vmr,o3pvmr,n2covmr,pnb,np)
167         endif
168       
169         do i=1,nlayer  ! loop over layers
170C
171C test if p lies outside range (p > 3.5 Pa)
172C changed to 1 Pa since transition will always be higher than this
173C
174            if(pyy(i) .gt. 1.0 .or. pyy(i) .lt. 4.0e-6) then
175               hr(i)=0.0
176               dtnlte(j,i)=0.0
177            else
178c
179c           if(pt(j,i).lt.1.0)print*,pt(j,i)
180               nt = pyy(i)/(1.381e-17*pt(j,i)) ! nt in cm-3
181               !Dynamic composition
182               if(nltemodel.eq.1) then
183                  co2(i)=pq(j,i,igcm_co2)*mmean(j,i)/mmol(igcm_co2)
184                  o3p(i)=pq(j,i,igcm_o)*mmean(j,i)/mmol(igcm_o)
185                  n2co(i)=pq(j,i,igcm_co)*mmean(j,i)/mmol(igcm_co) +
186     $                 pq(j,i,igcm_n2)*mmean(j,i)/mmol(igcm_n2)
187               endif
188
189               !Mixing ratio to density
190               co2(i)=co2(i)*nt                ! CO2 density in cm-3
191               o3p(i)=o3p(i)*nt                ! O3p density in cm-3
192               n2co(i)=n2co(i)*nt              ! N2+CO in cm-3
193c molecular populations
194               n1 = co2(i) * imr1
195               n2 = co2(i) * imr2
196               co2t = n1 + n2
197
198c intermediate collisional rates
199               tt = pt(j,i)*pt(j,i)
200
201               if (pt(j,i).le.175.0) then
202                  k19xca = 3.3e-15
203                  k19xcb = 7.6e-16
204               else
205                  k19xca = 4.2e-12 * exp( -2988.0/pt(j,i) + 303930.0/tt)
206                  k19xcb = 2.1e-12 * exp( -2659.0/pt(j,i) + 223052.0/tt)
207               endif
208               k19xca = k19xca * rfvt
209               k19xcb = k19xcb * rfvt
210               k19cap1 = k19xca * 2.0 * exp( -ee*nu1/pt(j,i) )
211               k19cap2 = k19xca * 2.0 * exp( -ee*nu2/pt(j,i) )
212               k19cbp1 = k19xcb * 2.0 * exp( -ee*nu1/pt(j,i) )
213               k19cbp2 = k19xcb * 2.0 * exp( -ee*nu2/pt(j,i) )
214               d19c = k19xca*co2t + k19xcb*n2co(i)
215               d19cp1 = k19cap1*co2t + k19cbp1*n2co(i)
216               d19cp2 = k19cap2*co2t + k19cbp2*n2co(i)
217                                !
218               k20xc = 3.e-12 * rfvto3p
219               k20cp1 = k20xc * 2.0 * exp( -ee/pt(j,i) * nu1 )
220               k20cp2 = k20xc * 2.0 * exp( -ee/pt(j,i) * nu2 )
221                                !
222               k21xc = 2.49e-11 * 0.5 * rfvv
223               k21cp2 = k21xc * exp( - ee/pt(j,i) * (nu2-nu1) )
224                                !
225               l1 = d19c + k20xc*o3p(i) + k21cp2*n2
226               p1 = ( d19cp1 + k20cp1*o3p(i) ) * n1
227               p12 = k21xc*n1
228                                !
229               l2 = d19c + k20xc*o3p(i) + k21xc*n1
230               p2 = ( d19cp2 + k20cp2*o3p(i) ) * n2
231               p21 = k21cp2*n2
232
233c radiative rates
234               ae1 = 1.3546 * 1.66 / 4.0 * escf1(i)
235               ae2 = ( 1.3452 + 1.1878 ) * 1.66 / 4.0 * escf2(i)
236               l1 = l1 + ae1
237               l2 = l2 + ae2
238
239c solving the system
240               c1 = gamma*nu1**3. * 0.5
241               c2 = gamma*nu2**3. * 0.5
242               a1 = c1 * p1 / (n1*l1)
243               a2 = c2 * p2 / (n2*l2)
244               a12 = (nu1/nu2)**3. * n2/n1 * p12/l1
245               a21 = (nu2/nu1)**3. * n1/n2 * p21/l2
246               el2 = (a2 + a21 * a1 ) / ( 1.0 - a21 * a12 )
247               el1 = a1 + a12 * el2
248               pl1 = el1 * n1 / c1
249               pl2 = el2 * n2 / c2
250
251c  heating rate
252               hr1 = - hplanck*vlight * nu1 * ae1 * pl1
253               hr2 = - hplanck*vlight * nu2 * ae2 * pl2
254               hr(i) = hr1 + hr2
255               dtnlte(j,i)=0.1*hr(i)*pt(j,i)/(4.4*pyy(i)) ! dtnlte in K s-1
256c              write(7,25)pxx(i),hr1,hr2,hr(i),qt
257c  25         format(' ',1p5e12.4)
258
259            endif
260
261         enddo  ! end loop over layers
262      enddo     ! end loop over grid points
263c     close(7)
264c
265        return
266        end
267
268c***********************************************************************
269
270      subroutine interp1(escout,p,nlayer,escin,pin,nl)
271C
272C subroutine to perform linear interpolation in pressure from 1D profile
273C escin(nl) sampled on pressure grid pin(nl) to profile
274C escout(nlayer) on pressure grid p(nlayer).
275C
276      real escout(nlayer),p(nlayer)
277      real escin(nl),pin(nl),wm,wp
278      integer nl,nlayer,n1,n,nm,np
279      do n1=1,nlayer
280         if(p(n1) .gt. 3.5 .or. p(n1) .lt. 4.0e-6) then
281            escout(n1) = 0.0
282         else
283            do n = 1,nl-1
284               if (p(n1).le.pin(n).and.p(n1).ge.pin(n+1)) then
285                  nm=n
286                  np=n+1
287                  wm=abs(pin(np)-p(n1))/(pin(nm)-pin(np))
288                  wp=1.0 - wm
289               endif
290            enddo
291            escout(n1) = escin(nm)*wm + escin(np)*wp
292         endif
293      enddo
294      return
295      end
296
297c***********************************************************************
298
299      subroutine interp3(esco1,esco2,esco3,p,nlayer,
300     1     esci1,esci2,esci3,pin,nl)
301C
302C subroutine to perform 3 simultaneous linear interpolations in pressure from
303C 1D profiles esci1-3(nl) sampled on pressure grid pin(nl) to 1D profiles
304C esco1-3(nlayer) on pressure grid p(ngrid,nlayer).
305C
306      real esco1(nlayer),esco2(nlayer),esco3(nlayer),p(nlayer)
307      real esci1(nl),    esci2(nl),    esci3(nl),    pin(nl),wm,wp
308      integer nl,nlayer,n1,n,nm,np
309      do n1=1,nlayer
310         if (p(n1).gt. 3.5 .or. p(n1) .lt. 4.0e-6) then
311            esco1(n1)=0.0
312            esco2(n1)=0.0
313            esco3(n1)=0.0
314         else
315            do n = 1,nl-1
316               if (p(n1).le.pin(n).and.p(n1).ge.pin(n+1)) then
317                  nm=n
318                  np=n+1
319                  wm=abs(pin(np)-p(n1))/(pin(nm)-pin(np))
320                  wp=1.0 - wm
321               endif
322            enddo
323            esco1(n1) = esci1(nm)*wm + esci1(np)*wp
324            esco2(n1) = esci2(nm)*wm + esci2(np)*wp
325            esco3(n1) = esci3(nm)*wm + esci3(np)*wp
326         endif
327      enddo
328      return
329      end
Note: See TracBrowser for help on using the repository browser.