source: trunk/LMDZ.MARS/libf/aeronomars/param_read.F @ 1621

Last change on this file since 1621 was 1381, checked in by emillour, 10 years ago

Mars GCM:

EM

File size: 9.6 KB
Line 
1      subroutine param_read
2
3      use param_v4_h, only: jfotsout,crscabsi2,
4     .    c1_16,c17_24,c25_29,c30_31,c32,c33,c34,c35,c36,
5     .    co2crsc195,co2crsc295,t0,
6     .    jabsifotsintpar,ninter,nz2,
7     .    efdisco2,efdiso2,efdish2o,
8     .    efdish2o2,efdish2,efdiso3,
9     .    efdiso,efdisn,efdish,
10     .    efdisno,efdisn2,efdisno2,
11     .    efdisco,efionco2,efionn2,
12     .    efionco,efiono3p,efionn,
13     .    efionno,efionh,
14     .    fluxtop,ct1,ct2,p1,p2
15
16      implicit none
17
18
19c     common variables and constants
20      include 'datafile.h'
21 
22 
23c     local variables
24
25      integer    i,j,k,inter                          !indexes
26      integer ierr
27      real       nada
28 
29     
30c*************************+PROGRAM STARTS**************************
31
32c     Reads tabulated functions
33
34      !Tabulated column amount
35      open(210, status = 'old',
36c    $file=trim(datafile)//'/EUVDAT/coln.dat',iostat=ierr)
37     $file=trim(datafile)//'/EUVDAT/param_v5/coln.dat',iostat=ierr)
38
39      IF (ierr.NE.0) THEN
40       write(*,*)'cant find directory EUVDAT containing param_v5 subdir'
41       write(*,*)'(in aeronomars/param_read.F)'
42       write(*,*)'It should be in :', trim(datafile),'/'
43       write(*,*)'1) You can change this directory address in '
44       write(*,*)'   callphys.def with datadir=/path/to/dir'
45       write(*,*)'2) If necessary, EUVDAT (and other datafiles)'
46       write(*,*)'   can be obtained online on:'
47       write(*,*)' http://www.lmd.jussieu.fr/~lmdz/planets/mars/datadir'
48       STOP
49      ENDIF
50 
51      !Tabulated photoabsorption coefficients
52      open(220,file=trim(datafile)//'/EUVDAT/param_v5/j2_an.dat')
53      open(230,file=trim(datafile)//'/EUVDAT/param_v5/j3_an.dat')
54      open(240,file=trim(datafile)//'/EUVDAT/param_v5/j1_an.dat')
55      open(250,file=trim(datafile)//'/EUVDAT/param_v5/j2_bn.dat')
56      open(260,file=trim(datafile)//'/EUVDAT/param_v5/j2_cn.dat')
57      open(300,file=trim(datafile)//'/EUVDAT//param_v5/j2_dn.dat')
58      open(270,file=trim(datafile)//'/EUVDAT//param_v5/j1_bn.dat')
59      open(280,file=trim(datafile)//'/EUVDAT//param_v5/j1_cn.dat')
60      open(290,file=trim(datafile)//'/EUVDAT//param_v5/j1_dn.dat')
61      open(150,file=trim(datafile)//'/EUVDAT//param_v5/j4n.dat')
62      open(160,file=trim(datafile)//'/EUVDAT//param_v5/j5n.dat')
63      open(170,file=trim(datafile)//'/EUVDAT//param_v5/j6n.dat')
64      open(180,file=trim(datafile)//'/EUVDAT//param_v5/j7n.dat')
65      open(390,file=trim(datafile)//'/EUVDAT//param_v5/j8_an.dat')
66      open(400,file=trim(datafile)//'/EUVDAT//param_v5/j8_bn.dat')
67      open(410,file=trim(datafile)//'/EUVDAT//param_v5/j9n.dat')
68      open(420,file=trim(datafile)//'/EUVDAT//param_v5/j10_an.dat')
69      open(430,file=trim(datafile)//'/EUVDAT//param_v5/j10_bn.dat')
70      open(440,file=trim(datafile)//'/EUVDAT//param_v5/j10_cn.dat')
71      open(450,file=trim(datafile)//'/EUVDAT//param_v5/j11_an.dat')
72      open(460,file=trim(datafile)//'/EUVDAT//param_v5/j11_bn.dat')
73      open(470,file=trim(datafile)//'/EUVDAT//param_v5/j11_cn.dat')
74      open(480,file=trim(datafile)//'/EUVDAT//param_v5/j12n.dat')
75      open(490,file=trim(datafile)//'/EUVDAT//param_v5/j13_an.dat')
76      open(500,file=trim(datafile)//'/EUVDAT//param_v5/j13_bn.dat')
77      open(510,file=trim(datafile)//'/EUVDAT//param_v5/j13_cn.dat')
78
79     
80      do i=210,300,10
81         read(i,*)
82         read(i,*)
83      end do
84
85      do i=150,180,10
86         read(i,*)
87         read(i,*)
88      end do
89
90      do i=390,510,10
91         read(i,*)
92         read(i,*)
93      enddo
94
95      do i=nz2,1,-1
96         read(210,*) (c1_16(i,j),j=1,16),c17_24(i),c25_29(i),c30_31(i),
97     $        c32(i),c33(i),c34(i),c35(i),c36(i)
98      end do
99
100      do i=nz2,1,-1
101         read(220,*) (jabsifotsintpar(i,2,j),j=1,16)
102      end do
103     
104      do i=nz2,1,-1
105         read(230,*) (jabsifotsintpar(i,3,j),j=1,16)
106      end do
107
108      do i=nz2,1,-1
109         read(240,*) (jabsifotsintpar(i,1,j),j=1,16)
110      end do
111
112      do i=nz2,1,-1
113         read(250,*) (jabsifotsintpar(i,2,j),j=17,24)
114      end do
115
116
117      do i=nz2,1,-1
118         read(260,*) (jabsifotsintpar(i,2,j),j=25,31)
119      end do
120
121      do i=nz2,1,-1
122         read(270,*) (jabsifotsintpar(i,1,j),j=17,24)
123      end do
124
125      do i=nz2,1,-1
126         read(280,*) (jabsifotsintpar(i,1,j),j=25,31)
127      end do
128
129      do i=nz2,1,-1
130         read(290,*) jabsifotsintpar(i,1,32)
131      end do
132
133      do i=nz2,1,-1
134         read(300,*) (jabsifotsintpar(i,2,j),j=32,34)
135      end do
136
137      do i=nz2,1,-1
138         read(160,*) (jabsifotsintpar(i,5,j),j=1,15)
139      end do
140
141      do i=nz2,1,-1
142         read(150,*) (jabsifotsintpar(i,4,j),j=25,31)
143      end do
144
145      do i=nz2,1,-1
146         read(170,*) (jabsifotsintpar(i,6,j),j=25,35)
147      end do
148
149      do i=nz2,1,-1
150         read(180,*) (jabsifotsintpar(i,7,j),j=34,36)
151      end do
152
153      do i=nz2,1,-1
154         read(390,*) (jabsifotsintpar(i,8,j),j=2,16)
155      enddo
156
157      do i=nz2,1,-1
158         read(400,*) (jabsifotsintpar(i,8,j),j=17,24)
159      enddo
160
161      do i=nz2,1,-1
162         read(410,*) (jabsifotsintpar(i,9,j),j=1,16)
163      enddo
164
165      do i=nz2,1,-1
166         read(420,*) (jabsifotsintpar(i,10,j),j=2,16)
167      enddo
168
169      do i=nz2,1,-1
170         read(430,*) (jabsifotsintpar(i,10,j),j=17,24)
171      enddo
172
173      do i=nz2,1,-1
174         read(440,*) (jabsifotsintpar(i,10,j),j=25,32)
175      enddo
176
177      do i=nz2,1,-1
178         read(450,*) (jabsifotsintpar(i,11,j),j=2,16)
179      enddo
180
181      do i=nz2,1,-1
182         read(460,*) (jabsifotsintpar(i,11,j),j=17,24)
183      enddo
184
185      do i=nz2,1,-1
186         read(470,*) (jabsifotsintpar(i,11,j),j=25,29)
187      enddo
188     
189      do i=nz2,1,-1
190         read(480,*) (jabsifotsintpar(i,12,j),j=2,16)
191      enddo
192
193      do i=nz2,1,-1
194         read(490,*) (jabsifotsintpar(i,13,j),j=2,16)
195      enddo
196     
197      do i=nz2,1,-1
198         read(500,*) (jabsifotsintpar(i,13,j),j=17,24)
199      enddo
200     
201      do i=nz2,1,-1
202         read(510,*) (jabsifotsintpar(i,13,j),j=25,36)
203      enddo
204
205      do i=210,300,10
206         close(i)
207      end do
208
209      do i=150,180,10
210         close(i)
211      end do
212
213      do i=390,510,10
214         close(i)
215      enddo
216
217
218c     set t0
219
220      do i=1,nz2
221         t0(i)=195.
222      end do
223
224
225      do i=1,ninter
226         fluxtop(i)=1.
227      end do
228
229      !Parameters for the variation of the solar flux with 11 years cycle
230      open(100,file=trim(datafile)//'/EUVDAT/param_v5/varflujo.dat')
231      read(100,*)
232      do i=1,24
233         read(100,*) inter,ct1(i),p1(i),ct2(i),p2(i),nada
234      end do
235      close(100)
236
237c     dissociation and ionization efficiencies
238
239      do inter=1,ninter
240         efdisco2(inter)=0.
241         efdiso2(inter)=0.
242         efdish2(inter)=0.
243         efdish2o(inter)=0.
244         efdish2o2(inter)=0.
245         efdiso3(inter)=0.
246         efdisco(inter)=0.
247         efdisn2(inter)=0.
248         efdisno(inter)=0.
249         efdisno2(inter)=0.
250         efionco2(inter,1)=0.
251         efionco2(inter,2)=0.
252         efionco2(inter,3)=0.
253         efionco2(inter,4)=0.
254         efiono3p(inter)=0.
255         efionn2(inter,1)=0.
256         efionn2(inter,2)=0.
257         efionco(inter,1)=0.
258         efionco(inter,2)=0.
259         efionco(inter,3)=0.
260         efionn(inter)=0.
261         efionh(inter)=0.
262         efionno(inter)=0.
263      enddo
264
265
266c     CO2, O2, NO
267
268      open(120,file=trim(datafile)//'/EUVDAT/param_v5/efdis_inter.dat')
269      read(120,*)
270!      do i=1,21
271!         read(120,*)inter,efdisco2(inter),efdiso2(inter),efdisno(inter)
272      do inter=8,28
273         read(120,*)i,efdisco2(inter),efdiso2(inter),efdisno(inter)
274      enddo
275      do inter=29,ninter
276         efdisco2(inter)=1.
277         efdiso2(inter)=1.
278         efdisno(inter)=1.
279      enddo
280
281
282c     N2
283
284      efdisn2(15)=0.1
285      do inter=16,ninter
286         efdisn2(inter)=1.
287      enddo
288
289
290c     CO
291
292      efdisco(16)=0.5
293      do inter=17,ninter
294         efdisco(inter)=1.
295      enddo
296
297     
298c     O, N, H
299
300      do inter=1,ninter
301         efdiso(inter)=0.
302         efdisn(inter)=0.
303         efdish(inter)=0.
304      enddo
305
306
307c     H2O, H2O2, O3, NO2
308
309      do inter=25,31
310         efdish2o(inter)=1.
311      enddo
312      do inter=25,35
313         efdish2o2(inter)=1.
314      enddo
315      do inter=34,36
316         efdiso3(inter)=1.
317      enddo
318      do inter=27,36
319         efdisno2(inter)=1.
320      enddo
321      do inter=1,15
322         efdish2(inter)=1.
323      enddo
324         
325      !4 possible channels for CO2 ionization
326      do inter=14,16
327         efionco2(inter,1)=1.-efdisco2(inter)
328      enddo
329      efionco2(13,1)=0.805*(1.-efdisco2(13))
330      efionco2(13,2)=0.195*(1.-efdisco2(13))
331      do inter=11,12
332         efionco2(inter,3)=1.-efdisco2(inter)
333      enddo
334      efionco2(10,3)=0.9*(1.-efdisco2(10))
335      efionco2(10,4)=0.1*(1.-efdisco2(10))
336      do inter=2,9
337         efionco2(inter,4)=1.-efdisco2(inter)
338      enddo
339
340      !For O(3p) total ionization under 91.1 nm
341      do inter=1,16
342         efiono3p(inter)=1.d0
343      enddo
344
345      !2 channels for N2 ionization
346      do inter=9,15
347         efionn2(inter,1)=1.-efdisn2(inter)
348      enddo
349      do inter=2,8
350         efionn2(inter,2)=1.-efdisn2(inter)
351      enddo
352     
353      !3 channels for CO ionization
354      do inter=11,16
355         efionco(inter,1)=1.-efdisco(inter)
356      enddo
357      efionco(10,1)=0.87*(1.-efdisco(10))
358      efionco(10,2)=0.13*(1.-efdisco(10))
359      do inter=8,9
360         efionco(inter,2)=1.-efdisco(inter)
361      enddo
362      efionco(7,2)=0.1*(1.-efdisco(7))
363      efionco(7,3)=0.9*(1.-efdisco(7))
364      do inter=2,6
365         efionco(inter,3)=1.-efdisco(inter)
366      enddo
367
368      !Total ionization under 85 nm for N
369      do inter=1,16
370         efionn(inter)=1.
371      enddo
372
373      !NO
374      do inter=2,28
375         efionno(inter)=1.-efdisno(inter)
376      enddo
377
378      !Total ionization under 90 nm for H
379      do inter=3,16
380         efionh(inter)=1.
381      enddo
382
383
384      return
385
386
387      end
388
Note: See TracBrowser for help on using the repository browser.