source: trunk/LMDZ.MARS/util/lslin.F90 @ 2450

Last change on this file since 2450 was 2450, checked in by abierjon, 4 years ago

Mars GCM:
Changes in lslin.F90 :

  • make sure the interpolation takes into account potential missing values in the input file fields, especially when the average/binning option is not used (ticket #68)
  • enable the averaging option even when the Ls timestep is set automatically /!\ previous lslin.def must be updated in consequence !
  • code clean-up, especially remove duplicates on controle field and altitude units
  • replace all the tests "if (ndim.ge.3) then" by a unique test at the beginning of the var loop
  • change "title" attribute into "long_name" by default (ticket #48)

+ add an indicative updated lslin.def
+ update util/README

AB

  • Property svn:executable set to *
File size: 43.5 KB
Line 
1
2
3program lslin
4
5! Program to interpolate data in Solar Longitude linear time coordinate
6! (usable with grads) from Netcdf diagfi or concatnc  files
7! author Y. Wanherdrick, K. Dassas 2004
8! Modified by F.Forget 04/2005
9! More modifications by Ehouarn Millour 12/2005
10! Modified by Ehouarn Millour 10/2007 (changed evaluation of 'start_var'
11! from hard-coded values to a computed value)
12! Read controle field, if available TN, October 2013
13! Allow to chose a starting Ls and to average within Ls timestep a
14! (great to compare with TES and MCS data. F. Forget december 2013)
15
16
17implicit none
18
19include "netcdf.inc" ! NetCDF definitions
20
21character (len=50) :: infile
22! infile(): input file name
23character (len=50), dimension(:), allocatable :: var
24! var(): names of the variables
25character (len=50) :: long_name,units
26! long_name(),units(): [netcdf] "long_name" and "units" attributes
27character (len=50) :: altlong_name,altunits,altpositive
28! altlong_name(): [netcdf] altitude "long_name" attribute
29! altunits(): [netcdf] altitude "units" attribute
30! altpositive(): [netcdf] altitude "positive" attribute
31
32character (len=100) :: outfile
33! outfile(): output file name
34character (len=100) :: vartmp
35! vartmp(): used for temporary storage of various strings
36character (len=1) :: answer
37! answer: user's answers at multiple questions (especially, answer='y'/'Y')
38character (len=1) :: average
39! average: the character 'y' to average within the Ls timestep
40integer :: nid,varid,ierr,miss,validr
41! nid: [netcdf] ID # of input file
42! varid: [netcdf] ID # of a variable
43! ierr: [netcdf] error code (returned by subroutines)
44integer :: nout,varidout
45! nout: [netcdf] ID # of output file
46! varidout: [netcdf] ID # of a variable (to write in the output file)
47integer :: i,j,k,l,x,y,n
48! counters for various loops
49integer :: start_var
50! starting index/ID # from which physical variables are to be found
51integer :: reptime ! Ehouarn: integer or real ?
52! rep_time: starting date/time of the run (given by user)
53integer :: day_ini ! Ehouarn: integer or real ?
54! day_ini: starting date/time of the run stored in input file
55real, dimension(:), allocatable:: lat,lon,alt,time,lsgcm,timels,ctl
56! lat(): latitude coordinates (read from input file)
57! lon(): longitude coordinates (read from input file)
58! alt(): altitude coordinates (read from input file)
59! time(): time coordinates (in "sol", read from input file)
60! lsgcm(): time coordinate (in unevenly spaced "Ls")
61! timels(): new time coordinates (evenly spaced "Ls"; written to output file)
62! ctl(): array, stores controle array
63integer :: latlen,lonlen,altlen,timelen,Nls,Sls,ctllen
64! latlen: # of elements of lat() array
65! lonlen: # of elements of lon() array
66! altvar: # of elements of alt() array
67! timelen: # of elements of time() and lsgcm() arrays
68! Nls: # of elements of timels() array
69integer :: nbvarfile,nbvar,ndim !,nbfile
70! nbvar: # of time-dependent variables
71! nbvarfile: total # of variables (in netcdf file)
72! ndim: [netcdf] # of dimensions (3 or 4) of a variable
73integer :: latdim,londim,altdim,timedim,ctldim
74! latdim: [netcdf] "latitude" dim ID
75! londim: [netcdf] "longitude" dim ID
76! altdim: [netcdf] "altdim" dim ID
77! timedim: [netcdf] "timedim" dim ID
78! ctldim: [netcdf] "controle" dim ID
79integer :: latvar,lonvar,altvar,timevar,ctlvar
80! latvar: [netcdf] ID of "latitude" variable
81! lonvar: [netcdf] ID of "longitude" variable
82! altvar: [netcdf] ID of "altitude" variable
83! timevar: [netcdf] ID of "Time" variable
84integer :: latdimout,londimout,altdimout,timedimout,timevarout
85! latdimout: [netcdf] output latitude (dimension) ID
86! londimout: [netcdf] output longitude (dimension) ID
87! altdimout: [netcdf] output altitude (dimension) ID
88! timedimout: [netcdf] output time (dimension) ID
89! timevarout: [netcdf] ID of output "Time" variable
90integer, dimension(4) :: corner,edges,dim
91! corner(4): [netcdf] start indexes (where block of data will be written)
92! edges(4): [netcdf] lenght (along dimensions) of block of data to write
93! dim(4): [netcdf] lat, long, alt and time dimensions
94real, dimension(:,:,:,:), allocatable :: var3d,var3dls
95! var3d(,,,): 4D array to store a variable (on initial lat/lon/alt/sol grid)
96! var3dls(,,,): 4D array to store a variable (on new lat/lon/alt/Ls grid)
97real, dimension(:), allocatable :: var3dxy
98! var3dxy(): to store the temporal evolution of a variable (at fixed lat/lon/alt)
99real :: deltatsol,deltals,resultat,ls0
100! deltatsol: time step (in sols) of input file data
101! deltals: time step (in Ls) for the data sent to output
102! ls0: first Ls value for the data sent to output
103! resultat: to temporarily store the result of the interpolation
104character (len=3) :: mon
105! mon(3): to store (and write to file) the 3 first letters of a month
106real :: date
107! date: used to compute/build a fake date (for grads output)
108real :: missing
109! to handle missing value when reading /writing files
110real, dimension(2) :: valid_range
111! valid range
112
113!==============================================================================
114! 1. Initialisation step
115!==============================================================================
116
117!==============================================================================
118! 1.1. Get input file name and 'type' (to initialize start_var and reptime)
119!==============================================================================
120
121write(*,*) "which file do you want to use?"
122read(*,'(a50)') infile
123
124!==============================================================================
125! 1.2. Open input file and read/list the variables it contains
126!==============================================================================
127
128write(*,*) "opening "//trim(infile)//"..."
129ierr = NF_OPEN(infile,NF_NOWRITE,nid)
130if (ierr.NE.NF_NOERR) then
131   write(*,*) 'Failed to open file '//infile
132   write(*,*) NF_STRERROR(ierr)
133   stop ""
134endif
135
136! Compute 'start_var', the index from which variables are lon-lat-time
137! and/or lon-lat-alt-time
138! WARNING: We assume here that the ordering of variables in the NetCDF
139! file is such that 0D, 1D and 2D variables are placed BEFORE 3D and 4D
140! variables
141
142i=1 ! dummy initialization to enter loop below
143start_var=0 ! initialization
144do while (i.lt.3)
145  start_var=start_var+1
146  ! request number of dims of variable number 'start_var'
147  ierr=NF_INQ_VARNDIMS(nid,start_var,i)
148  if (ierr.ne.NF_NOERR) then
149    write(*,*) "Failed to get number of dims for variable number:",start_var
150    write(*,*) NF_STRERROR(ierr)
151    stop ""
152  endif
153enddo
154write(*,*) "start_var=",start_var
155
156
157ierr=NF_INQ_NVARS(nid,nbvarfile)
158! nbvarfile is now equal to the (total) number of variables in input file
159                                       
160! Get the variables' names from the input file (and store them in var())
161write(*,*) ""
162nbvar=nbvarfile-(start_var-1)
163allocate(var(nbvar))
164do i=1,nbvar
165   ierr=NF_INQ_VARNAME(nid,i+(start_var-1),var(i))
166   write(*,'(a9,1x,i2,1x,a1,1x,a50)') "variable ",i,":",var(i)
167enddo
168
169!==============================================================================
170! 1.3. Output file name
171!==============================================================================
172outfile=infile(1:len_trim(infile)-3)//"_Ls.nc"
173write(*,*) ""
174write(*,*) "Output file : "//trim(outfile)
175
176
177!==============================================================================
178! 2. Work: read input, build new time coordinate and write it to output
179!==============================================================================
180
181!==============================================================================
182! 2.1. Read (and check) latitude, longitude and altitude from input file
183!==============================================================================
184
185   ierr=NF_INQ_DIMID(nid,"latitude",latdim)
186   ierr=NF_INQ_VARID(nid,"latitude",latvar)
187   if (ierr.NE.NF_NOERR) then
188      write(*,*) 'ERROR: Field <latitude> is missing'
189      write(*,*) NF_STRERROR(ierr)
190      stop "" 
191   endif
192   ierr=NF_INQ_DIMLEN(nid,latdim,latlen)
193!  write(*,*) "latlen: ",latlen
194
195   ierr=NF_INQ_DIMID(nid,"longitude",londim)
196   ierr=NF_INQ_VARID(nid,"longitude",lonvar)
197   if (ierr.NE.NF_NOERR) then
198      write(*,*) 'ERROR: Field <longitude> is missing'
199      write(*,*) NF_STRERROR(ierr)
200      stop ""
201   endif
202   ierr=NF_INQ_DIMLEN(nid,londim,lonlen)
203!  write(*,*) "lonlen: ",lonlen
204
205   ierr=NF_INQ_DIMID(nid,"altitude",altdim)
206   ierr=NF_INQ_VARID(nid,"altitude",altvar)
207   if (ierr.NE.NF_NOERR) then
208      write(*,*) 'ERROR: Field <altitude> is missing'
209      write(*,*) NF_STRERROR(ierr)
210!     stop ""
211       altlen = 1
212   else
213      ierr=NF_INQ_DIMLEN(nid,altdim,altlen)
214     
215    ! Get altitude attributes to handle files with any altitude type
216      ierr=nf_get_att_text(nid,altvar,'long_name',altlong_name)
217      ierr=nf_get_att_text(nid,altvar,'units',altunits)
218      ierr=nf_get_att_text(nid,altvar,'positive',altpositive)
219   endif
220!   write(*,*) "altlen: ",altlen
221
222
223! Allocate lat(), lon() and alt()
224   allocate(lat(latlen))
225   allocate(lon(lonlen))
226   allocate(alt(altlen))
227
228! Read lat(),lon() and alt() from input file
229   ierr = NF_GET_VAR_REAL(nid,latvar,lat)
230   ierr = NF_GET_VAR_REAL(nid,lonvar,lon)
231   ierr = NF_GET_VAR_REAL(nid,altvar,alt)
232   if (ierr.NE.NF_NOERR) then
233       if (altlen.eq.1) alt(1)=0.
234   end if
235
236
237   ierr=NF_INQ_DIMID(nid,"index",ctldim)
238   if (ierr.NE.NF_NOERR) then
239      write(*,*) ' Dimension <index> is missing in file '//trim(infile)
240      ctllen=0
241      !stop "" 
242   endif
243   ierr=NF_INQ_VARID(nid,"controle",ctlvar)
244   if (ierr.NE.NF_NOERR) then
245      write(*,*) 'Field <controle> is missing in file '//trim(infile)
246      ctllen=0
247      !stop ""
248   else
249      ierr=NF_INQ_DIMLEN(nid,ctldim,ctllen)
250   endif
251
252   allocate(ctl(ctllen),stat=ierr)
253   if (ierr.ne.0) then
254     write(*,*) "Error: failed to allocate ctl(ctllen)"
255     stop
256   endif
257
258   if (ctllen .ne. 0) then
259     ierr = NF_GET_VAR_REAL(nid,ctlvar,ctl)
260     if (ierr.ne.0) then
261       write(*,*) "Error: failed to load controle"
262       write(*,*) NF_STRERROR(ierr)
263       stop
264     endif
265   endif ! of if (ctllen .ne. 0)
266
267
268
269!==============================================================================
270! 2.2. Create (and initialize) latitude, longitude and altitude in output file
271!==============================================================================
272
273  ! Initialize output file's lat,lon,alt and time dimensions
274   call initiate(outfile,lat,lon,alt,ctl,latlen,lonlen,altlen,ctllen,&
275         altlong_name,altunits,altpositive,&
276         nout,latdimout,londimout,altdimout,timedimout,timevarout)
277
278  ! Initialize output file's aps,bps and phisinit variables
279   call init2(nid,lonlen,latlen,altlen,altdim,&
280              nout,londimout,latdimout,altdimout)
281
282   write(*,*)""
283!==============================================================================
284! 2.3. Read time from input file
285!==============================================================================
286
287ierr=NF_INQ_DIMID(nid,"Time",timedim)
288if (ierr.NE.NF_NOERR) then
289   write(*,*) 'ERROR: Dimension <Time> is missing'
290   write(*,*) NF_STRERROR(ierr)
291   stop ""
292endif
293ierr=NF_INQ_VARID(nid,"Time",timevar)
294if (ierr.NE.NF_NOERR) then
295   write(*,*) 'ERROR: Field <Time> is missing'
296   write(*,*) NF_STRERROR(ierr)
297   stop ""
298endif
299
300ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
301write(*,*) "timelen: ",timelen
302
303allocate(time(timelen))
304allocate(lsgcm(timelen))
305
306ierr = NF_GET_VAR_REAL(nid,timevar,time)
307
308!==============================================================================
309! 2.4. Initialize day_ini (starting day of the run)
310!==============================================================================
311
312if (ctllen .ne. 0) then ! controle is present in input file
313  write(*,*) "Do you want to specify the beginning day of the file (y/n) ?"
314  write(*,*) "If 'n', the value stored in the controle field will be used."
315  read(*,*) answer
316else ! controle is not present in the input file
317  write(*,*) "There is no field controle in file "//trim(infile)
318  write(*,*) "You have to specify the beginning day of the file."
319  answer="y"
320endif
321
322if ((answer=="y").or.(answer=="Y")) then
323   write(*,*) "Beginning day of the file?"
324   read(*,*) reptime
325endif
326
327write(*,*) 'answer : ',answer
328if ((answer/="y").and.(answer/="Y")) then
329   ! the starting date of the run
330   ! is stored in the "controle" array
331   
332   day_ini = nint(ctl(4))                                               
333   day_ini = modulo(day_ini,669)                                                   
334   write(*,*) 'day_ini', day_ini
335else
336   day_ini= reptime
337   write(*,*) 'day_ini', day_ini
338endif
339
340!==============================================================================
341! 2.5. Build timels() (i.e.: time, but in evenly spaced "Ls" time steps)
342!==============================================================================
343
344! compute time step, in sols, of input dataset
345deltatsol=time(2)-time(1)
346write(*,*) 'deltatsol',deltatsol
347
348! compute time/dates in "Ls"; result stored in lsgcm()
349do i=1,timelen
350   call sol2ls(day_ini+time(i),lsgcm(i))
351  if (lsgcm(i).lt.lsgcm(1)) lsgcm(i)=lsgcm(i) + 360.
352enddo
353
354write(*,*) 'Input data LS range :', lsgcm(1),lsgcm(timelen)
355
356! generate the time step of the new Ls axis for the output
357write(*,*) "Automatic Ls timestep (y/n)?"
358read(*,*) answer
359if ((answer=="y").or.(answer=="Y")) then
360!   *********************************************
361!   Trick of the trade:
362!   Use the value of deltatsol to determine
363!   a suitable value for deltals
364!   *********************************************
365    deltals=1./12.
366    if (0.6*deltatsol.ge.1/6.) deltals=1./6.
367    if (0.6*deltatsol.ge.0.25) deltals=0.25
368    if (0.6*deltatsol.ge.0.5) deltals=0.5
369    if (0.6*deltatsol.ge.0.75) deltals=0.75
370    if (0.6*deltatsol.ge.1) deltals=1.
371    if (0.6*deltatsol.ge.2) deltals=2.
372    if (0.6*deltatsol.ge.3) deltals=3.
373    if (0.6*deltatsol.ge.5) deltals=5.
374    ls0=lsgcm(1) ! First Ls date
375else
376    write(*,*) "New timestep in Ls (deg) :"
377    read(*,*) deltals
378    write(*,*) "First Ls (deg) :"
379    read(*,*) ls0
380    if (ls0.lt.lsgcm(1)) then
381       write(*,*) 'with this file, the earliest Ls is ',lsgcm(1),'let s use it'
382       ls0=lsgcm(1)
383    end if
384    write(*,*) "First Ls date (deg) = ", ls0  ! FF2013
385endif
386
387NLs=int(int((lsgcm(timelen)-ls0)/deltals) +1)   ! FF2013
388allocate(timels(Nls))
389
390! Build timels()
391timels(1) = 0.01*nint(100.*ls0) ! 2 decimals
392do k=2,Nls
393   timels(k) = timels(k-1) + deltals
394end do
395
396write(*,*) 'New timestep in Ls (deg) ', deltals
397write(*,*) 'Output data LS range : ', timels(1),timels(Nls)
398
399! create averaged bins or interpolate at exact Ls ?
400write(*,*) "Average data within the Ls timestep (y/n) or interpolate ? "
401read(*,*) average
402write(*,*)""
403
404!==============================================================================
405! 2.6. Write timels() to output file
406!==============================================================================
407
408do k=1,Nls
409  ierr= NF_PUT_VARA_REAL(nout,timevarout,k,1,timels(k))
410enddo
411
412!==============================================================================
413! 3. Read variables, interpolate them along the new time coordinate
414!    and send the result to output
415!==============================================================================
416
417do j=1,nbvar ! loop on the variables to read/interpolate/write
418
419!==============================================================================
420! 3.1 Check that variable exists, and get some of its attributes
421!==============================================================================
422   write(*,*) "variable ",var(j)
423   ! Get the variable's ID
424   ierr=NF_INQ_VARID(nid,var(j),varid)
425   if (ierr.NE.NF_NOERR) then
426      write(*,*) 'ERROR: Field <',var(j),'> not found'
427      write(*,*) NF_STRERROR(ierr)
428      stop "Better stop now..."
429   endif
430   
431   ! Get the value of 'ndim' for this varriable
432   ierr=NF_INQ_VARNDIMS(nid,varid,ndim)
433   write(*,*) 'ndim', ndim
434   
435   ! Check that it is a 3D or 4D variable
436   if (ndim.lt.3) then
437     write(*,*) "Error:",trim(var(j))," is neither a 3D nor a 4D variable"
438     write(*,*) "We'll skip that variable..."
439     CYCLE ! go directly to the next variable
440   endif
441
442!==============================================================================
443! 3.2 Prepare a few things in order to interpolate/write
444!==============================================================================
445
446   if (ndim==3) then
447      allocate(var3d(lonlen,latlen,timelen,1))
448      allocate(var3dls(lonlen,latlen,Nls,1))
449      allocate(var3dxy(timelen))
450
451      dim(1)=londimout
452      dim(2)=latdimout
453      dim(3)=timedimout
454
455      corner(1)=1
456      corner(2)=1
457      corner(3)=1
458      corner(4)=1
459
460      edges(1)=lonlen
461      edges(2)=latlen
462      edges(3)=Nls
463      edges(4)=1
464
465   else if (ndim==4) then
466      allocate(var3d(lonlen,latlen,altlen,timelen))
467      allocate(var3dls(lonlen,latlen,altlen,Nls))
468      allocate(var3dxy(timelen))
469
470      dim(1)=londimout
471      dim(2)=latdimout
472      dim(3)=altdimout
473      dim(4)=timedimout
474
475      corner(1)=1
476      corner(2)=1
477      corner(3)=1
478      corner(4)=1
479
480      edges(1)=lonlen
481      edges(2)=latlen
482      edges(3)=altlen
483      edges(4)=Nls
484   endif
485
486!==============================================================================
487! 3.3 Write this variable's definition and attributes to the output file
488!==============================================================================
489   units=" "
490   long_name=" "
491   ierr=nf_get_att_text(nid,varid,"long_name",long_name)
492   if (ierr.ne.nf_noerr) then
493   ! if no attribute "long_name", try "title"
494     ierr=nf_get_att_text(nid,varid,"title",long_name)
495   endif
496   ierr=nf_get_att_text(nid,varid,"units",units)
497
498   call def_var(nout,var(j),long_name,units,ndim,dim,varidout,ierr)
499
500!==============================================================================
501! 3.4 Read variable
502!==============================================================================
503
504   ierr = NF_GET_VAR_REAL(nid,varid,var3d)
505   miss=NF_GET_ATT_REAL(nid,varid,"missing_value",missing)
506   validr=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range)
507
508!==============================================================================
509! 3.6 interpolate or average
510!==============================================================================
511
512! 2D variable  (lon, lat, time)   
513! interpolation of var at timels
514   if (ndim==3) then
515      do x=1,lonlen
516       do y=1,latlen
517!        write(*,*) 'd: x, y', x, y
518        do l=1,timelen
519          var3dxy(l)=var3d(x,y,l,1)
520        enddo
521        do n=1,Nls
522          if(average.eq.'y') then
523            resultat=0.
524            Sls=0 ! (gcm data counter within each Ls timestep)
525            do l=1,timelen
526               if((lsgcm(l).ge.(timels(n)-deltals/2.)).and.   &
527               (lsgcm(l).lt.(timels(n)+deltals/2.))) then
528                  if(var3dxy(l) .ne. missing) then
529                    Sls= Sls +1
530                    resultat = resultat + var3dxy(l)
531                  end if
532               end if
533               if (Sls.ne.0) then
534                  var3dls(x,y,n,1)=resultat/float(Sls)
535               else
536                  var3dls(x,y,n,1)=missing
537               endif
538            enddo
539          else   ! average = 'n'
540            call interpolf(timels(n),resultat,missing,lsgcm,var3dxy,timelen)
541            var3dls(x,y,n,1)=resultat
542          endif
543        enddo
544       enddo
545      enddo
546! 3D variable (lon, lat, alt, time)     
547! interpolation of var at timels
548   else if (ndim==4) then
549      do x=1,lonlen
550       do y=1,latlen
551        do k=1,altlen
552        do l=1,timelen
553        var3dxy(l)=var3d(x,y,k,l)
554        enddo
555        do n=1,Nls
556          if(average.eq.'y') then
557            resultat=0.
558            Sls=0 ! (gcm data counter within each Ls timestep)
559            do l=1,timelen
560               if((lsgcm(l).ge.(timels(n)-deltals/2.)).and.   &
561               (lsgcm(l).lt.(timels(n)+deltals/2.))) then
562                  if(var3dxy(l) .ne. missing) then
563                    Sls= Sls +1
564                    resultat = resultat + var3dxy(l)
565                  end if
566               end if
567               if (Sls.ne.0) then
568                  var3dls(x,y,k,n)=resultat/float(Sls)
569               else
570                  var3dls(x,y,k,n)=missing
571               endif
572            enddo
573          else   ! average = 'n'
574             call interpolf(timels(n),resultat,missing,lsgcm,var3dxy,timelen)
575             var3dls(x,y,k,n)=resultat
576          endif
577         enddo
578        enddo
579       enddo
580      enddo
581   endif
582
583!==============================================================================
584! 3.7 Write variable to output file
585!==============================================================================
586
587   ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3dls)
588
589   if (ierr.ne.NF_NOERR) then
590     write(*,*) 'PUT_VAR ERROR: ',NF_STRERROR(ierr)
591     stop ""
592   endif
593
594! In case there is a "missing value" attribute and "valid range"
595   if (miss.eq.NF_NOERR) then
596     call missing_value(nout,varidout,missing) 
597   endif
598
599   deallocate(var3d)
600   deallocate(var3dls)
601   deallocate(var3dxy)
602
603enddo ! of do j=1,nbvar loop
604
605deallocate(time)
606deallocate(lsgcm)
607
608! close input and output files
609ierr=nf_close(nid)
610ierr=NF_CLOSE(nout)
611
612
613!==============================================================================
614! 4. Build a companion file 'lslin.ctl', so that output file can be
615!    processed by Grads
616!==============================================================================
617
618!  ----------------------------------------------------
619!  Writing ctl file to directly read Ls coordinate in grads
620!  (because of bug in grads that refuse to read date like 0089 in .nc files)
621write(*,*)""
622write(*,*) "Writing a controle file to directly read Ls coordinate with Grads..."
623
624 open(33,file=infile(1:len_trim(infile)-3)//"_Ls.ctl")
625 date= (timels(1)-int(timels(1)))*365.
626 mon='jan'
627 if(date.ge.32) mon='feb'
628 if(date.ge.60) mon='mar'
629 if(date.ge.91) mon='apr'
630 if(date.ge.121) mon='may'
631 if(date.ge.152) mon='jun'
632 if(date.ge.182) mon='jul'
633 if(date.ge.213) mon='aug'
634 if(date.ge.244) mon='sep'
635 if(date.ge.274) mon='oct'
636 if(date.ge.305) mon='nov'
637 if(date.ge.335) mon='dec'
638write(33,98) "^"//outfile
63998 format("DSET ",a)
640  write(33,99) Nls, 15,mon, int(timels(1)),nint(deltals*12),'mo'
64199 format("TDEF Time ",i5," LINEAR ", i2,a3,i4.4,1x,i2,a2)
642   
643deallocate(timels)
644   
645write(*,*)""   
646write(*,*) "Program completed !"
647
648contains
649
650!******************************************************************************
651
652subroutine initiate(filename,lat,lon,alt,ctl,latlen,lonlen,altlen,ctllen,&
653                     altlong_name,altunits,altpositive,&
654                     nout,latdimout,londimout,altdimout,timedimout,timevarout)
655!==============================================================================
656! Purpose:
657! Create and initialize a data file (NetCDF format)
658!==============================================================================
659! Remarks:
660! The NetCDF file (created in this subroutine) remains open
661!==============================================================================
662
663implicit none
664
665include "netcdf.inc" ! NetCDF definitions
666
667!==============================================================================
668! Arguments:
669!==============================================================================
670character (len=*), intent(in):: filename
671! filename(): the file's name
672integer,intent(in) ::latlen,lonlen,altlen,ctllen
673real, intent(in):: lat(latlen)
674! lat(): latitude
675real, intent(in):: lon(lonlen)
676! lon(): longitude
677real, intent(in):: alt(altlen)
678! alt(): altitude
679real, intent(in):: ctl(ctllen)
680! ctl(): controle
681character (len=*), intent(in) :: altlong_name
682! altlong_name(): [netcdf] altitude "long_name" attribute
683character (len=*), intent(in) :: altunits
684! altunits(): [netcdf] altitude "units" attribute
685character (len=*), intent(in) :: altpositive
686! altpositive(): [netcdf] altitude "positive" attribute
687integer, intent(out):: nout
688! nout: [netcdf] file ID
689integer, intent(out):: latdimout
690! latdimout: [netcdf] lat() (i.e.: latitude)  ID
691integer, intent(out):: londimout
692! londimout: [netcdf] lon()  ID
693integer, intent(out):: altdimout
694! altdimout: [netcdf] alt()  ID
695integer, intent(out):: timedimout
696! timedimout: [netcdf] "Time"  ID
697integer, intent(out):: timevarout
698! timevarout: [netcdf] "Time" (considered as a variable) ID
699
700!==============================================================================
701! Local variables:
702!==============================================================================
703!integer :: latdim,londim,altdim,timedim
704integer :: nvarid,ierr
705integer :: ctldimout
706! nvarid: [netcdf] ID of a variable
707! ierr: [netcdf]  return error code (from called subroutines)
708
709!==============================================================================
710! 1. Create (and open) output file
711!==============================================================================
712write(*,*) "creating "//trim(adjustl(filename))//'...'
713ierr = NF_CREATE(filename,IOR(NF_CLOBBER,NF_64BIT_OFFSET),nout)
714! NB: setting NF_CLOBBER mode means that it's OK to overwrite an existing file
715if (ierr.NE.NF_NOERR) then
716   WRITE(*,*)'ERROR: Impossible to create the file ',trim(filename)
717   write(*,*) NF_STRERROR(ierr)
718   stop ""
719endif
720
721!==============================================================================
722! 2. Define/write "dimensions" and get their IDs
723!==============================================================================
724
725ierr = NF_DEF_DIM(nout, "latitude", latlen, latdimout)
726if (ierr.NE.NF_NOERR) then
727   WRITE(*,*)'initiate: error failed to define dimension <latitude>'
728   write(*,*) NF_STRERROR(ierr)
729   stop ""
730endif
731ierr = NF_DEF_DIM(nout, "longitude", lonlen, londimout)
732if (ierr.NE.NF_NOERR) then
733   WRITE(*,*)'initiate: error failed to define dimension <longitude>'
734   write(*,*) NF_STRERROR(ierr)
735   stop ""
736endif
737ierr = NF_DEF_DIM(nout, "altitude", altlen, altdimout)
738if (ierr.NE.NF_NOERR) then
739   WRITE(*,*)'initiate: error failed to define dimension <altitude>'
740   write(*,*) NF_STRERROR(ierr)
741   stop ""
742endif
743if (size(ctl).ne.0) then
744  ierr = NF_DEF_DIM(nout, "index", ctllen, ctldimout)
745  if (ierr.NE.NF_NOERR) then
746    WRITE(*,*)'initiate: error failed to define dimension <index>'
747    write(*,*) NF_STRERROR(ierr)
748    stop ""
749  endif
750endif
751ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout)
752if (ierr.NE.NF_NOERR) then
753   WRITE(*,*)'initiate: error failed to define dimension <Time>'
754   write(*,*) NF_STRERROR(ierr)
755   stop ""
756endif
757
758! End netcdf define mode
759ierr = NF_ENDDEF(nout)
760if (ierr.NE.NF_NOERR) then
761   WRITE(*,*)'initiate: error failed to switch out of define mode'
762   write(*,*) NF_STRERROR(ierr)
763   stop ""
764endif
765
766!==============================================================================
767! 3. Write "Time" and its attributes
768!==============================================================================
769
770call def_var(nout,"Time","Ls (Solar Longitude)","degrees",1,&
771             (/timedimout/),timevarout,ierr)
772! Switch to netcdf define mode
773ierr = NF_REDEF (nout)
774ierr = NF_ENDDEF(nout)
775
776!==============================================================================
777! 4. Write "latitude" (data and attributes)
778!==============================================================================
779
780call def_var(nout,"latitude","latitude","degrees_north",1,&
781             (/latdimout/),nvarid,ierr)
782
783ierr = NF_PUT_VAR_REAL (nout,nvarid,lat)
784if (ierr.NE.NF_NOERR) then
785   WRITE(*,*)'initiate: error failed writing variable <latitude>'
786   write(*,*) NF_STRERROR(ierr)
787   stop ""
788endif
789
790!==============================================================================
791! 4. Write "longitude" (data and attributes)
792!==============================================================================
793
794call def_var(nout,"longitude","East longitude","degrees_east",1,&
795             (/londimout/),nvarid,ierr)
796
797ierr = NF_PUT_VAR_REAL (nout,nvarid,lon)
798if (ierr.NE.NF_NOERR) then
799   WRITE(*,*)'initiate: error failed writing variable <longitude>'
800   write(*,*) NF_STRERROR(ierr)
801   stop ""
802endif
803
804!==============================================================================
805! 5. Write "altitude" (data and attributes)
806!==============================================================================
807
808! Switch to netcdf define mode
809ierr = NF_REDEF (nout)
810
811ierr = NF_DEF_VAR (nout,"altitude",NF_FLOAT,1,altdimout,nvarid)
812
813ierr = NF_PUT_ATT_TEXT (nout,nvarid,'long_name',len_trim(adjustl(altlong_name)),adjustl(altlong_name))
814ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',len_trim(adjustl(altunits)),adjustl(altunits))
815ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',len_trim(adjustl(altpositive)),adjustl(altpositive))
816
817! End netcdf define mode
818ierr = NF_ENDDEF(nout)
819
820ierr = NF_PUT_VAR_REAL (nout,nvarid,alt)
821if (ierr.NE.NF_NOERR) then
822   WRITE(*,*)'initiate: error failed writing variable <altitude>'
823   write(*,*) NF_STRERROR(ierr)
824   stop ""
825endif
826
827!==============================================================================
828! 6. Write "controle" (data and attributes)
829!==============================================================================
830
831if (size(ctl).ne.0) then
832   ! Switch to netcdf define mode
833   ierr = NF_REDEF (nout)
834
835   ierr = NF_DEF_VAR (nout,"controle",NF_FLOAT,1,ctldimout,nvarid)
836
837   ierr = NF_PUT_ATT_TEXT (nout,nvarid,"long_name",18,"Control parameters")
838
839   ! End netcdf define mode
840   ierr = NF_ENDDEF(nout)
841
842   ierr = NF_PUT_VAR_REAL (nout,nvarid,ctl)
843   if (ierr.NE.NF_NOERR) then
844      WRITE(*,*)'initiate: error failed writing variable <controle>'
845      write(*,*) NF_STRERROR(ierr)
846      stop ""
847   endif
848endif
849
850end Subroutine initiate
851
852!******************************************************************************
853subroutine init2(infid,lonlen,latlen,altlen,altdimid, &
854                 outfid,londimout,latdimout,altdimout)
855!==============================================================================
856! Purpose:
857! Copy aps(), bps() and phisinit() from input file to output file
858!==============================================================================
859! Remarks:
860! The NetCDF files must be open
861!==============================================================================
862
863implicit none
864
865include "netcdf.inc" ! NetCDF definitions
866
867!==============================================================================
868! Arguments:
869!==============================================================================
870integer, intent(in) :: infid  ! NetCDF output file ID
871integer, intent(in) :: lonlen ! # of grid points along longitude
872integer, intent(in) :: latlen ! # of grid points along latitude
873integer, intent(in) :: altlen ! # of grid points along altitude
874integer, intent(in) :: altdimid ! ID of altitude dimension
875integer, intent(in) :: outfid ! NetCDF output file ID
876integer, intent(in) :: londimout ! longitude dimension ID
877integer, intent(in) :: latdimout ! latitude dimension ID
878integer, intent(in) :: altdimout ! altitude dimension ID
879!==============================================================================
880! Local variables:
881!==============================================================================
882real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates
883real,dimension(:,:),allocatable :: phisinit ! Ground geopotential
884integer :: apsid,bpsid,phisinitid
885integer :: ierr
886integer :: tmpdimid ! temporary dimension ID
887integer :: tmpvarid ! temporary variable ID
888integer :: tmplen ! temporary variable length
889logical :: phis, aps_ok, bps_ok ! are "phisinit" "aps" "bps" available ?
890
891
892!==============================================================================
893! 1. Read data from input file
894!==============================================================================
895
896! hybrid coordinate aps
897  allocate(aps(altlen),stat=ierr)
898  if (ierr.ne.0) then
899    write(*,*) "init2: failed to allocate aps(altlen)"
900    stop
901  endif
902
903ierr=NF_INQ_VARID(infid,"aps",tmpvarid)
904if (ierr.ne.NF_NOERR) then
905  write(*,*) "oops Failed to get aps ID. OK"
906  aps_ok=.false.
907else
908  ! Check that aps() is the right size (it most likely won't be if
909  ! zrecast has be used to generate the input file)
910  ierr=NF_INQ_VARDIMID(infid,tmpvarid,tmpdimid)
911  ierr=NF_INQ_DIMLEN(infid,tmpdimid,tmplen)
912  if (tmplen.ne.altlen) then
913    write(*,*) "tmplen,altlen", tmpdimid, altdimid
914    write(*,*) "init2: wrong dimension size for aps(), skipping it ..."
915    aps_ok=.false.
916  else
917    ierr=NF_GET_VAR_REAL(infid,tmpvarid,aps)
918    if (ierr.ne.NF_NOERR) then
919     stop "init2 error: Failed reading aps"
920     aps_ok=.false.
921    endif
922    aps_ok=.true.
923  endif
924endif
925
926! hybrid coordinate bps
927  allocate(bps(altlen),stat=ierr)
928  if (ierr.ne.0) then
929    write(*,*) "init2: failed to allocate bps(altlen)"
930    stop
931  endif
932
933ierr=NF_INQ_VARID(infid,"bps",tmpvarid)
934if (ierr.ne.NF_NOERR) then
935  write(*,*) "oops: Failed to get bps ID. OK"
936  bps_ok=.false.
937else
938  ! Check that bps() is the right size
939  ierr=NF_INQ_VARDIMID(infid,tmpvarid,tmpdimid)
940  ierr=NF_INQ_DIMLEN(infid,tmpdimid,tmplen)
941  if (tmplen.ne.altlen) then
942    write(*,*) "init2: wrong dimension size for bps(), skipping it ..."
943    bps_ok=.false.
944  else
945    ierr=NF_GET_VAR_REAL(infid,tmpvarid,bps)
946    if (ierr.ne.NF_NOERR) then
947      stop "init2 Error: Failed reading bps"
948      bps_ok=.false.
949    endif
950    bps_ok=.true.
951  endif
952endif
953
954! ground geopotential phisinit
955allocate(phisinit(lonlen,latlen),stat=ierr)
956if (ierr.ne.0) then
957  write(*,*) "init2: failed to allocate phisinit(lonlen,latlen)"
958  stop
959endif
960ierr=NF_INQ_VARID(infid,"phisinit",tmpvarid)
961if (ierr.ne.NF_NOERR) then
962  write(*,*) "Failed to get phisinit ID. OK"
963  phisinit = 0.
964  phis = .false.
965else
966  ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
967  if (ierr.ne.NF_NOERR) then
968    stop "init2 Error: Failed reading phisinit"
969  endif
970  phis = .true.
971endif
972
973
974!==============================================================================
975! 2. Write
976!==============================================================================
977
978!==============================================================================
979! 2.2. Hybrid coordinates aps() and bps()
980!==============================================================================
981
982IF (aps_ok) then
983  ! define aps
984  call def_var(outfid,"aps","hybrid pressure at midlayers"," ",1,&
985               (/altdimout/),apsid,ierr)
986  if (ierr.ne.NF_NOERR) then
987    stop "Error: Failed to def_var aps"
988  endif
989
990  ! write aps
991  ierr=NF_PUT_VAR_REAL(outfid,apsid,aps)
992  if (ierr.ne.NF_NOERR) then
993    stop "Error: Failed to write aps"
994  endif
995ENDIF ! of IF (aps_ok)
996
997IF (bps_ok) then
998  ! define bps
999  call def_var(outfid,"bps","hybrid sigma at midlayers"," ",1,&
1000               (/altdimout/),bpsid,ierr)
1001  if (ierr.ne.NF_NOERR) then
1002    stop "Error: Failed to def_var bps"
1003  endif
1004
1005  ! write bps
1006  ierr=NF_PUT_VAR_REAL(outfid,bpsid,bps)
1007  if (ierr.ne.NF_NOERR) then
1008    stop "Error: Failed to write bps"
1009  endif
1010ENDIF ! of IF (bps_ok)
1011
1012!==============================================================================
1013! 2.2. phisinit()
1014!==============================================================================
1015
1016IF (phis) THEN
1017
1018  !define phisinit
1019  call def_var(outfid,"phisinit","Ground level geopotential"," ",2,&
1020              (/londimout,latdimout/),phisinitid,ierr)
1021  if (ierr.ne.NF_NOERR) then
1022     stop "Error: Failed to def_var phisinit"
1023  endif
1024
1025  ! write phisinit
1026  ierr=NF_PUT_VAR_REAL(outfid,phisinitid,phisinit)
1027  if (ierr.ne.NF_NOERR) then
1028    stop "Error: Failed to write phisinit"
1029  endif
1030
1031ENDIF ! of IF (phis)
1032
1033
1034! Cleanup
1035deallocate(aps)
1036deallocate(bps)
1037deallocate(phisinit)
1038
1039end subroutine init2
1040
1041!******************************************************************************
1042subroutine def_var(nid,name,long_name,units,nbdim,dim,nvarid,ierr)
1043!==============================================================================
1044! Purpose: Write a variable (i.e: add a variable to a dataset)
1045! called "name"; along with its attributes "long_name", "units"...
1046! to a file (following the NetCDF format)
1047!==============================================================================
1048! Remarks:
1049! The NetCDF file must be open
1050!==============================================================================
1051
1052implicit none
1053
1054include "netcdf.inc" ! NetCDF definitions
1055
1056!==============================================================================
1057! Arguments:
1058!==============================================================================
1059integer, intent(in) :: nid
1060! nid: [netcdf] file ID #
1061character (len=*), intent(in) :: name
1062! name(): [netcdf] variable's name
1063character (len=*), intent(in) :: long_name
1064! long_name(): [netcdf] variable's "long_name" attribute
1065character (len=*), intent(in) :: units
1066! unit(): [netcdf] variable's "units" attribute
1067integer, intent(in) :: nbdim
1068! nbdim: number of dimensions of the variable
1069integer, dimension(nbdim), intent(in) :: dim
1070! dim(nbdim): [netcdf] dimension(s) ID(s)
1071integer, intent(out) :: nvarid
1072! nvarid: [netcdf] ID # of the variable
1073integer, intent(out) :: ierr
1074! ierr: [netcdf] subroutines returned error code
1075
1076! Switch to netcdf define mode
1077ierr=NF_REDEF(nid)
1078
1079! Insert the definition of the variable
1080ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid)
1081
1082! Write the attributes
1083ierr=NF_PUT_ATT_TEXT(nid,nvarid,"long_name",len_trim(adjustl(long_name)),adjustl(long_name))
1084ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units))
1085
1086! End netcdf define mode
1087ierr=NF_ENDDEF(nid)
1088
1089end subroutine def_var
1090
1091!******************************************************************************
1092subroutine  missing_value(nout,nvarid,missing)
1093!==============================================================================
1094! Purpose:
1095! Write "valid_range" and "missing_value" attributes (of a given
1096! variable) to a netcdf file
1097!==============================================================================
1098! Remarks:
1099! NetCDF file must be open
1100! Variable (nvarid) ID must be set
1101!==============================================================================
1102
1103implicit none
1104
1105include "netcdf.inc"  ! NetCDF definitions
1106
1107!==============================================================================
1108! Arguments:
1109!==============================================================================
1110integer, intent(in) :: nout
1111! nout: [netcdf] file ID #
1112integer, intent(in) :: nvarid
1113! varid: [netcdf] variable ID #
1114!real, dimension(2), intent(in) :: valid_range
1115! valid_range(2): [netcdf] "valid_range" attribute (min and max)
1116real, intent(in) :: missing
1117! missing: [netcdf] "missing_value" attribute
1118
1119!==============================================================================
1120! Local variables:
1121!==============================================================================
1122integer :: ierr
1123! ierr: [netcdf] subroutine returned error code
1124!      INTEGER nout,nvarid,ierr
1125!      REAL missing
1126!      REAL valid_range(2)
1127
1128! Switch to netcdf dataset define mode
1129ierr = NF_REDEF (nout)
1130if (ierr.ne.NF_NOERR) then
1131   print*,'missing_value: '
1132   print*, NF_STRERROR(ierr)
1133endif
1134
1135
1136!********* valid range not used in Lslin ****************
1137! Write valid_range() attribute
1138!ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range)
1139!
1140!if (ierr.ne.NF_NOERR) then
1141!   print*,'missing_value: valid_range attribution failed'
1142!   print*, NF_STRERROR(ierr)
1143!   stop ""
1144!endif
1145!*********************************************************
1146
1147! Write "missing_value" attribute
1148ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,missing)
1149if (ierr.NE.NF_NOERR) then
1150   print*, 'missing_value: missing value attribution failed'
1151   print*, NF_STRERROR(ierr)
1152   stop ""
1153endif
1154
1155! End netcdf dataset define mode
1156ierr = NF_ENDDEF(nout)
1157
1158end subroutine  missing_value
1159
1160!******************************************************************************
1161subroutine sol2ls(sol,Ls)
1162!==============================================================================
1163! Purpose:
1164! Convert a date/time, given in sol (martian day),
1165! into solar longitude date/time, in Ls (in degrees),
1166! where sol=0 is (by definition) the northern hemisphere
1167!  spring equinox (where Ls=0).
1168!==============================================================================
1169! Notes:
1170! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year,
1171! "Ls" will be increased by N*360
1172! Won't work as expected if sol is negative (then again,
1173! why would that ever happen?)
1174!==============================================================================
1175
1176implicit none
1177
1178!==============================================================================
1179! Arguments:
1180!==============================================================================
1181real,intent(in) :: sol
1182real,intent(out) :: Ls
1183
1184!==============================================================================
1185! Local variables:
1186!==============================================================================
1187real year_day,peri_day,timeperi,e_elips,twopi,degrad
1188data year_day /669./            ! # of sols in a martian year
1189data peri_day /485.0/           
1190data timeperi /1.9082314/! True anomaly at vernal equinox = 2*PI-Lsperi
1191data e_elips  /0.093358/
1192data twopi       /6.2831853/    ! 2.*pi
1193data degrad   /57.2957795/      ! pi/180
1194
1195real zanom,xref,zx0,zdx,zteta,zz
1196
1197integer count_years
1198integer iter
1199
1200!==============================================================================
1201! 1. Compute Ls
1202!==============================================================================
1203
1204zz=(sol-peri_day)/year_day
1205zanom=twopi*(zz-nint(zz))
1206xref=abs(zanom)
1207
1208!  The equation zx0 - e * sin (zx0) = xref, solved by Newton
1209zx0=xref+e_elips*sin(xref)
1210do iter=1,20 ! typically, 2 or 3 iterations are enough
1211   zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0))
1212   zx0=zx0+zdx
1213   if(abs(zdx).le.(1.e-7)) then
1214!      write(*,*)'iter:',iter,'     |zdx|:',abs(zdx)
1215      exit
1216   endif
1217enddo
1218
1219if(zanom.lt.0.) zx0=-zx0
1220
1221zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
1222Ls=zteta-timeperi
1223
1224if(Ls.lt.0.) then
1225   Ls=Ls+twopi
1226else
1227   if(Ls.gt.twopi) then
1228      Ls=Ls-twopi
1229   endif
1230endif
1231
1232Ls=degrad*Ls
1233! Ls is now in degrees
1234
1235!==============================================================================
1236! 1. Account for (eventual) years included in input date/time sol
1237!==============================================================================
1238
1239count_years=0 ! initialize
1240zz=sol  ! use "zz" to store (and work on) the value of sol
1241do while (zz.ge.year_day)
1242   count_years=count_years+1
1243   zz=zz-year_day
1244enddo
1245
1246! Add 360 degrees to Ls for every year
1247if (count_years.ne.0) then
1248   Ls=Ls+360.*count_years
1249endif
1250
1251
1252end subroutine sol2ls
1253
1254!******************************************************************************
1255subroutine interpolf(x,y,missing,xd,yd,nd)
1256!==============================================================================
1257! Purpose:
1258! Yield y=f(x), where xd() end yd() are arrays of known values,
1259! using linear interpolation
1260! If x is not included in the interval spaned by xd(), then y is set
1261! to a default value 'missing'
1262! Note:
1263! Array xd() should contain ordered (either increasing or decreasing) abscissas
1264!==============================================================================
1265implicit none
1266!==============================================================================
1267! Arguments:
1268!==============================================================================
1269real,intent(in) :: x ! abscissa at which interpolation is to be done
1270real,intent(in) :: missing ! missing value (if no interpolation is performed)
1271integer :: nd ! size of arrays
1272real,dimension(nd),intent(in) :: xd ! array of known absissas
1273real,dimension(nd),intent(in) :: yd ! array of correponding values
1274
1275real,intent(out) :: y ! interpolated value
1276!==============================================================================
1277! Local variables:
1278!==============================================================================
1279integer :: i
1280
1281! default: set y to 'missing'
1282y=missing
1283
1284   do i=1,nd-1
1285     if (((x.ge.xd(i)).and.(x.le.xd(i+1))).or.&
1286          ((x.le.xd(i)).and.(x.ge.xd(i+1)))) then
1287        if ((yd(i).eq.missing).or.(yd(i+1).eq.missing)) then
1288          ! cannot perform the interpolation if an encompasing value
1289          ! is already set to 'missing'
1290        else
1291          !linear interpolation based on encompassing values
1292          y=yd(i)+(x-xd(i))*(yd(i+1)-yd(i))/(xd(i+1)-xd(i))
1293        endif
1294        exit
1295     endif
1296   enddo
1297
1298
1299end subroutine interpolf
1300
1301!******************************************************************************
1302
1303end program lslin
Note: See TracBrowser for help on using the repository browser.