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

Last change on this file since 2883 was 2577, checked in by emillour, 3 years ago

Mars GCM utilities:
Fixes in the utilities for the picky gfortran 10+ compiler.
JL+EM

  • Property svn:executable set to *
File size: 43.6 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      write(*,*) "Better stop now..."
429      stop
430   endif
431   
432   ! Get the value of 'ndim' for this varriable
433   ierr=NF_INQ_VARNDIMS(nid,varid,ndim)
434   write(*,*) 'ndim', ndim
435   
436   ! Check that it is a 3D or 4D variable
437   if (ndim.lt.3) then
438     write(*,*) "Error:",trim(var(j))," is neither a 3D nor a 4D variable"
439     write(*,*) "We'll skip that variable..."
440     CYCLE ! go directly to the next variable
441   endif
442
443!==============================================================================
444! 3.2 Prepare a few things in order to interpolate/write
445!==============================================================================
446
447   if (ndim==3) then
448      allocate(var3d(lonlen,latlen,timelen,1))
449      allocate(var3dls(lonlen,latlen,Nls,1))
450      allocate(var3dxy(timelen))
451
452      dim(1)=londimout
453      dim(2)=latdimout
454      dim(3)=timedimout
455
456      corner(1)=1
457      corner(2)=1
458      corner(3)=1
459      corner(4)=1
460
461      edges(1)=lonlen
462      edges(2)=latlen
463      edges(3)=Nls
464      edges(4)=1
465
466   else if (ndim==4) then
467      allocate(var3d(lonlen,latlen,altlen,timelen))
468      allocate(var3dls(lonlen,latlen,altlen,Nls))
469      allocate(var3dxy(timelen))
470
471      dim(1)=londimout
472      dim(2)=latdimout
473      dim(3)=altdimout
474      dim(4)=timedimout
475
476      corner(1)=1
477      corner(2)=1
478      corner(3)=1
479      corner(4)=1
480
481      edges(1)=lonlen
482      edges(2)=latlen
483      edges(3)=altlen
484      edges(4)=Nls
485   endif
486
487!==============================================================================
488! 3.3 Write this variable's definition and attributes to the output file
489!==============================================================================
490   units=" "
491   long_name=" "
492   ierr=nf_get_att_text(nid,varid,"long_name",long_name)
493   if (ierr.ne.nf_noerr) then
494   ! if no attribute "long_name", try "title"
495     ierr=nf_get_att_text(nid,varid,"title",long_name)
496   endif
497   ierr=nf_get_att_text(nid,varid,"units",units)
498
499   call def_var(nout,var(j),long_name,units,ndim,dim,varidout,ierr)
500
501!==============================================================================
502! 3.4 Read variable
503!==============================================================================
504
505   ierr = NF_GET_VAR_REAL(nid,varid,var3d)
506   miss=NF_GET_ATT_REAL(nid,varid,"missing_value",[missing])
507   validr=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range)
508
509!==============================================================================
510! 3.6 interpolate or average
511!==============================================================================
512
513! 2D variable  (lon, lat, time)   
514! interpolation of var at timels
515   if (ndim==3) then
516      do x=1,lonlen
517       do y=1,latlen
518!        write(*,*) 'd: x, y', x, y
519        do l=1,timelen
520          var3dxy(l)=var3d(x,y,l,1)
521        enddo
522        do n=1,Nls
523          if(average.eq.'y') then
524            resultat=0.
525            Sls=0 ! (gcm data counter within each Ls timestep)
526            do l=1,timelen
527               if((lsgcm(l).ge.(timels(n)-deltals/2.)).and.   &
528               (lsgcm(l).lt.(timels(n)+deltals/2.))) then
529                  if(var3dxy(l) .ne. missing) then
530                    Sls= Sls +1
531                    resultat = resultat + var3dxy(l)
532                  end if
533               end if
534               if (Sls.ne.0) then
535                  var3dls(x,y,n,1)=resultat/float(Sls)
536               else
537                  var3dls(x,y,n,1)=missing
538               endif
539            enddo
540          else   ! average = 'n'
541            call interpolf(timels(n),resultat,missing,lsgcm,var3dxy,timelen)
542            var3dls(x,y,n,1)=resultat
543          endif
544        enddo
545       enddo
546      enddo
547! 3D variable (lon, lat, alt, time)     
548! interpolation of var at timels
549   else if (ndim==4) then
550      do x=1,lonlen
551       do y=1,latlen
552        do k=1,altlen
553        do l=1,timelen
554        var3dxy(l)=var3d(x,y,k,l)
555        enddo
556        do n=1,Nls
557          if(average.eq.'y') then
558            resultat=0.
559            Sls=0 ! (gcm data counter within each Ls timestep)
560            do l=1,timelen
561               if((lsgcm(l).ge.(timels(n)-deltals/2.)).and.   &
562               (lsgcm(l).lt.(timels(n)+deltals/2.))) then
563                  if(var3dxy(l) .ne. missing) then
564                    Sls= Sls +1
565                    resultat = resultat + var3dxy(l)
566                  end if
567               end if
568               if (Sls.ne.0) then
569                  var3dls(x,y,k,n)=resultat/float(Sls)
570               else
571                  var3dls(x,y,k,n)=missing
572               endif
573            enddo
574          else   ! average = 'n'
575             call interpolf(timels(n),resultat,missing,lsgcm,var3dxy,timelen)
576             var3dls(x,y,k,n)=resultat
577          endif
578         enddo
579        enddo
580       enddo
581      enddo
582   endif
583
584!==============================================================================
585! 3.7 Write variable to output file
586!==============================================================================
587
588   ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3dls)
589
590   if (ierr.ne.NF_NOERR) then
591     write(*,*) 'PUT_VAR ERROR: ',NF_STRERROR(ierr)
592     stop
593   endif
594
595! In case there is a "missing value" attribute and "valid range"
596   if (miss.eq.NF_NOERR) then
597     call missing_value(nout,varidout,missing) 
598   endif
599
600   deallocate(var3d)
601   deallocate(var3dls)
602   deallocate(var3dxy)
603
604enddo ! of do j=1,nbvar loop
605
606deallocate(time)
607deallocate(lsgcm)
608
609! close input and output files
610ierr=nf_close(nid)
611ierr=NF_CLOSE(nout)
612
613
614!==============================================================================
615! 4. Build a companion file 'lslin.ctl', so that output file can be
616!    processed by Grads
617!==============================================================================
618
619!  ----------------------------------------------------
620!  Writing ctl file to directly read Ls coordinate in grads
621!  (because of bug in grads that refuse to read date like 0089 in .nc files)
622write(*,*)""
623write(*,*) "Writing a controle file to directly read Ls coordinate with Grads..."
624
625 open(33,file=infile(1:len_trim(infile)-3)//"_Ls.ctl")
626 date= (timels(1)-int(timels(1)))*365.
627 mon='jan'
628 if(date.ge.32) mon='feb'
629 if(date.ge.60) mon='mar'
630 if(date.ge.91) mon='apr'
631 if(date.ge.121) mon='may'
632 if(date.ge.152) mon='jun'
633 if(date.ge.182) mon='jul'
634 if(date.ge.213) mon='aug'
635 if(date.ge.244) mon='sep'
636 if(date.ge.274) mon='oct'
637 if(date.ge.305) mon='nov'
638 if(date.ge.335) mon='dec'
639write(33,98) "^"//outfile
64098 format("DSET ",a)
641  write(33,99) Nls, 15,mon, int(timels(1)),nint(deltals*12),'mo'
64299 format("TDEF Time ",i5," LINEAR ", i2,a3,i4.4,1x,i2,a2)
643   
644deallocate(timels)
645   
646write(*,*)""   
647write(*,*) "Program completed !"
648
649contains
650
651!******************************************************************************
652
653subroutine initiate(filename,lat,lon,alt,ctl,latlen,lonlen,altlen,ctllen,&
654                     altlong_name,altunits,altpositive,&
655                     nout,latdimout,londimout,altdimout,timedimout,timevarout)
656!==============================================================================
657! Purpose:
658! Create and initialize a data file (NetCDF format)
659!==============================================================================
660! Remarks:
661! The NetCDF file (created in this subroutine) remains open
662!==============================================================================
663
664implicit none
665
666include "netcdf.inc" ! NetCDF definitions
667
668!==============================================================================
669! Arguments:
670!==============================================================================
671character (len=*), intent(in):: filename
672! filename(): the file's name
673integer,intent(in) ::latlen,lonlen,altlen,ctllen
674real, intent(in):: lat(latlen)
675! lat(): latitude
676real, intent(in):: lon(lonlen)
677! lon(): longitude
678real, intent(in):: alt(altlen)
679! alt(): altitude
680real, intent(in):: ctl(ctllen)
681! ctl(): controle
682character (len=*), intent(in) :: altlong_name
683! altlong_name(): [netcdf] altitude "long_name" attribute
684character (len=*), intent(in) :: altunits
685! altunits(): [netcdf] altitude "units" attribute
686character (len=*), intent(in) :: altpositive
687! altpositive(): [netcdf] altitude "positive" attribute
688integer, intent(out):: nout
689! nout: [netcdf] file ID
690integer, intent(out):: latdimout
691! latdimout: [netcdf] lat() (i.e.: latitude)  ID
692integer, intent(out):: londimout
693! londimout: [netcdf] lon()  ID
694integer, intent(out):: altdimout
695! altdimout: [netcdf] alt()  ID
696integer, intent(out):: timedimout
697! timedimout: [netcdf] "Time"  ID
698integer, intent(out):: timevarout
699! timevarout: [netcdf] "Time" (considered as a variable) ID
700
701!==============================================================================
702! Local variables:
703!==============================================================================
704!integer :: latdim,londim,altdim,timedim
705integer :: nvarid,ierr
706integer :: ctldimout
707! nvarid: [netcdf] ID of a variable
708! ierr: [netcdf]  return error code (from called subroutines)
709
710!==============================================================================
711! 1. Create (and open) output file
712!==============================================================================
713write(*,*) "creating "//trim(adjustl(filename))//'...'
714ierr = NF_CREATE(filename,IOR(NF_CLOBBER,NF_64BIT_OFFSET),nout)
715! NB: setting NF_CLOBBER mode means that it's OK to overwrite an existing file
716if (ierr.NE.NF_NOERR) then
717   WRITE(*,*)'ERROR: Impossible to create the file ',trim(filename)
718   write(*,*) NF_STRERROR(ierr)
719   stop
720endif
721
722!==============================================================================
723! 2. Define/write "dimensions" and get their IDs
724!==============================================================================
725
726ierr = NF_DEF_DIM(nout, "latitude", latlen, latdimout)
727if (ierr.NE.NF_NOERR) then
728   WRITE(*,*)'initiate: error failed to define dimension <latitude>'
729   write(*,*) NF_STRERROR(ierr)
730   stop
731endif
732ierr = NF_DEF_DIM(nout, "longitude", lonlen, londimout)
733if (ierr.NE.NF_NOERR) then
734   WRITE(*,*)'initiate: error failed to define dimension <longitude>'
735   write(*,*) NF_STRERROR(ierr)
736   stop
737endif
738ierr = NF_DEF_DIM(nout, "altitude", altlen, altdimout)
739if (ierr.NE.NF_NOERR) then
740   WRITE(*,*)'initiate: error failed to define dimension <altitude>'
741   write(*,*) NF_STRERROR(ierr)
742   stop
743endif
744if (size(ctl).ne.0) then
745  ierr = NF_DEF_DIM(nout, "index", ctllen, ctldimout)
746  if (ierr.NE.NF_NOERR) then
747    WRITE(*,*)'initiate: error failed to define dimension <index>'
748    write(*,*) NF_STRERROR(ierr)
749    stop
750  endif
751endif
752ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout)
753if (ierr.NE.NF_NOERR) then
754   WRITE(*,*)'initiate: error failed to define dimension <Time>'
755   write(*,*) NF_STRERROR(ierr)
756   stop
757endif
758
759! End netcdf define mode
760ierr = NF_ENDDEF(nout)
761if (ierr.NE.NF_NOERR) then
762   WRITE(*,*)'initiate: error failed to switch out of define mode'
763   write(*,*) NF_STRERROR(ierr)
764   stop
765endif
766
767!==============================================================================
768! 3. Write "Time" and its attributes
769!==============================================================================
770
771call def_var(nout,"Time","Ls (Solar Longitude)","degrees",1,&
772             (/timedimout/),timevarout,ierr)
773! Switch to netcdf define mode
774ierr = NF_REDEF (nout)
775ierr = NF_ENDDEF(nout)
776
777!==============================================================================
778! 4. Write "latitude" (data and attributes)
779!==============================================================================
780
781call def_var(nout,"latitude","latitude","degrees_north",1,&
782             (/latdimout/),nvarid,ierr)
783
784ierr = NF_PUT_VAR_REAL (nout,nvarid,lat)
785if (ierr.NE.NF_NOERR) then
786   WRITE(*,*)'initiate: error failed writing variable <latitude>'
787   write(*,*) NF_STRERROR(ierr)
788   stop
789endif
790
791!==============================================================================
792! 4. Write "longitude" (data and attributes)
793!==============================================================================
794
795call def_var(nout,"longitude","East longitude","degrees_east",1,&
796             (/londimout/),nvarid,ierr)
797
798ierr = NF_PUT_VAR_REAL (nout,nvarid,lon)
799if (ierr.NE.NF_NOERR) then
800   WRITE(*,*)'initiate: error failed writing variable <longitude>'
801   write(*,*) NF_STRERROR(ierr)
802   stop
803endif
804
805!==============================================================================
806! 5. Write "altitude" (data and attributes)
807!==============================================================================
808
809! Switch to netcdf define mode
810ierr = NF_REDEF (nout)
811
812ierr = NF_DEF_VAR (nout,"altitude",NF_FLOAT,1,[altdimout],nvarid)
813
814ierr = NF_PUT_ATT_TEXT (nout,nvarid,'long_name',len_trim(adjustl(altlong_name)),adjustl(altlong_name))
815ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',len_trim(adjustl(altunits)),adjustl(altunits))
816ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',len_trim(adjustl(altpositive)),adjustl(altpositive))
817
818! End netcdf define mode
819ierr = NF_ENDDEF(nout)
820
821ierr = NF_PUT_VAR_REAL (nout,nvarid,alt)
822if (ierr.NE.NF_NOERR) then
823   WRITE(*,*)'initiate: error failed writing variable <altitude>'
824   write(*,*) NF_STRERROR(ierr)
825   stop
826endif
827
828!==============================================================================
829! 6. Write "controle" (data and attributes)
830!==============================================================================
831
832if (size(ctl).ne.0) then
833   ! Switch to netcdf define mode
834   ierr = NF_REDEF (nout)
835
836   ierr = NF_DEF_VAR (nout,"controle",NF_FLOAT,1,[ctldimout],nvarid)
837
838   ierr = NF_PUT_ATT_TEXT (nout,nvarid,"long_name",18,"Control parameters")
839
840   ! End netcdf define mode
841   ierr = NF_ENDDEF(nout)
842
843   ierr = NF_PUT_VAR_REAL (nout,nvarid,ctl)
844   if (ierr.NE.NF_NOERR) then
845      WRITE(*,*)'initiate: error failed writing variable <controle>'
846      write(*,*) NF_STRERROR(ierr)
847      stop
848   endif
849endif
850
851end Subroutine initiate
852
853!******************************************************************************
854subroutine init2(infid,lonlen,latlen,altlen,altdimid, &
855                 outfid,londimout,latdimout,altdimout)
856!==============================================================================
857! Purpose:
858! Copy aps(), bps() and phisinit() from input file to output file
859!==============================================================================
860! Remarks:
861! The NetCDF files must be open
862!==============================================================================
863
864implicit none
865
866include "netcdf.inc" ! NetCDF definitions
867
868!==============================================================================
869! Arguments:
870!==============================================================================
871integer, intent(in) :: infid  ! NetCDF output file ID
872integer, intent(in) :: lonlen ! # of grid points along longitude
873integer, intent(in) :: latlen ! # of grid points along latitude
874integer, intent(in) :: altlen ! # of grid points along altitude
875integer, intent(in) :: altdimid ! ID of altitude dimension
876integer, intent(in) :: outfid ! NetCDF output file ID
877integer, intent(in) :: londimout ! longitude dimension ID
878integer, intent(in) :: latdimout ! latitude dimension ID
879integer, intent(in) :: altdimout ! altitude dimension ID
880!==============================================================================
881! Local variables:
882!==============================================================================
883real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates
884real,dimension(:,:),allocatable :: phisinit ! Ground geopotential
885integer :: apsid,bpsid,phisinitid
886integer :: ierr
887integer :: tmpdimid ! temporary dimension ID
888integer :: tmpvarid ! temporary variable ID
889integer :: tmplen ! temporary variable length
890logical :: phis, aps_ok, bps_ok ! are "phisinit" "aps" "bps" available ?
891
892
893!==============================================================================
894! 1. Read data from input file
895!==============================================================================
896
897! hybrid coordinate aps
898  allocate(aps(altlen),stat=ierr)
899  if (ierr.ne.0) then
900    write(*,*) "init2: failed to allocate aps(altlen)"
901    stop
902  endif
903
904ierr=NF_INQ_VARID(infid,"aps",tmpvarid)
905if (ierr.ne.NF_NOERR) then
906  write(*,*) "oops Failed to get aps ID. OK"
907  aps_ok=.false.
908else
909  ! Check that aps() is the right size (it most likely won't be if
910  ! zrecast has be used to generate the input file)
911  ierr=NF_INQ_VARDIMID(infid,tmpvarid,tmpdimid)
912  ierr=NF_INQ_DIMLEN(infid,tmpdimid,tmplen)
913  if (tmplen.ne.altlen) then
914    write(*,*) "tmplen,altlen", tmpdimid, altdimid
915    write(*,*) "init2: wrong dimension size for aps(), skipping it ..."
916    aps_ok=.false.
917  else
918    ierr=NF_GET_VAR_REAL(infid,tmpvarid,aps)
919    if (ierr.ne.NF_NOERR) then
920     write(*,*) "init2 error: Failed reading aps"
921     stop
922     aps_ok=.false.
923    endif
924    aps_ok=.true.
925  endif
926endif
927
928! hybrid coordinate bps
929  allocate(bps(altlen),stat=ierr)
930  if (ierr.ne.0) then
931    write(*,*) "init2: failed to allocate bps(altlen)"
932    stop
933  endif
934
935ierr=NF_INQ_VARID(infid,"bps",tmpvarid)
936if (ierr.ne.NF_NOERR) then
937  write(*,*) "oops: Failed to get bps ID. OK"
938  bps_ok=.false.
939else
940  ! Check that bps() is the right size
941  ierr=NF_INQ_VARDIMID(infid,tmpvarid,tmpdimid)
942  ierr=NF_INQ_DIMLEN(infid,tmpdimid,tmplen)
943  if (tmplen.ne.altlen) then
944    write(*,*) "init2: wrong dimension size for bps(), skipping it ..."
945    bps_ok=.false.
946  else
947    ierr=NF_GET_VAR_REAL(infid,tmpvarid,bps)
948    if (ierr.ne.NF_NOERR) then
949      write(*,*) "init2 Error: Failed reading bps"
950      stop
951      bps_ok=.false.
952    endif
953    bps_ok=.true.
954  endif
955endif
956
957! ground geopotential phisinit
958allocate(phisinit(lonlen,latlen),stat=ierr)
959if (ierr.ne.0) then
960  write(*,*) "init2: failed to allocate phisinit(lonlen,latlen)"
961  stop
962endif
963ierr=NF_INQ_VARID(infid,"phisinit",tmpvarid)
964if (ierr.ne.NF_NOERR) then
965  write(*,*) "Failed to get phisinit ID. OK"
966  phisinit = 0.
967  phis = .false.
968else
969  ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
970  if (ierr.ne.NF_NOERR) then
971    write(*,*) "init2 Error: Failed reading phisinit"
972    stop
973  endif
974  phis = .true.
975endif
976
977
978!==============================================================================
979! 2. Write
980!==============================================================================
981
982!==============================================================================
983! 2.2. Hybrid coordinates aps() and bps()
984!==============================================================================
985
986IF (aps_ok) then
987  ! define aps
988  call def_var(outfid,"aps","hybrid pressure at midlayers"," ",1,&
989               (/altdimout/),apsid,ierr)
990  if (ierr.ne.NF_NOERR) then
991    write(*,*) "Error: Failed to def_var aps"
992    stop
993  endif
994
995  ! write aps
996  ierr=NF_PUT_VAR_REAL(outfid,apsid,aps)
997  if (ierr.ne.NF_NOERR) then
998    write(*,*) "Error: Failed to write aps"
999    stop
1000  endif
1001ENDIF ! of IF (aps_ok)
1002
1003IF (bps_ok) then
1004  ! define bps
1005  call def_var(outfid,"bps","hybrid sigma at midlayers"," ",1,&
1006               (/altdimout/),bpsid,ierr)
1007  if (ierr.ne.NF_NOERR) then
1008    write(*,*) "Error: Failed to def_var bps"
1009    stop
1010  endif
1011
1012  ! write bps
1013  ierr=NF_PUT_VAR_REAL(outfid,bpsid,bps)
1014  if (ierr.ne.NF_NOERR) then
1015    write(*,*) "Error: Failed to write bps"
1016    stop
1017  endif
1018ENDIF ! of IF (bps_ok)
1019
1020!==============================================================================
1021! 2.2. phisinit()
1022!==============================================================================
1023
1024IF (phis) THEN
1025
1026  !define phisinit
1027  call def_var(outfid,"phisinit","Ground level geopotential"," ",2,&
1028              (/londimout,latdimout/),phisinitid,ierr)
1029  if (ierr.ne.NF_NOERR) then
1030     write(*,*) "Error: Failed to def_var phisinit"
1031     stop
1032  endif
1033
1034  ! write phisinit
1035  ierr=NF_PUT_VAR_REAL(outfid,phisinitid,phisinit)
1036  if (ierr.ne.NF_NOERR) then
1037    write(*,*) "Error: Failed to write phisinit"
1038    stop
1039  endif
1040
1041ENDIF ! of IF (phis)
1042
1043
1044! Cleanup
1045deallocate(aps)
1046deallocate(bps)
1047deallocate(phisinit)
1048
1049end subroutine init2
1050
1051!******************************************************************************
1052subroutine def_var(nid,name,long_name,units,nbdim,dim,nvarid,ierr)
1053!==============================================================================
1054! Purpose: Write a variable (i.e: add a variable to a dataset)
1055! called "name"; along with its attributes "long_name", "units"...
1056! to a file (following the NetCDF format)
1057!==============================================================================
1058! Remarks:
1059! The NetCDF file must be open
1060!==============================================================================
1061
1062implicit none
1063
1064include "netcdf.inc" ! NetCDF definitions
1065
1066!==============================================================================
1067! Arguments:
1068!==============================================================================
1069integer, intent(in) :: nid
1070! nid: [netcdf] file ID #
1071character (len=*), intent(in) :: name
1072! name(): [netcdf] variable's name
1073character (len=*), intent(in) :: long_name
1074! long_name(): [netcdf] variable's "long_name" attribute
1075character (len=*), intent(in) :: units
1076! unit(): [netcdf] variable's "units" attribute
1077integer, intent(in) :: nbdim
1078! nbdim: number of dimensions of the variable
1079integer, dimension(nbdim), intent(in) :: dim
1080! dim(nbdim): [netcdf] dimension(s) ID(s)
1081integer, intent(out) :: nvarid
1082! nvarid: [netcdf] ID # of the variable
1083integer, intent(out) :: ierr
1084! ierr: [netcdf] subroutines returned error code
1085
1086! Switch to netcdf define mode
1087ierr=NF_REDEF(nid)
1088
1089! Insert the definition of the variable
1090ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid)
1091
1092! Write the attributes
1093ierr=NF_PUT_ATT_TEXT(nid,nvarid,"long_name",len_trim(adjustl(long_name)),adjustl(long_name))
1094ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units))
1095
1096! End netcdf define mode
1097ierr=NF_ENDDEF(nid)
1098
1099end subroutine def_var
1100
1101!******************************************************************************
1102subroutine  missing_value(nout,nvarid,missing)
1103!==============================================================================
1104! Purpose:
1105! Write "valid_range" and "missing_value" attributes (of a given
1106! variable) to a netcdf file
1107!==============================================================================
1108! Remarks:
1109! NetCDF file must be open
1110! Variable (nvarid) ID must be set
1111!==============================================================================
1112
1113implicit none
1114
1115include "netcdf.inc"  ! NetCDF definitions
1116
1117!==============================================================================
1118! Arguments:
1119!==============================================================================
1120integer, intent(in) :: nout
1121! nout: [netcdf] file ID #
1122integer, intent(in) :: nvarid
1123! varid: [netcdf] variable ID #
1124!real, dimension(2), intent(in) :: valid_range
1125! valid_range(2): [netcdf] "valid_range" attribute (min and max)
1126real, intent(in) :: missing
1127! missing: [netcdf] "missing_value" attribute
1128
1129!==============================================================================
1130! Local variables:
1131!==============================================================================
1132integer :: ierr
1133! ierr: [netcdf] subroutine returned error code
1134!      INTEGER nout,nvarid,ierr
1135!      REAL missing
1136!      REAL valid_range(2)
1137
1138! Switch to netcdf dataset define mode
1139ierr = NF_REDEF (nout)
1140if (ierr.ne.NF_NOERR) then
1141   print*,'missing_value: '
1142   print*, NF_STRERROR(ierr)
1143endif
1144
1145
1146!********* valid range not used in Lslin ****************
1147! Write valid_range() attribute
1148!ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range)
1149!
1150!if (ierr.ne.NF_NOERR) then
1151!   print*,'missing_value: valid_range attribution failed'
1152!   print*, NF_STRERROR(ierr)
1153!   stop
1154!endif
1155!*********************************************************
1156
1157! Write "missing_value" attribute
1158ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,missing)
1159if (ierr.NE.NF_NOERR) then
1160   print*, 'missing_value: missing value attribution failed'
1161   print*, NF_STRERROR(ierr)
1162   stop
1163endif
1164
1165! End netcdf dataset define mode
1166ierr = NF_ENDDEF(nout)
1167
1168end subroutine  missing_value
1169
1170!******************************************************************************
1171subroutine sol2ls(sol,Ls)
1172!==============================================================================
1173! Purpose:
1174! Convert a date/time, given in sol (martian day),
1175! into solar longitude date/time, in Ls (in degrees),
1176! where sol=0 is (by definition) the northern hemisphere
1177!  spring equinox (where Ls=0).
1178!==============================================================================
1179! Notes:
1180! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year,
1181! "Ls" will be increased by N*360
1182! Won't work as expected if sol is negative (then again,
1183! why would that ever happen?)
1184!==============================================================================
1185
1186implicit none
1187
1188!==============================================================================
1189! Arguments:
1190!==============================================================================
1191real,intent(in) :: sol
1192real,intent(out) :: Ls
1193
1194!==============================================================================
1195! Local variables:
1196!==============================================================================
1197real year_day,peri_day,timeperi,e_elips,twopi,degrad
1198data year_day /669./            ! # of sols in a martian year
1199data peri_day /485.0/           
1200data timeperi /1.9082314/! True anomaly at vernal equinox = 2*PI-Lsperi
1201data e_elips  /0.093358/
1202data twopi       /6.2831853/    ! 2.*pi
1203data degrad   /57.2957795/      ! pi/180
1204
1205real zanom,xref,zx0,zdx,zteta,zz
1206
1207integer count_years
1208integer iter
1209
1210!==============================================================================
1211! 1. Compute Ls
1212!==============================================================================
1213
1214zz=(sol-peri_day)/year_day
1215zanom=twopi*(zz-nint(zz))
1216xref=abs(zanom)
1217
1218!  The equation zx0 - e * sin (zx0) = xref, solved by Newton
1219zx0=xref+e_elips*sin(xref)
1220do iter=1,20 ! typically, 2 or 3 iterations are enough
1221   zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0))
1222   zx0=zx0+zdx
1223   if(abs(zdx).le.(1.e-7)) then
1224!      write(*,*)'iter:',iter,'     |zdx|:',abs(zdx)
1225      exit
1226   endif
1227enddo
1228
1229if(zanom.lt.0.) zx0=-zx0
1230
1231zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
1232Ls=zteta-timeperi
1233
1234if(Ls.lt.0.) then
1235   Ls=Ls+twopi
1236else
1237   if(Ls.gt.twopi) then
1238      Ls=Ls-twopi
1239   endif
1240endif
1241
1242Ls=degrad*Ls
1243! Ls is now in degrees
1244
1245!==============================================================================
1246! 1. Account for (eventual) years included in input date/time sol
1247!==============================================================================
1248
1249count_years=0 ! initialize
1250zz=sol  ! use "zz" to store (and work on) the value of sol
1251do while (zz.ge.year_day)
1252   count_years=count_years+1
1253   zz=zz-year_day
1254enddo
1255
1256! Add 360 degrees to Ls for every year
1257if (count_years.ne.0) then
1258   Ls=Ls+360.*count_years
1259endif
1260
1261
1262end subroutine sol2ls
1263
1264!******************************************************************************
1265subroutine interpolf(x,y,missing,xd,yd,nd)
1266!==============================================================================
1267! Purpose:
1268! Yield y=f(x), where xd() end yd() are arrays of known values,
1269! using linear interpolation
1270! If x is not included in the interval spaned by xd(), then y is set
1271! to a default value 'missing'
1272! Note:
1273! Array xd() should contain ordered (either increasing or decreasing) abscissas
1274!==============================================================================
1275implicit none
1276!==============================================================================
1277! Arguments:
1278!==============================================================================
1279real,intent(in) :: x ! abscissa at which interpolation is to be done
1280real,intent(in) :: missing ! missing value (if no interpolation is performed)
1281integer :: nd ! size of arrays
1282real,dimension(nd),intent(in) :: xd ! array of known absissas
1283real,dimension(nd),intent(in) :: yd ! array of correponding values
1284
1285real,intent(out) :: y ! interpolated value
1286!==============================================================================
1287! Local variables:
1288!==============================================================================
1289integer :: i
1290
1291! default: set y to 'missing'
1292y=missing
1293
1294   do i=1,nd-1
1295     if (((x.ge.xd(i)).and.(x.le.xd(i+1))).or.&
1296          ((x.le.xd(i)).and.(x.ge.xd(i+1)))) then
1297        if ((yd(i).eq.missing).or.(yd(i+1).eq.missing)) then
1298          ! cannot perform the interpolation if an encompasing value
1299          ! is already set to 'missing'
1300        else
1301          !linear interpolation based on encompassing values
1302          y=yd(i)+(x-xd(i))*(yd(i+1)-yd(i))/(xd(i+1)-xd(i))
1303        endif
1304        exit
1305     endif
1306   enddo
1307
1308
1309end subroutine interpolf
1310
1311!******************************************************************************
1312
1313end program lslin
Note: See TracBrowser for help on using the repository browser.