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

Last change on this file since 560 was 552, checked in by emillour, 13 years ago

Mars GCM:

bug fix in "aeronomars/moldiff_red.F90" (wrong mmol array size in included

subroutine) + adapted extreme limit test for H to altitudes of ~4000km
(compatble with 50 layer model).

bug fix in "nirco2abs.F" index of "co2" and "o" tracers were hard coded as

1 and 3!

updates from FGG of euvheat.F, callkeys.h and inifis.F to have the

"euveff" parameter read from callphys.def

EM

File size: 19.5 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   Modified 11/2011 Addition of methane (Franck Lefevre)
20c
21c   Arguments:
22c   ----------
23c
24c    pq(iip1,jjp1,llm,nqmx)  Advected fields, ie chemical species here
25c    sig = sigma                vertical coordinate (interface of layers)
26c    nq: number of tracers to initialize (used to evaluate if only
27c        chemistry species are to be initialized or if water vapour
28c        and water ice should also be initialized.
29c        NB: number of chemistry tracers is ncomp (set in chimiedata.h)
30c=======================================================================
31
32c    Declarations :
33c    --------------
34
35#include "dimensions.h"
36#include "dimphys.h"
37#include "paramet.h"
38#include "chimiedata.h"
39#include "tracer.h"
40#include "comcstfi.h"
41#include "comdiurn.h"
42#include "callkeys.h"
43#include "temps.h"
44#include "datafile.h"
45
46c Arguments :
47c -----------
48
49      real    ps(iip1,jjp1)
50      real    pq(iip1,jjp1,llm,nqmx)     ! Advected fields, ie chemical species
51      real    sig(llm+1)                 ! vertical coordinate (interface of layers)
52      integer nq            !  =nqmx
53                            ! or =nqmx-1 if H2O kept
54                            ! or =nqmx-2 if H2O and ice kept
55      REAL latfi(ngridmx),lonfi(ngridmx),airefi(ngridmx)
56      integer thermo ! flag for thermosphere initialisation only
57      real :: qsurf(ngridmx,nqmx) ! surface tracers
58
59c Local variables :
60c -----------------
61
62      integer  iq,i,j,l,n, nspecini,inil,iqmax,iiq
63      INTEGER aa(nqmx)
64      REAL qtot(iip1,jjp1,llm)
65      real densconc,zgcm,zfile(252),vmrint,nt, mmean
66      real nxf(252),zfilemin(252),sigfile(252)
67      real vmrf(252,14),nf(252,14)
68      real tfile(252),zzfile(252)
69      real ni(14),niold(14)
70      real mmolf(14)    ! molecular mass in amu (species in files)
71      data mmolf/44.,40.,28.,32.,28.,16.,2.,   ! majors
72     &           1.,17.,33.,18.,34.,16.,48./   ! minors
73      character*3 tmp ! temporary variable
74      integer ierr
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-2
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_ch4=0
263      igcm_n2=0
264      igcm_ar=0
265      igcm_ar_n2=0
266
267      ! 1. find dust tracers
268      count=0
269      if (dustbin.gt.0) then
270        do iq=1,nqmx
271          txt=" "
272          write(txt,'(a4,i2.2)')'dust',count+1
273          if (noms(iq).eq.txt) then
274            count=count+1
275            igcm_dustbin(count)=iq
276            mmol(iq)=100.
277          endif
278        enddo !do iq=1,nqmx
279      endif ! of if (dustbin.gt.0)
280      if (doubleq) then
281        do iq=1,nqmx
282          if (noms(iq).eq."dust_mass") then
283            igcm_dust_mass=iq
284            count=count+1
285          endif
286          if (noms(iq).eq."dust_number") then
287            igcm_dust_number=iq
288            count=count+1
289          endif
290        enddo
291      endif ! of if (doubleq)
292      ! 2. find chemistry and water tracers
293      nbqchem=0
294      do iq=1,nqmx
295        if (noms(iq).eq."co2") then
296          igcm_co2=iq
297          mmol(igcm_co2)=44.
298          count=count+1
299          nbqchem=nbqchem+1
300        endif
301        if (noms(iq).eq."co") then
302          igcm_co=iq
303          mmol(igcm_co)=28.
304          count=count+1
305          nbqchem=nbqchem+1
306        endif
307        if (noms(iq).eq."o") then
308          igcm_o=iq
309          mmol(igcm_o)=16.
310          count=count+1
311          nbqchem=nbqchem+1
312        endif
313        if (noms(iq).eq."o1d") then
314          igcm_o1d=iq
315          mmol(igcm_o1d)=16.
316          count=count+1
317          nbqchem=nbqchem+1
318        endif
319        if (noms(iq).eq."o2") then
320          igcm_o2=iq
321          mmol(igcm_o2)=32.
322          count=count+1
323          nbqchem=nbqchem+1
324        endif
325        if (noms(iq).eq."o3") then
326          igcm_o3=iq
327          mmol(igcm_o3)=48.
328          count=count+1
329          nbqchem=nbqchem+1
330        endif
331        if (noms(iq).eq."h") then
332          igcm_h=iq
333          mmol(igcm_h)=1.
334          count=count+1
335          nbqchem=nbqchem+1
336        endif
337        if (noms(iq).eq."h2") then
338          igcm_h2=iq
339          mmol(igcm_h2)=2.
340          count=count+1
341          nbqchem=nbqchem+1
342        endif
343        if (noms(iq).eq."oh") then
344          igcm_oh=iq
345          mmol(igcm_oh)=17.
346          count=count+1
347          nbqchem=nbqchem+1
348        endif
349        if (noms(iq).eq."ho2") then
350          igcm_ho2=iq
351          mmol(igcm_ho2)=33.
352          count=count+1
353          nbqchem=nbqchem+1
354        endif
355        if (noms(iq).eq."h2o2") then
356          igcm_h2o2=iq
357          mmol(igcm_h2o2)=34.
358          count=count+1
359          nbqchem=nbqchem+1
360        endif
361        if (noms(iq).eq."ch4") then
362          igcm_ch4=iq
363          mmol(igcm_ch4)=16.
364          count=count+1
365          nbqchem=nbqchem+1
366        endif
367        if (noms(iq).eq."n2") then
368          igcm_n2=iq
369          mmol(igcm_n2)=28.
370          count=count+1
371          nbqchem=nbqchem+1
372        endif
373        if (noms(iq).eq."ar") then
374          igcm_ar=iq
375          mmol(igcm_ar)=40.
376          count=count+1
377          nbqchem=nbqchem+1
378        endif
379        if (noms(iq).eq."h2o_vap") then
380          igcm_h2o_vap=iq
381          mmol(igcm_h2o_vap)=18.
382          count=count+1
383          nbqchem=nbqchem+1
384        endif
385        if (noms(iq).eq."h2o_ice") then
386          igcm_h2o_ice=iq
387          mmol(igcm_h2o_ice)=18.
388          count=count+1
389          nbqchem=nbqchem+1
390        endif
391        ! Other stuff: e.g. for simulations using co2 + neutral gaz
392        if (noms(iq).eq."Ar_N2") then
393          igcm_ar_n2=iq
394          mmol(igcm_ar_n2)=30.
395          count=count+1
396        endif
397      enddo ! of do iq=1,nqmx
398     
399      ! check that we identified all tracers:
400      if (count.ne.nqmx) then
401        write(*,*) "inichim_newstart: found only ",count," tracers"
402        write(*,*) "               expected ",nqmx
403        do iq=1,count
404          write(*,*)'      ',iq,' ',trim(noms(iq))
405        enddo
406        stop
407      else
408       write(*,*)"inichim_newstart: found all expected tracers, namely:"
409        do iq=1,nqmx
410          write(*,*)'      ',iq,' ',trim(noms(iq))
411        enddo
412      endif
413
414      ! build local chemistry tracer index array nqchem(:)
415      ! as in the old days, water vapour is last and water ice,
416      ! if present, is just before water vapour
417      do iq=1,16 ! loop so as to avoid out-of-bounds on array
418      if (iq==1) nqchem(i)=igcm_co2
419      if (iq==2) nqchem(i)=igcm_co
420      if (iq==3) nqchem(i)=igcm_o
421      if (iq==4) nqchem(i)=igcm_o1d
422      if (iq==5) nqchem(i)=igcm_o2
423      if (iq==6) nqchem(i)=igcm_o3
424      if (iq==7) nqchem(i)=igcm_h
425      if (iq==8) nqchem(i)=igcm_h2
426      if (iq==9) nqchem(i)=igcm_oh
427      if (iq==10) nqchem(i)=igcm_ho2
428      if (iq==11) nqchem(i)=igcm_h2o2
429      if (iq==12) nqchem(i)=igcm_ch4
430      if (iq==13) nqchem(i)=igcm_n2
431      if (iq==14) nqchem(i)=igcm_ar
432      if (iq==15) nqchem(i)=igcm_h2o_ice
433      if (iq==16) nqchem(i)=igcm_h2o_vap
434      enddo
435
436! Load in chemistry data for initialization:
437c----------------------------------------------------------------------
438c Order of Major species in file:
439c
440c    1=CO2
441c    2=AR
442c    3=N2 
443c    4=O2 
444c    5=CO 
445c    6=O   
446c    7=H2
447c
448c Order of Minor species in files:
449c
450c    1=H 
451c    2=OH
452c    3=HO2
453c    4=H2O
454c    5=H2O2
455c    6=O1D
456c    7=O3
457c
458c----------------------------------------------------------------------
459c Major species:
460        aa(igcm_co2) = 1
461        aa(igcm_ar)  = 2
462        aa(igcm_n2)  = 3
463        aa(igcm_o2)  = 4
464        aa(igcm_co)  = 5
465        aa(igcm_o)   = 6
466        aa(igcm_h2)  = 7
467c Minor species:
468        aa(igcm_h)       = 8
469        aa(igcm_oh)      = 9
470        aa(igcm_ho2)     = 10
471        aa(igcm_h2o_vap) = 11
472        aa(igcm_h2o2)    = 12
473        aa(igcm_o1d)     = 13
474        aa(igcm_o3)      = 14
475c----------------------------------------------------------------------
476      open(210, iostat=ierr,file=
477     & trim(datafile)//'/atmosfera_LMD_may.dat')
478      if (ierr.ne.0) then
479        write(*,*)'Error : cannot open file atmosfera_LMD_may.dat '
480        write(*,*)'(in aeronomars/inichim.F)'
481        write(*,*)'It should be in :', trim(datafile),'/'
482        write(*,*)'1) You can change this path in callphys.def with'
483        write(*,*)'   datadir=/path/to/datafiles/'
484        write(*,*)'2) If necessary atmosfera_LMD_may.dat (and others)'
485        write(*,*)'   can be obtained online on:'
486        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
487         STOP
488      endif
489      open(220, iostat=ierr,file=
490     & trim(datafile)//'/atmosfera_LMD_min.dat')
491      if (ierr.ne.0) then
492        write(*,*)'Error : cannot open file atmosfera_LMD_min.dat '
493        write(*,*)'(in aeronomars/inichim.F)'
494        write(*,*)'It should be in :', trim(datafile),'/'
495        write(*,*)'1) You can change this path in callphys.def with'
496        write(*,*)'   datadir=/path/to/datafiles/'
497        write(*,*)'2) If necessary atmosfera_LMD_min.dat (and others)'
498        write(*,*)'   can be obtained online on:'
499        write(*,*)' http://www.lmd.jussieu.fr/~forget/datagcm/datafile'
500         STOP
501      endif
502
503      read(210,*) tmp
504      read(220,*) tmp
505        do l=1,252
506          read(210,112) tfile(l),zfile(l),nxf(l),(vmrf(l,n),n=1,7)
507          zfile(l)=zfile(l)*100.              !pressure (Pa)
508          sigfile(l)=zfile(l)/zfile(1)       
509        enddo
510        do l=1,252
511          read(220,113)zfilemin(l),(vmrf(l,n),n=8,14)
512          zfilemin(l)=zfilemin(l)*1000.       !height (m)
513          do n=1,14
514            nf(l,n)=vmrf(l,n)*nxf(l)
515          enddo
516        enddo
517        close(210)
518        close(220)
519
520c flag thermo set to 1 or 0 in newstart
521c inil=33 for initialisation of thermosphere only       
522c inil=1 for initialisation of all layers       
523        if (thermo.eq.1) then
524        inil=33
525        else
526        inil=1
527        endif
528
529! Initialization of tracers
530
531      do i=1,iip1
532       do j=1,jjp1
533        do l=inil,llm
534
535c          zgcm=sig(l)
536          zgcm=sig(l)*ps(i,j)
537          densconc=0.
538          nt=0.
539
540          do n=1,14
541            call intrplf(zgcm,vmrint,zfile,nf(1,n),252)
542c            call intrplf(zgcm,vmrint,sigfile,nf(1,n),252)
543            ni(n)=vmrint
544
545            densconc=densconc+ni(n)*mmolf(n)
546            nt=nt+ni(n)
547          enddo
548
549          mmean=densconc/nt                     ! in amu
550
551          if (nq.ne.nbqchem) then ! don't initialize water vapour
552            do n=1,nq-1
553            pq(i,j,l,nqchem(n))=
554     &                 ni(aa(nqchem(n)))/nt*mmol(nqchem(n))/mmean
555            if(i.eq.iip1) pq(i,j,l,nqchem(n))=pq(1,j,l,nqchem(n))
556            enddo
557            if (water .and. (nq.gt.nbqchem-2)) then
558              pq(i,j,l,nqchem(nq)) = 0.
559              if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
560            else
561              pq(i,j,l,nqchem(nq))=
562     &              ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean
563              if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
564            endif
565          endif ! of if (nq.ne.nbqchem)
566         
567          if (nq.eq.nbqchem) then ! also initialize water vapour
568            do n=1,nq-2
569              pq(i,j,l,nqchem(n))=
570     &                   ni(aa(nqchem(n)))/nt*mmol(nqchem(n))/mmean
571              if(i.eq.iip1) pq(i,j,l,nqchem(n))=pq(1,j,l,nqchem(n))
572            enddo
573            if (water) then
574              pq(i,j,l,nqchem(nq-1)) = 0.
575              if(i.eq.iip1) then
576                pq(i,j,l,nqchem(nq-1))=pq(1,j,l,nqchem(nq-1))
577              endif
578              pq(i,j,l,nqchem(nq))=
579     &                   ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean
580              if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
581            else
582              do n=nq-1,nq
583                pq(i,j,l,nqchem(nq))=
584     &                   ni(aa(nqchem(nq)))/nt*mmol(nqchem(nq))/mmean
585                if(i.eq.iip1) pq(i,j,l,nqchem(nq))=pq(1,j,l,nqchem(nq))
586              enddo
587            endif ! of if (water)
588          endif ! of if (nq.eq.nqmx)
589        enddo     !nlayer
590       enddo      !ngrid_j
591      enddo       !ngrid_i
592
593cccccccccccccccccccccccccccccc     
594c  make sure that sum of q = 1     
595c dominent species is = 1 - sum(all other species)     
596cccccccccccccccccccccccccccccc     
597!      iqmax=nqchem_min
598       iqmax=1
599     
600!      if ((nqmx-nqchem_min).gt.10) then
601      if (nbqchem.gt.10) then
602       do l=1,llm
603        do j=1,jjp1
604          do i=1,iip1
605!           do iq=nqchem_min,nqmx
606            do iiq=1,nbqchem
607              iq=nqchem(iiq)
608              if ( (pq(i,j,l,iq).gt.pq(i,j,l,iqmax)).and.
609     .           (noms(iq).ne."ice") )  then
610              iqmax=iq
611              endif
612            enddo
613            pq(i,j,l,iqmax)=1.
614            qtot(i,j,l)=0
615!           do iq=nqchem_min,nqmx
616            do iiq=1,nbqchem
617             iq=nqchem(iiq)
618             if ( (iq.ne.iqmax).and.
619     .           (noms(iq).ne."ice") ) then       
620               pq(i,j,l,iqmax)=pq(i,j,l,iqmax)-pq(i,j,l,iq)       
621             endif
622            enddo !iq
623!            do iq=nqchem_min,nqmx
624            do iiq=1,nbqchem
625             iq=nqchem(iiq)
626              if (noms(iq).ne."ice") then
627                qtot(i,j,l)=qtot(i,j,l)+pq(i,j,l,iq)
628              endif
629c            if (i.eq.1.and.j.eq.1.and.l.eq.1) write(*,*) 'qtot(i,j,l)',
630c     $    qtot(i,j,l)
631            enddo !iq
632            if (i.eq.1.and.j.eq.1) write(*,*) 'inichim_newstart: ',
633     $        'qtot(i,j,l)=',qtot(i,j,l)
634          enddo !i   
635         enddo !j   
636       enddo !l 
637      endif ! of if (nbqchem.gt.10)
638ccccccccccccccccccccccccccccccc
639
64066      format(i2,6(1x,e11.5))
641112     format(7x,f7.3,3x,e12.6,3x,e12.6,7(3x,e12.6))
642113     format(1x,f6.2,7(3x,e12.6))
643
644      end
Note: See TracBrowser for help on using the repository browser.