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

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

LMDZ_MARS RV : Open_MP; Reading files in parallel for nlte parametrisation
(callnlte = .true.)

File size: 7.2 KB
RevLine 
[757]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
[1918]19      use datafile_mod, only: datadir
[2606]20      USE mod_phys_lmdz_para, ONLY: is_master
21      USE mod_phys_lmdz_transfert_para, ONLY: bcast
[1918]22     
[757]23      implicit none
24
25      include   'nlte_paramdef.h'
26      include   'nlte_commons.h'
27
28
29c***************
30
31c     local variables
32
33      integer   i, k, lun1, lun2
34      real*8          xx
35      character isotcode*2
36
37
38c     formats
39 132  format (i2)
40
41c**********************************************************************
42
43c     *** Groups old 1-d model subroutines SETTINGS and LeeESCTVCISO_dlvr11
44c     *** Both were called in old NLTEdlvr11_SETUP ***
45
46c     *** Old SETTINGS ***
47
48      lun1 = 1
49      lun2 = 2
50
[2606]51      if(is_master) then
52
[757]53      do k=1,nisot
54         write (isotcode,132) indexisot(k)
55         open (lun1,
[1918]56     $        file=trim(datadir)//'/NLTEDAT/enelow'
[757]57     $        //isotcode//'.dat',status='old')
58         open (lun2,
[1918]59     $        file=trim(datadir)//'/NLTEDAT/deltanu'
[757]60     $        //isotcode//'.dat',status='old')
61         read (lun1,*)
62         read (lun2,*)
63         read (lun1,*) (elow(k,i), i=1,nb)
64         read (lun2,*) (deltanu(k,i), i=1,nb)
65         close (lun1)
66         close (lun2)
67      end do
68
[2606]69      endif !      if(is_master)
70
71      call bcast(elow)
72      call bcast(deltanu)
73
[757]74      a1_010_000 = 1.3546d00
75      a2_010_000 = 1.3452d00
76      a3_010_000 = 1.1878d00
77      a4_010_000 = 1.2455d00
78      a1_020_010 = 4.35d0
79     
80
81c     *** Old LeeESCTVCISO_dlvr11 ***
82
[2606]83      if(is_master) then
84
[1918]85      open( 11, file=trim(datadir)//
[757]86     $     '/NLTEDAT/parametp_Tstar_IAA1204.dat' )
87      read (11, *)
88      do i=1,nztabul
89         read (11,*) lnpnbtab(i), tstar11tab(i),
90     $        tstar21tab(i), tstar31tab(i), tstar41tab(i)
91      enddo
92      close (11)
93
[1918]94      open( 12, file=trim(datadir)//
[757]95     $     '/NLTEDAT/parametp_VC_IAA1204.dat' )
96      read (12, *)
97      do i=1,nztabul
98         read (12,*) xx, vc210tab(i), vc310tab(i), vc410tab(i)
99      enddo
100      close (12)
[2606]101
102
103      endif !is_master
104
105      call bcast(lnpnbtab)
106      call bcast(tstar11tab)
107      call bcast(tstar21tab)
108      call bcast(tstar31tab)
109      call bcast(tstar41tab)
110      call bcast(xx)
111      call bcast(vc210tab)
112      call bcast(vc310tab)
113      call bcast(vc410tab)
114
[757]115      xx=xx
116
117
118
119      call LeeHISTOGRMS
120
121c     end subroutine
122
123      return
124      end
125
126
127
128c***********************************************************************
129      subroutine LeeHISTOGRMS
130c***********************************************************************
131
[1918]132      use datafile_mod, only: datadir
133
[757]134      implicit none
135
136      include   'nlte_paramdef.h'
137      include   'nlte_commons.h'
138
139
140c     local variables and constants
141      integer   ihist
142
143
144c***********************************************************************
145
146                                ! Banda fundamental
147                                !
[1918]148      hisfile = trim(datadir)//
[757]149     $     '/NLTEDAT/hid26-1.dat'
150      ihist = 1
151      call rhist_03 (ihist)
152
153
154                                ! First Hot bands       
155                                !     
[1918]156      hisfile = trim(datadir)//
[757]157     $     '/NLTEDAT/hid26-2.dat'
158      ihist = 2
159      call rhist_03 (ihist)
160
[1918]161      hisfile = trim(datadir)//
[757]162     $     '/NLTEDAT/hid26-3.dat'
163      ihist = 3
164      call rhist_03 (ihist)
165
[1918]166      hisfile = trim(datadir)//
[757]167     $     '/NLTEDAT/hid26-4.dat'
168      ihist = 4
169      call rhist_03 (ihist)
170
171
172
173
174      return
175      end
176
177
178c     *** Old GETK_dlvr11.f ***
179
180c***********************************************************************
181
182      subroutine GETK_dlvr11 (tt)
183
184c***********************************************************************
185
186      implicit none
187
188      include 'nlte_paramdef.h'
189      include 'nlte_commons.h'
190
191c     arguments
192      real              tt      ! i. temperature
193
194!     ! local variables:
195      real*8 k20x, k20xb, k20xc
196      real*8 k19xca,k19xcb,k19xcc
197      real*8 k19xba,k19xbb,k19xbc
198      real*8 k21x,k21xa,k21xb,k21xc
199      real*8 anu, factor , tdt
200      integer   i
201
202c***********************************************************************
203
204      tdt = dble(tt)
205
206                                !! k19 & k20
207
208      k20x = 3.d-12
209      k20xc = k20x * rf20
210      k20xb = 2.d0 * k20xc
211
212      k19xca = 4.2d-12 * exp( -2988.d0/tdt + 303930.d0/tdt**2.d0 )
213      if (tt.le.175.) k19xca = 3.3d-15
214      k19xcb = 2.1d-12 * exp( -2659.d0/tdt + 223052.d0/tdt**2.d0 )
215      if (tt.le.175.) k19xcb = 7.6d-16
216      k19xca = k19xca * rf19
217      k19xcb = k19xcb * rf19
218      k19xcc = k19xcb
219
220      factor = 2.5d0
221      k19xba = factor * k19xca
222      k19xbb = factor * k19xcb
223      k19xbc = factor * k19xcc
224
225      do i = 1, nisot
226
227         k19ba(i) = k19xba
228         k19ca(i) = k19xca
229         k19bb(i) = k19xbb
230         k19cb(i) = k19xcb
231         k19bc(i) = k19xbc
232         k19cc(i) = k19xcc
233
234         k20b(i) = k20xb
235         k20c(i) = k20xc
236
237         anu = dble( nu(i,2)-nu(i,1) )
238
239         k19bap(i) = k19ba(i) * 2.d0 * exp( -ee*anu/tdt )
240         k19bbp(i) = k19bb(i) * 2.d0 * exp( -ee*anu/tdt )
241         k19bcp(i) = k19bc(i) * 2.d0 * exp( -ee*anu/tdt )
242
243         k20bp(i) = k20b(i)*4.d0/2.d0 * exp( -ee/tdt * anu )
244
245         anu = dble( nu(i,1) )
246
247         k19cap(i) = k19ca(i) * 2.d0 * exp( -ee*anu/tdt )
248         k19cbp(i) = k19cb(i) * 2.d0 * exp( -ee*anu/tdt )
249         k19ccp(i) = k19cc(i) * 2.d0 * exp( -ee*anu/tdt )
250
251         k20cp(i) = k20c(i)*2.d0/1.d0 * exp( -ee/tdt * anu )
252
253      end do
254
255
256                                !! k21 &  k23k21c &  k24k21c & k34k21c
257
258      k21x = 2.49d-11
259      k21xb = k21x
260      k21xa = 3.d0/2.d0 * k21xb
261      k21xc = k21xb / 2.d0
262
263      k21xa = k21xa * rf21a
264      k21xb = k21xb * rf21b
265      k21xc = k21xc * rf21c
266
267      do i = 1, nisot
268         k21b(i) = k21xb
269         k21c(i) = k21xc
270         k21bp(i) = k21b(i) *
271     @        exp( -ee/tdt* dble(nu(i,2)-nu(i,1)-nu(1,1)) )
272         k21cp(i) = k21c(i) *
273     @        exp( -ee/tdt * dble(nu(i,1)-nu(1,1)) )
274      end do
275
276      k23k21c = k21xc
277      k24k21c = k21xc
278      k34k21c = k21xc
279      k23k21cp = k23k21c*2.d0/2.d0 *
280     @     exp( -ee/tdt* dble(nu(2,1)-nu(3,1)) )
281      k24k21cp = k24k21c*2.d0/2.d0 *
282     @     exp( -ee/tdt* dble(nu(2,1)-nu(4,1)) )
283      k34k21cp = k34k21c*2.d0/2.d0 *
284     @     exp( -ee/tdt* dble(nu(3,1)-nu(4,1)) )
285
286
287                                !! k33
288
289      k33c = k21x * rf33bc
290      do i=2,nisot
291         k33cp(i) = k33c *
292     @        exp( -ee/tdt * dble(nu(1,2)-nu(1,1)-nu(i,1)) )
293      end do
294
295
296      return
297      end
298
299
300
301
302
303
304
Note: See TracBrowser for help on using the repository browser.