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

Last change on this file since 2800 was 2586, checked in by romain.vande, 3 years ago

LMDZ_MARS Implementation of Open_MP in the physic.
Works with radiative transfert

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