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

Last change on this file was 3466, checked in by emillour, 6 weeks ago

Mars PCM:
More tidying in aeronomars:

  • remove unused "inv.F" and remove "dtridgl.F" which is not used here and is a duplicate of the "dtridgl" routine in phymars/swr_toon.F
  • turn aeronomars routines to modules, for those which aren't in modules yet.

EM

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