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

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