source: trunk/LMDZ.VENUS/libf/phyvenus/nlte_setup.F @ 1661

Last change on this file since 1661 was 1442, checked in by slebonnois, 9 years ago

SL: update of the Venus GCM, + corrections on routines used for newstart/start2archive for Titan and Venus, + some modifications on tools

File size: 6.9 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      implicit none
20
21#include "nlte_paramdef.h"
22#include "nlte_commons.h"
23
24
25c***************
26
27c     local variables
28
29      integer   i, k, lun1, lun2
30      real*8          xx
31      character isotcode*2
32      character (len=100) :: datafile="HIGHATM"
33
34c     formats
35 132  format (i2)
36
37c**********************************************************************
38
39c     *** Groups old 1-d model subroutines SETTINGS and LeeESCTVCISO_dlvr11
40c     *** Both were called in old NLTEdlvr11_SETUP ***
41
42c     *** Old SETTINGS ***
43
44      lun1 = 1
45      lun2 = 2
46
47      do k=1,nisot
48         write (isotcode,132) indexisot(k)
49         open (lun1,
50     $        file=trim(datafile)//'/enelow'
51     $        //isotcode//'.dat',status='old')
52         open (lun2,
53     $        file=trim(datafile)//'/deltanu'
54     $        //isotcode//'.dat',status='old')
55         read (lun1,*)
56         read (lun2,*)
57         read (lun1,*) (elow(k,i), i=1,nb)
58         read (lun2,*) (deltanu(k,i), i=1,nb)
59         close (lun1)
60         close (lun2)
61      end do
62
63      a1_010_000 = 1.3546d00
64      a2_010_000 = 1.3452d00
65      a3_010_000 = 1.1878d00
66      a4_010_000 = 1.2455d00
67      a1_020_010 = 4.35d0
68     
69
70c     *** Old LeeESCTVCISO_dlvr11 ***
71
72      open( 11, file=trim(datafile)//
73     $     '/parametp_Tstar_IAA1204.dat' )
74      read (11, *)
75      do i=1,nztabul
76         read (11,*) lnpnbtab(i), tstar11tab(i),
77     $        tstar21tab(i), tstar31tab(i), tstar41tab(i)
78      enddo
79      close (11)
80
81      open( 12, file=trim(datafile)//
82     $     '/parametp_VC_IAA1204.dat' )
83      read (12, *)
84      do i=1,nztabul
85         read (12,*) xx, vc210tab(i), vc310tab(i), vc410tab(i)
86      enddo
87      close (12)
88      xx=xx
89
90
91
92      call LeeHISTOGRMS
93
94c     end subroutine
95
96      return
97      end
98
99
100
101c***********************************************************************
102      subroutine LeeHISTOGRMS
103c***********************************************************************
104
105      implicit none
106
107      include   'nlte_paramdef.h'
108      include   'nlte_commons.h'
109
110
111c     local variables and constants
112      integer   ihist
113
114      character (len=100) :: datafile="HIGHATM"
115c***********************************************************************
116
117                                ! Banda fundamental
118                                !
119      hisfile = trim(datafile)//
120     $     '/hid26-1.dat'
121      ihist = 1
122      call rhist_03 (ihist)
123
124
125                                ! First Hot bands       
126                                !     
127      hisfile = trim(datafile)//
128     $     '/hid26-2.dat'
129      ihist = 2
130      call rhist_03 (ihist)
131
132      hisfile = trim(datafile)//
133     $     '/hid26-3.dat'
134      ihist = 3
135      call rhist_03 (ihist)
136
137      hisfile = trim(datafile)//
138     $     '/hid26-4.dat'
139      ihist = 4
140      call rhist_03 (ihist)
141
142
143
144
145      return
146      end
147
148
149c     *** Old GETK_dlvr11.f ***
150
151c***********************************************************************
152
153      subroutine GETK_dlvr11 (tt)
154
155c***********************************************************************
156
157      implicit none
158
159      include 'nlte_paramdef.h'
160      include 'nlte_commons.h'
161
162c     arguments
163      real              tt      ! i. temperature
164
165!     ! local variables:
166      real*8 k20x, k20xb, k20xc
167      real*8 k19xca,k19xcb,k19xcc
168      real*8 k19xba,k19xbb,k19xbc
169      real*8 k21x,k21xa,k21xb,k21xc
170      real*8 anu, factor , tdt
171      integer   i
172
173c***********************************************************************
174
175      tdt = dble(tt)
176
177                                !! k19 & k20
178
179      k20x = 3.d-12
180c  TEST GG: double the values of Kvv as recently found by Sharma et al.2014
181c      k20x = 6.d-12   
182c  TEST GG: use the minimum value of the experimental bracket's values [1-6]
183c      k20x = 1.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
234CCC TEST GG
235c      k21x =  2.49d-11*0.5
236C      k21x =  2.49d-11*2
237
238      k21xb = k21x
239      k21xa = 3.d0/2.d0 * k21xb
240      k21xc = k21xb / 2.d0
241
242      k21xa = k21xa * rf21a
243      k21xb = k21xb * rf21b
244      k21xc = k21xc * rf21c
245
246      do i = 1, nisot
247         k21b(i) = k21xb
248         k21c(i) = k21xc
249         k21bp(i) = k21b(i) *
250     @        exp( -ee/tdt* dble(nu(i,2)-nu(i,1)-nu(1,1)) )
251         k21cp(i) = k21c(i) *
252     @        exp( -ee/tdt * dble(nu(i,1)-nu(1,1)) )
253      end do
254
255      k23k21c = k21xc
256      k24k21c = k21xc
257      k34k21c = k21xc
258      k23k21cp = k23k21c*2.d0/2.d0 *
259     @     exp( -ee/tdt* dble(nu(2,1)-nu(3,1)) )
260      k24k21cp = k24k21c*2.d0/2.d0 *
261     @     exp( -ee/tdt* dble(nu(2,1)-nu(4,1)) )
262      k34k21cp = k34k21c*2.d0/2.d0 *
263     @     exp( -ee/tdt* dble(nu(3,1)-nu(4,1)) )
264
265
266                                !! k33
267
268      k33c = k21x * rf33bc
269      do i=2,nisot
270         k33cp(i) = k33c *
271     @        exp( -ee/tdt * dble(nu(1,2)-nu(1,1)-nu(i,1)) )
272      end do
273
274
275      return
276      end
277
278
279
280
281
282
283
Note: See TracBrowser for help on using the repository browser.