source: trunk/LMDZ.VENUS/libf/phyvenus/nltecool.F @ 3461

Last change on this file since 3461 was 1530, checked in by emillour, 9 years ago

Venus and Titan GCMs:
Updates in the physics to keep up with updates in LMDZ5 (up to
LMDZ5 trunk, rev 2350) concerning dynamics/physics separation:

  • Adapted makelmdz and makelmdz_fcm script to stop if trying to compile 1d model or newstart or start2archive in parallel.
  • got rid of references to "dimensions.h" in physics. Within physics packages, use nbp_lon (=iim), nbp_lat (=jjmp1) and nbp_lev (=llm) from module mod_grid_phy_lmdz (in phy_common) instead. Only partially done for Titan, because of many hard-coded commons; a necessary first step will be to clean these up (using modules).

EM

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