source: trunk/LMDZ.MARS/libf/phymars/nlte_setup.F @ 2266

Last change on this file since 2266 was 1918, checked in by emillour, 7 years ago

Mars GCM:
Code cleanup:

  • remove "comorbit.h" since it is no longer used.
  • turn "datafile.h" into module datafile_mod.F90 (and rename variable "datafile" as "datadir" since it stores the path to the datafile directory).

EM

File size: 6.7 KB
Line 
1c***********************************************************************
2
3      subroutine nlte_setup
4
5c     malv    Oct 09          Adapt mz1d_onlyTCR_MUCHASveces.f to "V09"
6c     malv    Sep 07          Add LU deccomp & repetition option to test CPU
7c     malv    Jan 07          Add new vertical fine-grid for NLTE
8c     apr 06  malv            Read date,effuv from Driver. T fixed at zbott.
9c     2003    fgg             Double precission in UV, Photoq, Conduct & Diff
10c     oct 02  malv            V02: New scheme to allow for continuity eq.
11c     dec 01  malv            See changes/progress of the code in mz1d.actual
12c     nov 01  malv            adapt for parameterizations of tcr y shr
13c     nov 98  malv            add chemical & photochem. processes
14c     jul 98  malv            transic hiperb con zs fuera de la region
15c     equil hidrostatico. smoothing en cr y sh
16c     jan 98    malv            first version
17c***********************************************************************
18
19      use datafile_mod, only: datadir
20     
21      implicit none
22
23      include   'nlte_paramdef.h'
24      include   'nlte_commons.h'
25
26
27c***************
28
29c     local variables
30
31      integer   i, k, lun1, lun2
32      real*8          xx
33      character isotcode*2
34
35
36c     formats
37 132  format (i2)
38
39c**********************************************************************
40
41c     *** Groups old 1-d model subroutines SETTINGS and LeeESCTVCISO_dlvr11
42c     *** Both were called in old NLTEdlvr11_SETUP ***
43
44c     *** Old SETTINGS ***
45
46      lun1 = 1
47      lun2 = 2
48
49      do k=1,nisot
50         write (isotcode,132) indexisot(k)
51         open (lun1,
52     $        file=trim(datadir)//'/NLTEDAT/enelow'
53     $        //isotcode//'.dat',status='old')
54         open (lun2,
55     $        file=trim(datadir)//'/NLTEDAT/deltanu'
56     $        //isotcode//'.dat',status='old')
57         read (lun1,*)
58         read (lun2,*)
59         read (lun1,*) (elow(k,i), i=1,nb)
60         read (lun2,*) (deltanu(k,i), i=1,nb)
61         close (lun1)
62         close (lun2)
63      end do
64
65      a1_010_000 = 1.3546d00
66      a2_010_000 = 1.3452d00
67      a3_010_000 = 1.1878d00
68      a4_010_000 = 1.2455d00
69      a1_020_010 = 4.35d0
70     
71
72c     *** Old LeeESCTVCISO_dlvr11 ***
73
74      open( 11, file=trim(datadir)//
75     $     '/NLTEDAT/parametp_Tstar_IAA1204.dat' )
76      read (11, *)
77      do i=1,nztabul
78         read (11,*) lnpnbtab(i), tstar11tab(i),
79     $        tstar21tab(i), tstar31tab(i), tstar41tab(i)
80      enddo
81      close (11)
82
83      open( 12, file=trim(datadir)//
84     $     '/NLTEDAT/parametp_VC_IAA1204.dat' )
85      read (12, *)
86      do i=1,nztabul
87         read (12,*) xx, vc210tab(i), vc310tab(i), vc410tab(i)
88      enddo
89      close (12)
90      xx=xx
91
92
93
94      call LeeHISTOGRMS
95
96c     end subroutine
97
98      return
99      end
100
101
102
103c***********************************************************************
104      subroutine LeeHISTOGRMS
105c***********************************************************************
106
107      use datafile_mod, only: datadir
108
109      implicit none
110
111      include   'nlte_paramdef.h'
112      include   'nlte_commons.h'
113
114
115c     local variables and constants
116      integer   ihist
117
118
119c***********************************************************************
120
121                                ! Banda fundamental
122                                !
123      hisfile = trim(datadir)//
124     $     '/NLTEDAT/hid26-1.dat'
125      ihist = 1
126      call rhist_03 (ihist)
127
128
129                                ! First Hot bands       
130                                !     
131      hisfile = trim(datadir)//
132     $     '/NLTEDAT/hid26-2.dat'
133      ihist = 2
134      call rhist_03 (ihist)
135
136      hisfile = trim(datadir)//
137     $     '/NLTEDAT/hid26-3.dat'
138      ihist = 3
139      call rhist_03 (ihist)
140
141      hisfile = trim(datadir)//
142     $     '/NLTEDAT/hid26-4.dat'
143      ihist = 4
144      call rhist_03 (ihist)
145
146
147
148
149      return
150      end
151
152
153c     *** Old GETK_dlvr11.f ***
154
155c***********************************************************************
156
157      subroutine GETK_dlvr11 (tt)
158
159c***********************************************************************
160
161      implicit none
162
163      include 'nlte_paramdef.h'
164      include 'nlte_commons.h'
165
166c     arguments
167      real              tt      ! i. temperature
168
169!     ! local variables:
170      real*8 k20x, k20xb, k20xc
171      real*8 k19xca,k19xcb,k19xcc
172      real*8 k19xba,k19xbb,k19xbc
173      real*8 k21x,k21xa,k21xb,k21xc
174      real*8 anu, factor , tdt
175      integer   i
176
177c***********************************************************************
178
179      tdt = dble(tt)
180
181                                !! k19 & k20
182
183      k20x = 3.d-12
184      k20xc = k20x * rf20
185      k20xb = 2.d0 * k20xc
186
187      k19xca = 4.2d-12 * exp( -2988.d0/tdt + 303930.d0/tdt**2.d0 )
188      if (tt.le.175.) k19xca = 3.3d-15
189      k19xcb = 2.1d-12 * exp( -2659.d0/tdt + 223052.d0/tdt**2.d0 )
190      if (tt.le.175.) k19xcb = 7.6d-16
191      k19xca = k19xca * rf19
192      k19xcb = k19xcb * rf19
193      k19xcc = k19xcb
194
195      factor = 2.5d0
196      k19xba = factor * k19xca
197      k19xbb = factor * k19xcb
198      k19xbc = factor * k19xcc
199
200      do i = 1, nisot
201
202         k19ba(i) = k19xba
203         k19ca(i) = k19xca
204         k19bb(i) = k19xbb
205         k19cb(i) = k19xcb
206         k19bc(i) = k19xbc
207         k19cc(i) = k19xcc
208
209         k20b(i) = k20xb
210         k20c(i) = k20xc
211
212         anu = dble( nu(i,2)-nu(i,1) )
213
214         k19bap(i) = k19ba(i) * 2.d0 * exp( -ee*anu/tdt )
215         k19bbp(i) = k19bb(i) * 2.d0 * exp( -ee*anu/tdt )
216         k19bcp(i) = k19bc(i) * 2.d0 * exp( -ee*anu/tdt )
217
218         k20bp(i) = k20b(i)*4.d0/2.d0 * exp( -ee/tdt * anu )
219
220         anu = dble( nu(i,1) )
221
222         k19cap(i) = k19ca(i) * 2.d0 * exp( -ee*anu/tdt )
223         k19cbp(i) = k19cb(i) * 2.d0 * exp( -ee*anu/tdt )
224         k19ccp(i) = k19cc(i) * 2.d0 * exp( -ee*anu/tdt )
225
226         k20cp(i) = k20c(i)*2.d0/1.d0 * exp( -ee/tdt * anu )
227
228      end do
229
230
231                                !! k21 &  k23k21c &  k24k21c & k34k21c
232
233      k21x = 2.49d-11
234      k21xb = k21x
235      k21xa = 3.d0/2.d0 * k21xb
236      k21xc = k21xb / 2.d0
237
238      k21xa = k21xa * rf21a
239      k21xb = k21xb * rf21b
240      k21xc = k21xc * rf21c
241
242      do i = 1, nisot
243         k21b(i) = k21xb
244         k21c(i) = k21xc
245         k21bp(i) = k21b(i) *
246     @        exp( -ee/tdt* dble(nu(i,2)-nu(i,1)-nu(1,1)) )
247         k21cp(i) = k21c(i) *
248     @        exp( -ee/tdt * dble(nu(i,1)-nu(1,1)) )
249      end do
250
251      k23k21c = k21xc
252      k24k21c = k21xc
253      k34k21c = k21xc
254      k23k21cp = k23k21c*2.d0/2.d0 *
255     @     exp( -ee/tdt* dble(nu(2,1)-nu(3,1)) )
256      k24k21cp = k24k21c*2.d0/2.d0 *
257     @     exp( -ee/tdt* dble(nu(2,1)-nu(4,1)) )
258      k34k21cp = k34k21c*2.d0/2.d0 *
259     @     exp( -ee/tdt* dble(nu(3,1)-nu(4,1)) )
260
261
262                                !! k33
263
264      k33c = k21x * rf33bc
265      do i=2,nisot
266         k33cp(i) = k33c *
267     @        exp( -ee/tdt * dble(nu(1,2)-nu(1,1)-nu(i,1)) )
268      end do
269
270
271      return
272      end
273
274
275
276
277
278
279
Note: See TracBrowser for help on using the repository browser.