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

Last change on this file since 4022 was 4008, checked in by aslmd, 3 weeks ago

MESOSCALE: use precompiling flags to hide instructions related to parallel computations (we consider physics as being like a 1D model without any attached dynamical core when we compile).

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