source: trunk/LMDZ.MARS/libf/aeronomars/inichim_newstart.F @ 290

Last change on this file since 290 was 164, checked in by emillour, 13 years ago

Mars GCM:

Updates and corrections (to enable compiling/running in debug mode with

ifort)

  • removed option "-free-form" from makegcm_ifort and set mod_loc_dir="." so that module files (produced in local directory by ifort) are moved to LIBO
  • updated initracer.F, physdem1.F, physiq.F, inichim_newstart.F to avoid referencing out-of-bound array indexes (even if unused)
  • cosmetic updates on inwrite.F, datareadnc.F
  • updated newstart.F to initialize and use 'datadir' when looking for files
  • corrected bug on interpolation of sub-surface temperatures in lect_start_archive.F

EM

File size: 19.3 KB
Line 
1      subroutine inichim_newstart(pq,ps,sig,nq,latfi,lonfi,airefi,
2     $ thermo,qsurf)
3
4      implicit none
5
6c=======================================================================
7c
8c   subject:
9c   --------
10c
11c  Initialization of the composition for use in the building of a new start file
12c  used by program newstart.F
13c  also used by program testphys1d.F
14c
15c   Author:   Sebastien Lebonnois (08/11/2002)
16c   -------
17c Revised 07/2003 by Monica Angelats-i-Coll to use input files
18c Modified 10/2008 Identify tracers by their names Ehouarn Millour
19c
20c   Arguments:
21c   ----------
22c
23c    pq(iip1,jjp1,llm,nqmx)  Advected fields, ie chemical species here
24c    sig = sigma                vertical coordinate (interface of layers)
25c    nq: number of tracers to initialize (used to evaluate if only
26c        chemistry species are to be initialized or if water vapour
27c        and water ice should also be initialized.
28c        NB: number of chemistry tracers is ncomp (set in chimiedata.h)
29c=======================================================================
30
31c    Declarations :
32c    --------------
33
34#include "dimensions.h"
35#include "dimphys.h"
36#include "paramet.h"
37#include "chimiedata.h"
38#include "tracer.h"
39#include "comcstfi.h"
40#include "comdiurn.h"
41#include "callkeys.h"
42#include "temps.h"
43#include "datafile.h"
44
45c Arguments :
46c -----------
47
48      real    ps(iip1,jjp1)
49      real    pq(iip1,jjp1,llm,nqmx)     ! Advected fields, ie chemical species
50      real    sig(llm+1)                 ! vertical coordinate (interface of layers)
51      integer nq            !  =nqmx
52                            ! or =nqmx-1 if H2O kept
53                            ! or =nqmx-2 if H2O and ice kept
54      REAL latfi(ngridmx),lonfi(ngridmx),airefi(ngridmx)
55      integer thermo ! flag for thermosphere initialisation only
56      real :: qsurf(ngridmx,nqmx) ! surface tracers
57
58c Local variables :
59c -----------------
60
61      integer  iq,i,j,l,n, nspecini,inil,iqmax,iiq
62      INTEGER aa(nqmx)
63      REAL qtot(iip1,jjp1,llm)
64      real densconc,zgcm,zfile(252),vmrint,nt, mmean
65      real nxf(252),zfilemin(252),sigfile(252)
66      real vmrf(252,14),nf(252,14)
67      real tfile(252),zzfile(252)
68      real ni(14),niold(14)
69      real mmolf(14)    ! molecular mass in amu (species in files)
70      data mmolf/44.,40.,28.,32.,28.,16.,2.,   ! majors
71     &           1.,17.,33.,18.,34.,16.,48./   ! minors
72      character*3 tmp ! temporary variable
73      integer ierr,lnblnk
74      external lnblnk
75
76      logical :: oldnames ! =.true. if old tracer naming convention (q01,...)
77      character(len=20) :: txt ! to store some text
78      integer :: nqchem_start,count
79      integer :: nqchem(nqmx) ! local chemistry tracer index array
80      integer :: nbqchem ! total number of chemical species (water included)
81
82c-----------------------------------------------------------------------
83c    Input/Output
84c    ------------
85! read in 'callphys.def'
86       call inichim_readcallphys()
87
88      ! check if tracers have 'old' names
89      oldnames=.false.
90      count=0
91      do iq=1,nqmx
92        txt=" "
93        write(txt,'(a1,i2.2)') 'q',iq
94        if (txt.eq.noms(iq)) then
95          count=count+1
96        endif
97      enddo ! of do iq=1,nqmx
98
99      if (count.eq.nqmx) then
100        write(*,*) "inichim_newstart: tracers seem to follow old ",
101     &             "naming convention (q01,q02,...)"
102        write(*,*) "   => will work for now ... "
103        write(*,*) "      but you should rename them"
104        oldnames=.true.
105      endif
106
107      ! copy/set tracer names
108      if (oldnames) then ! old names (derived from iq & option values)
109        ! 1. dust:
110        if (dustbin.ne.0) then ! transported dust
111          do iq=1,dustbin
112            txt=" "
113            write(txt,'(a4,i2.2)') 'dust',iq
114            noms(iq)=txt
115            mmol(iq)=100.
116          enddo
117          ! special case if doubleq
118          if (doubleq) then
119            if (dustbin.ne.2) then
120              write(*,*) 'initracer: error expected dustbin=2'
121            else
122!              noms(1)='dust_mass'   ! dust mass mixing ratio
123!              noms(2)='dust_number' ! dust number mixing ratio
124! same as above, but avoid explict possible out-of-bounds on array
125               noms(1)='dust_mass'   ! dust mass mixing ratio
126               do iq=2,2
127               noms(iq)='dust_number' ! dust number mixing ratio
128               enddo
129            endif
130          endif
131        endif
132        ! 2. water & ice
133        if (water) then
134          noms(nqmx)='h2o_vap'
135          mmol(nqmx)=18.
136!            noms(nqmx-1)='h2o_ice'
137!            mmol(nqmx-1)=18.
138          !"loop" to avoid potential out-of-bounds in array
139          do iq=nqmx-1,nqmx-1
140            noms(iq)='h2o_ice'
141            mmol(iq)=18.
142          enddo
143        endif
144        ! 3. Chemistry
145        if (photochem .or. callthermos) then
146          if (doubleq) then
147            nqchem_start=3
148          else
149            nqchem_start=dustbin+1
150          end if
151        endif ! of if (photochem .or. callthermos)
152        do iq=nqchem_start,nqchem_start+ncomp-1
153          noms(iq)=nomchem(iq-nqchem_start+1)
154          mmol(iq)=mmolchem(iq-nqchem_start+1)
155        enddo
156        ! 4. Other tracers ????
157        if ((dustbin.eq.0).and.(.not.water)) then
158          noms(1)='co2'
159          mmol(1)=44
160          if (nqmx.eq.2) then
161            noms(nqmx)='Ar_N2'
162            mmol(nqmx)=30
163          endif
164        endif
165      else
166        ! noms(:) already contain tracer names
167        ! mmol(:) array is initialized later (see below)
168      endif ! of if (oldnames)
169
170      ! special modification when using 'old' tracers:
171      ! qsurf(nqmx) was h2o ice whereas q(nqmx) was water vapour
172      ! and (if iceparty) q(nqmx-1) was null whereas q(nqmx-1) was water ice
173      if (oldnames.and.water) then
174        write(*,*)'inichim_newstart: moving surface water ice to index '
175     &            ,nqmx-1
176        ! "loop" to avoid potential out-of-bounds on arrays
177        do iq=nqmx-1,nqmx-1
178        qsurf(1:ngridmx,iq)=qsurf(1:ngridmx,iq+1)
179        enddo
180        qsurf(1:ngridmx,nqmx)=0
181      endif
182
183c Dimension test:
184c ---------------
185
186       if (water) then
187!         if ((nqchem_min+ncomp+1).ne.nqmx) then
188         if ((dustbin+ncomp+2).ne.nqmx) then
189!            print*,'********* Dimension problem! ********'
190!            print*,"nqchem_min+ncomp+1).ne.nqmx"
191            print*,'inichim_newstart: tracer dimension problem:'
192            print*,"(dustbin+ncomp+2).ne.nqmx"
193            print*,"ncomp: ",ncomp
194!            print*,"nqchem_min: ",nqchem_min
195            print*,"nqmx: ",nqmx
196            print*,'Change ncomp in chimiedata.h'
197            STOP
198         endif
199       else
200!         if ((nqchem_min+ncomp).ne.nqmx) then
201          if ((dustbin+ncomp+1).ne.nqmx) then
202!            print*,'********* Dimension problem! ********'
203!            print*,"nqchem_min+ncomp).ne.nqmx"
204            print*,'initracer: tracer dimension problem:'
205            print*,"(dustbin+ncomp+1).ne.nqmx"
206            print*,"ncomp: ",ncomp
207!            print*,"nqchem_min: ",nqchem_min
208            print*,"nqmx: ",nqmx
209            print*,'Change ncomp in chimiedata.h'
210            STOP
211         endif
212       endif
213
214c noms and mmol vectors:
215c ----------------------
216!
217!       if (iceparty) then
218!          do iq=nqchem_min,nqmx-2
219!            noms(iq) = nomchem(iq-nqchem_min+1)
220!            mmol(iq) = mmolchem(iq-nqchem_min+1)
221!          enddo
222!          noms(nqmx-1) = 'ice'
223!          mmol(nqmx-1) = 18.
224!          noms(nqmx)   = 'h2o'
225!          mmol(nqmx)   = 18.
226!       else
227!          do iq=nqchem_min,nqmx-1
228!            noms(iq) = nomchem(iq-nqchem_min+1)
229!            mmol(iq) = mmolchem(iq-nqchem_min+1)
230!          enddo
231!          noms(nqmx)   = 'h2o'
232!          mmol(nqmx)   = 18.
233!       endif
234!       if (nqchem_min.gt.1) then
235!          do iq=1,nqchem_min-1
236!            noms(iq) = 'dust'
237!            mmol(iq) = 100.
238!          enddo
239!       endif
240
241
242! Identify tracers by their names: (and set corresponding values of mmol)
243      ! 0. initialize tracer indexes to zero:
244      do iq=1,nqmx
245        igcm_dustbin(iq)=0
246      enddo
247      igcm_dust_mass=0
248      igcm_dust_number=0
249      igcm_h2o_vap=0
250      igcm_h2o_ice=0
251      igcm_co2=0
252      igcm_co=0
253      igcm_o=0
254      igcm_o1d=0
255      igcm_o2=0
256      igcm_o3=0
257      igcm_h=0
258      igcm_h2=0
259      igcm_oh=0
260      igcm_ho2=0
261      igcm_h2o2=0
262      igcm_n2=0
263      igcm_ar=0
264      igcm_ar_n2=0
265
266      ! 1. find dust tracers
267      count=0
268      if (dustbin.gt.0) then
269        do iq=1,nqmx
270          txt=" "
271          write(txt,'(a4,i2.2)')'dust',count+1
272          if (noms(iq).eq.txt) then
273            count=count+1
274            igcm_dustbin(count)=iq
275            mmol(iq)=100.
276          endif
277        enddo !do iq=1,nqmx
278      endif ! of if (dustbin.gt.0)
279      if (doubleq) then
280        do iq=1,nqmx
281          if (noms(iq).eq."dust_mass") then
282            igcm_dust_mass=iq
283            count=count+1
284          endif
285          if (noms(iq).eq."dust_number") then
286            igcm_dust_number=iq
287            count=count+1
288          endif
289        enddo
290      endif ! of if (doubleq)
291      ! 2. find chemistry and water tracers
292      nbqchem=0
293      do iq=1,nqmx
294        if (noms(iq).eq."co2") then
295          igcm_co2=iq
296          mmol(igcm_co2)=44.
297          count=count+1
298          nbqchem=nbqchem+1
299        endif
300        if (noms(iq).eq."co") then
301          igcm_co=iq
302          mmol(igcm_co)=28.
303          count=count+1
304          nbqchem=nbqchem+1
305        endif
306        if (noms(iq).eq."o") then
307          igcm_o=iq
308          mmol(igcm_o)=16.
309          count=count+1
310          nbqchem=nbqchem+1
311        endif
312        if (noms(iq).eq."o1d") then
313          igcm_o1d=iq
314          mmol(igcm_o1d)=16.
315          count=count+1
316          nbqchem=nbqchem+1
317        endif
318        if (noms(iq).eq."o2") then
319          igcm_o2=iq
320          mmol(igcm_o2)=32.
321          count=count+1
322          nbqchem=nbqchem+1
323        endif
324        if (noms(iq).eq."o3") then
325          igcm_o3=iq
326          mmol(igcm_o3)=48.
327          count=count+1
328          nbqchem=nbqchem+1
329        endif
330        if (noms(iq).eq."h") then
331          igcm_h=iq
332          mmol(igcm_h)=1.
333          count=count+1
334          nbqchem=nbqchem+1
335        endif
336        if (noms(iq).eq."h2") then
337          igcm_h2=iq
338          mmol(igcm_h2)=2.
339          count=count+1
340          nbqchem=nbqchem+1
341        endif
342        if (noms(iq).eq."oh") then
343          igcm_oh=iq
344          mmol(igcm_oh)=17.
345          count=count+1
346          nbqchem=nbqchem+1
347        endif
348        if (noms(iq).eq."ho2") then
349          igcm_ho2=iq
350          mmol(igcm_ho2)=33.
351          count=count+1
352          nbqchem=nbqchem+1
353        endif
354        if (noms(iq).eq."h2o2") then
355          igcm_h2o2=iq
356          mmol(igcm_h2o2)=34.
357          count=count+1
358          nbqchem=nbqchem+1
359        endif
360        if (noms(iq).eq."n2") then
361          igcm_n2=iq
362          mmol(igcm_n2)=28.
363          count=count+1
364          nbqchem=nbqchem+1
365        endif
366        if (noms(iq).eq."ar") then
367          igcm_ar=iq
368          mmol(igcm_ar)=40.
369          count=count+1
370          nbqchem=nbqchem+1
371        endif
372        if (noms(iq).eq."h2o_vap") then
373          igcm_h2o_vap=iq
374          mmol(igcm_h2o_vap)=18.
375          count=count+1
376          nbqchem=nbqchem+1
377        endif
378        if (noms(iq).eq."h2o_ice") then
379          igcm_h2o_ice=iq
380          mmol(igcm_h2o_ice)=18.
381          count=count+1
382          nbqchem=nbqchem+1
383        endif
384        ! Other stuff: e.g. for simulations using co2 + neutral gaz
385        if (noms(iq).eq."Ar_N2") then
386          igcm_ar_n2=iq
387          mmol(igcm_ar_n2)=30.
388          count=count+1
389        endif
390      enddo ! of do iq=1,nqmx
391     
392      ! check that we identified all tracers:
393      if (count.ne.nqmx) then
394        write(*,*) "inichim_newstart: found only ",count," tracers"
395        write(*,*) "               expected ",nqmx
396        do iq=1,count
397          write(*,*)'      ',iq,' ',trim(noms(iq))
398        enddo
399        stop
400      else
401       write(*,*)"inichim_newstart: found all expected tracers, namely:"
402        do iq=1,nqmx
403          write(*,*)'      ',iq,' ',trim(noms(iq))
404        enddo
405      endif
406
407      ! build local chemistry tracer index array nqchem(:)
408      ! as in the old days, water vapour is last and water ice,
409      ! if present, is just before water vapour
410      do iq=1,15 ! loop so as to avoid out-of-bounds on array
411      if (iq==1) nqchem(i)=igcm_co2
412      if (iq==2) nqchem(i)=igcm_co
413      if (iq==3) nqchem(i)=igcm_o
414      if (iq==4) nqchem(i)=igcm_o1d
415      if (iq==5) nqchem(i)=igcm_o2
416      if (iq==6) nqchem(i)=igcm_o3
417      if (iq==7) nqchem(i)=igcm_h
418      if (iq==8) nqchem(i)=igcm_h2
419      if (iq==9) nqchem(i)=igcm_oh
420      if (iq==10) nqchem(i)=igcm_ho2
421      if (iq==11) nqchem(i)=igcm_h2o2
422      if (iq==12) nqchem(i)=igcm_n2
423      if (iq==13) nqchem(i)=igcm_ar
424      if (iq==14) nqchem(i)=igcm_h2o_ice
425      if (iq==15) nqchem(i)=igcm_h2o_vap
426      enddo
427
428! Load in chemistry data for initialization:
429c----------------------------------------------------------------------
430c Order of Major species in file:
431c
432c    1=CO2
433c    2=AR
434c    3=N2 
435c    4=O2 
436c    5=CO 
437c    6=O   
438c    7=H2
439c
440c Order of Minor species in files:
441c
442c    1=H 
443c    2=OH
444c    3=HO2
445c    4=H2O
446c    5=H2O2
447c    6=O1D
448c    7=O3
449c
450c----------------------------------------------------------------------
451c Major species:
452        aa(igcm_co2) = 1
453        aa(igcm_ar)  = 2
454        aa(igcm_n2)  = 3
455        aa(igcm_o2)  = 4
456        aa(igcm_co)  = 5
457        aa(igcm_o)   = 6
458        aa(igcm_h2)  = 7
459c Minor species:
460        aa(igcm_h)       = 8
461        aa(igcm_oh)      = 9
462        aa(igcm_ho2)     = 10
463        aa(igcm_h2o_vap) = 11
464        aa(igcm_h2o2)    = 12
465        aa(igcm_o1d)     = 13
466        aa(igcm_o3)      = 14
467c----------------------------------------------------------------------
468      open(210, iostat=ierr,file=
469     & datafile(1:lnblnk(datafile))//'/atmosfera_LMD_may.dat')
470      if (ierr.ne.0) then
471        write(*,*)'Error : cannot open file atmosfera_LMD_may.dat '
472        write(*,*)'(in aeronomars/inichim.F)'
473        write(*,*)'It should be in :', datafile(1:lnblnk(datafile)),'/'
474        write(*,*)'1) You can change this directory address in '
475        write(*,*)'   file phymars/datafile.h'
476        write(*,*)'2) If necessary atmosfera_LMD_may.dat (and others)'
477        write(*,*)'   can be obtained online on:'
478        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
479         STOP
480      endif
481      open(220, iostat=ierr,file=
482     & datafile(1:lnblnk(datafile))//'/atmosfera_LMD_min.dat')
483      if (ierr.ne.0) then
484        write(*,*)'Error : cannot open file atmosfera_LMD_min.dat '
485        write(*,*)'(in aeronomars/inichim.F)'
486        write(*,*)'It should be in :', datafile(1:lnblnk(datafile)),'/'
487        write(*,*)'1) You can change this directory address in '
488        write(*,*)'   file phymars/datafile.h'
489        write(*,*)'2) If necessary atmosfera_LMD_min.dat (and others)'
490        write(*,*)'   can be obtained online on:'
491        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
492         STOP
493      endif
494
495      read(210,*) tmp
496      read(220,*) tmp
497        do l=1,252
498          read(210,112) tfile(l),zfile(l),nxf(l),(vmrf(l,n),n=1,7)
499          zfile(l)=zfile(l)*100.              !pressure (Pa)
500          sigfile(l)=zfile(l)/zfile(1)       
501        enddo
502        do l=1,252
503          read(220,113)zfilemin(l),(vmrf(l,n),n=8,14)
504          zfilemin(l)=zfilemin(l)*1000.       !height (m)
505          do n=1,14
506            nf(l,n)=vmrf(l,n)*nxf(l)
507          enddo
508        enddo
509        close(210)
510        close(220)
511
512c flag thermo set to 1 or 0 in newstart
513c inil=33 for initialisation of thermosphere only       
514c inil=1 for initialisation of all layers       
515        if (thermo.eq.1) then
516        inil=33
517        else
518        inil=1
519        endif
520
521! Initialization of tracers
522
523      do i=1,iip1
524       do j=1,jjp1
525        do l=inil,llm
526
527c          zgcm=sig(l)
528          zgcm=sig(l)*ps(i,j)
529          densconc=0.
530          nt=0.
531
532          do n=1,14
533            call intrplf(zgcm,vmrint,zfile,nf(1,n),252)
534c            call intrplf(zgcm,vmrint,sigfile,nf(1,n),252)
535            ni(n)=vmrint
536
537            densconc=densconc+ni(n)*mmolf(n)
538            nt=nt+ni(n)
539          enddo
540
541          mmean=densconc/nt                     ! in amu
542
543          if (nq.ne.nbqchem) then ! don't initialize water vapour
544            do n=1,nq-1
545            pq(i,j,l,nqchem(n))=
546     &                 ni(aa(nqchem(n)))/nt*mmol(nqchem(n))/mmean
547            if(i.eq.iip1) pq(i,j,l,nqchem(n))=pq(1,j,l,nqchem(n))
548            enddo
549            if (water .and. (nq.gt.nbqchem-2)) then
550              pq(i,j,l,nqchem(nq)) = 0.
551              if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
552            else
553              pq(i,j,l,nqchem(nq))=
554     &              ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean
555              if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
556            endif
557          endif ! of if (nq.ne.nbqchem)
558         
559          if (nq.eq.nbqchem) then ! also initialize water vapour
560            do n=1,nq-2
561              pq(i,j,l,nqchem(n))=
562     &                   ni(aa(nqchem(n)))/nt*mmol(nqchem(n))/mmean
563              if(i.eq.iip1) pq(i,j,l,nqchem(n))=pq(1,j,l,nqchem(n))
564            enddo
565            if (water) then
566              pq(i,j,l,nqchem(nq-1)) = 0.
567              if(i.eq.iip1) then
568                pq(i,j,l,nqchem(nq-1))=pq(1,j,l,nqchem(nq-1))
569              endif
570              pq(i,j,l,nqchem(nq))=
571     &                   ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean
572              if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
573            else
574              do n=nq-1,nq
575                pq(i,j,l,nqchem(nq))=
576     &                   ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean
577                if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
578              enddo
579            endif ! of if (water)
580          endif ! of if (nq.eq.nqmx)
581        enddo     !nlayer
582       enddo      !ngrid_j
583      enddo       !ngrid_i
584
585cccccccccccccccccccccccccccccc     
586c  make sure that sum of q = 1     
587c dominent species is = 1 - sum(all other species)     
588cccccccccccccccccccccccccccccc     
589!      iqmax=nqchem_min
590       iqmax=1
591     
592!      if ((nqmx-nqchem_min).gt.10) then
593      if (nbqchem.gt.10) then
594       do l=1,llm
595        do j=1,jjp1
596          do i=1,iip1
597!           do iq=nqchem_min,nqmx
598            do iiq=1,nbqchem
599              iq=nqchem(iiq)
600              if ( (pq(i,j,l,iq).gt.pq(i,j,l,iqmax)).and.
601     .           (noms(iq).ne."ice") )  then
602              iqmax=iq
603              endif
604            enddo
605            pq(i,j,l,iqmax)=1.
606            qtot(i,j,l)=0
607!           do iq=nqchem_min,nqmx
608            do iiq=1,nbqchem
609             iq=nqchem(iiq)
610             if ( (iq.ne.iqmax).and.
611     .           (noms(iq).ne."ice") ) then       
612               pq(i,j,l,iqmax)=pq(i,j,l,iqmax)-pq(i,j,l,iq)       
613             endif
614            enddo !iq
615!            do iq=nqchem_min,nqmx
616            do iiq=1,nbqchem
617             iq=nqchem(iiq)
618              if (noms(iq).ne."ice") then
619                qtot(i,j,l)=qtot(i,j,l)+pq(i,j,l,iq)
620              endif
621c            if (i.eq.1.and.j.eq.1.and.l.eq.1) write(*,*) 'qtot(i,j,l)',
622c     $    qtot(i,j,l)
623            enddo !iq
624            if (i.eq.1.and.j.eq.1) write(*,*) 'inichim_newstart: ',
625     $        'qtot(i,j,l)=',qtot(i,j,l)
626          enddo !i   
627         enddo !j   
628       enddo !l 
629      endif ! of if (nbqchem.gt.10)
630ccccccccccccccccccccccccccccccc
631
63266      format(i2,6(1x,e11.5))
633112     format(7x,f7.3,3x,e12.6,3x,e12.6,7(3x,e12.6))
634113     format(1x,f6.2,7(3x,e12.6))
635
636      end
Note: See TracBrowser for help on using the repository browser.