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

Last change on this file since 1409 was 1409, checked in by emillour, 10 years ago

Mars GCM:
Update of lslin utility: Added possibility to bin data (instead of interpolating) over the time intervals.
FF

  • Property svn:executable set to *
File size: 34.2 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) :: title,units
26! title(),units(): [netcdf] "title" and "units" attributes
27character (len=50) :: units_alt
28! var(): units of altitude which can be 'km' or 'Pa' (after zrecast)
29character (len=100) :: outfile
30! outfile(): output file name
31character (len=100) :: vartmp
32! vartmp(): used for temporary storage of various strings
33character (len=1) :: answer
34! answer: the character 'y' (or 'Y'), if input file is of 'concatnc' type
35character (len=1) :: average
36! average: the character 'y' to average within the Ls timestep
37integer :: nid,varid,ierr,miss
38! nid: [netcdf] ID # of input file
39! varid: [netcdf] ID # of a variable
40! ierr: [netcdf] error code (returned by subroutines)
41integer :: nout,varidout
42! nout: [netcdf] ID # of output file
43! varidout: [netcdf] ID # of a variable (to write in the output file)
44integer :: i,j,k,l,x,y,n
45! counters for various loops
46integer :: start_var
47! starting index/ID # from which physical variables are to be found
48integer :: reptime ! Ehouarn: integer or real ?
49! rep_time: starting date/time of the run (given by user)
50integer :: day_ini ! Ehouarn: integer or real ?
51! day_ini: starting date/time of the run stored in input file
52integer, parameter :: length=100
53! length: (default) lenght of tab_cntrl()
54real, dimension(length) :: tab_cntrl
55! tab_cntrl(length): array, stored in the input file,
56!                which contains various control parameters.
57real, dimension(:), allocatable:: lat,lon,alt,time,lsgcm,timels
58! lat(): latitude coordinates (read from input file)
59! lon(): longitude coordinates (read from input file)
60! alt(): altitude coordinates (read from input file)
61! time(): time coordinates (in "sol", read from input file)
62! lsgcm(): time coordinate (in unevenly spaced "Ls")
63! timels(): new time coordinates (evenly spaced "Ls"; written to output file)
64integer :: latlen,lonlen,altlen,timelen,Nls,Sls
65! latlen: # of elements of lat() array
66! lonlen: # of elements of lon() array
67! altvar: # of elements of alt() array
68! timelen: # of elements of time() and lsgcm() arrays
69! Nls: # of elements of timels() array
70integer :: nbvarfile,nbvar,ndim !,nbfile
71! nbvar: # of time-dependent variables
72! nbvarfile: total # of variables (in netcdf file)
73! ndim: [netcdf] # of dimensions (3 or 4) of a variable
74integer :: latdim,londim,altdim,timedim
75! latdim: [netcdf] "latitude" dim ID
76! londim: [netcdf] "longitude" dim ID
77! altdim: [netcdf] "altdim" dim ID
78! timedim: [netcdf] "timedim" dim ID
79integer :: latvar,lonvar,altvar,timevar
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!write(*,*) "Is it a concatnc file? (y/n)?"
125write(*,*) "Do you want to specify the beginning day of the file"
126write(*,*) "in case the controle field is not present ? (y/n)?"
127read(*,*) answer
128if ((answer=="y").or.(answer=="Y")) then
129   write(*,*) "Beginning day of the file?"
130   read(*,*) reptime
131!   start_var=8 ! 'concatnc' type of file
132!else
133!   start_var=16 ! 'diagfi' type of file
134! N.B.: start_var is now computed; see below
135endif
136
137
138!==============================================================================
139! 1.2. Open input file and read/list the variables it contains
140!==============================================================================
141
142write(*,*) "opening "//trim(infile)//"..."
143ierr = NF_OPEN(infile,NF_NOWRITE,nid)
144if (ierr.NE.NF_NOERR) then
145   write(*,*) 'Failed to open file '//infile
146   write(*,*) NF_STRERROR(ierr)
147   stop ""
148endif
149
150! Compute 'start_var', the index from which variables are lon-lat-time
151! and/or lon-lat-alt-time
152! WARNING: We assume here that the ordering of variables in the NetCDF
153! file is such that 0D, 1D and 2D variables are placed BEFORE 3D and 4D
154! variables
155
156i=1 ! dummy initialization to enter loop below
157start_var=0 ! initialization
158do while (i.lt.3)
159  start_var=start_var+1
160  ! request number of dims of variable number 'start_var'
161  ierr=NF_INQ_VARNDIMS(nid,start_var,i)
162  if (ierr.ne.NF_NOERR) then
163    write(*,*) "Failed to get number of dims for variable number:",start_var
164    write(*,*) NF_STRERROR(ierr)
165    stop ""
166  endif
167enddo
168write(*,*) "start_var=",start_var
169
170
171ierr=NF_INQ_NVARS(nid,nbvarfile)
172! nbvarfile is now equal to the (total) number of variables in input file
173
174allocate(var(nbvarfile-(start_var-1)))
175
176! Yield the list of variables (to the screen)
177write(*,*)
178do i=start_var,nbvarfile
179   ierr=NF_INQ_VARNAME(nid,i,vartmp)
180   write(*,*) trim(vartmp)
181enddo
182!nbvar=0
183
184! Ehouarn: Redundant (and obsolete) lines ?
185!  nbvar=nbvarfile-(start_var-1)
186!  do j=start_var,nbvarfile
187!  ierr=nf_inq_varname(nid,j,var(j-(start_var-1)))
188!  enddo
189                                       
190! Get the variables' names from the input file (and store them in var())
191nbvar=nbvarfile-(start_var-1)
192do i=1,nbvar
193   ierr=NF_INQ_VARNAME(nid,i+(start_var-1),var(i))
194   write(*,'(a9,1x,i2,1x,a1,1x,a50)') "variable ",i,":",var(i)
195enddo
196
197!==============================================================================
198! 1.3. Output file name
199!==============================================================================
200! outfile="lslin.nc"
201outfile=infile(1:len_trim(infile)-3)//"_Ls.nc"
202 write(*,*) outfile
203
204
205!==============================================================================
206! 2. Work: read input, build new time coordinate and write it to output
207!==============================================================================
208
209!==============================================================================
210! 2.1. Read (and check) latitude, longitude and altitude from input file
211!==============================================================================
212
213   ierr=NF_INQ_DIMID(nid,"latitude",latdim)
214   ierr=NF_INQ_VARID(nid,"latitude",latvar)
215   if (ierr.NE.NF_NOERR) then
216      write(*,*) 'ERROR: Field <latitude> is missing'
217      write(*,*) NF_STRERROR(ierr)
218      stop "" 
219   endif
220   ierr=NF_INQ_DIMLEN(nid,latdim,latlen)
221!  write(*,*) "latlen: ",latlen
222
223   ierr=NF_INQ_DIMID(nid,"longitude",londim)
224   ierr=NF_INQ_VARID(nid,"longitude",lonvar)
225   if (ierr.NE.NF_NOERR) then
226      write(*,*) 'ERROR: Field <longitude> is missing'
227      write(*,*) NF_STRERROR(ierr)
228      stop ""
229   endif
230   ierr=NF_INQ_DIMLEN(nid,londim,lonlen)
231!  write(*,*) "lonlen: ",lonlen
232
233   ierr=NF_INQ_DIMID(nid,"altitude",altdim)
234   ierr=NF_INQ_VARID(nid,"altitude",altvar)
235   if (ierr.NE.NF_NOERR) then
236      write(*,*) 'ERROR: Field <altitude> is missing'
237      write(*,*) NF_STRERROR(ierr)
238!     stop ""
239       altlen = 1
240   else
241      ierr=NF_INQ_DIMLEN(nid,altdim,altlen)
242      units_alt="                                                    "
243      ierr=NF_GET_ATT_TEXT(nid,varid,"units",units_alt)
244   endif
245   write(*,*) "altlen: ",altlen
246
247! Allocate lat(), lon() and alt()
248      allocate(lat(latlen))
249      allocate(lon(lonlen))
250      allocate(alt(altlen))
251
252! Read lat(),lon() and alt() from input file
253
254
255
256
257
258      ierr = NF_GET_VAR_REAL(nid,latvar,lat)
259      ierr = NF_GET_VAR_REAL(nid,lonvar,lon)
260      ierr = NF_GET_VAR_REAL(nid,altvar,alt)
261      if (ierr.NE.NF_NOERR) then
262          if (altlen.eq.1) alt(1)=0.
263      end if
264
265
266!==============================================================================
267! 2.2. Create (and initialize latitude, longitude and altitude in) output file
268!==============================================================================
269
270call initiate(outfile,lat,lon,alt,nout,&
271           latdimout,londimout,altdimout,timedimout,timevarout,units_alt)
272
273!==============================================================================
274! 2.3. Read time from input file
275!==============================================================================
276
277ierr=NF_INQ_DIMID(nid,"Time",timedim)
278if (ierr.NE.NF_NOERR) then
279   write(*,*) 'ERROR: Dimension <Time> is missing'
280   write(*,*) NF_STRERROR(ierr)
281   stop ""
282endif
283ierr=NF_INQ_VARID(nid,"Time",timevar)
284if (ierr.NE.NF_NOERR) then
285   write(*,*) 'ERROR: Field <Time> is missing'
286   write(*,*) NF_STRERROR(ierr)
287   stop ""
288endif
289
290ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
291write(*,*) "timelen: ",timelen
292
293allocate(time(timelen))
294allocate(lsgcm(timelen))
295
296
297
298
299ierr = NF_GET_VAR_REAL(nid,timevar,time)
300
301
302!==============================================================================
303! 2.4. Initialize day_ini (starting day of the run)
304!==============================================================================
305
306
307write(*,*) 'answer',answer
308if ((answer/="y").and.(answer/="Y")) then
309   ! input file is of 'concatnc' type; the starting date of the run
310   ! is stored in the "controle" array
311   ierr = NF_INQ_VARID (nid, "controle", varid)
312   if (ierr .NE. NF_NOERR) then
313      write(*,*) 'ERROR: <controle> variable is missing'
314      write(*,*) NF_STRERROR(ierr)
315      stop ""
316   endif
317
318
319
320
321   ierr = NF_GET_VAR_REAL(nid, varid, tab_cntrl)
322
323
324   if (ierr .NE. NF_NOERR) then
325      write(*,*) 'ERROR: Failed to read <controle>'
326      write(*,*) NF_STRERROR(ierr)
327     stop ""
328   endif
329
330   day_ini = tab_cntrl(4)                                               
331   day_ini = modulo(day_ini,669)                                                   
332   write(*,*) 'day_ini', day_ini
333else
334   day_ini= reptime
335   write(*,*) 'day_ini', day_ini
336endif
337
338!==============================================================================
339! 2.5. Build timels() (i.e.: time, but in evenly spaced "Ls" time steps)
340!==============================================================================
341
342! compute time step, in sols, of input dataset
343deltatsol=time(2)-time(1)
344write(*,*) 'deltatsol',deltatsol
345
346! compute time/dates in "Ls"; result stored in lsgcm()
347do i=1,timelen
348   call sol2ls(day_ini+time(i),lsgcm(i))
349  if (lsgcm(i).lt.lsgcm(1)) lsgcm(i)=lsgcm(i) + 360.
350enddo
351
352write(*,*) 'Input data LS range :', lsgcm(1),lsgcm(timelen)
353
354!******************************************
355write(*,*) "Automatic Ls timetsep (y/n)?"
356read(*,*) answer
357if ((answer=="y").or.(answer=="Y")) then
358!   *********************************************
359!   Trick of the trade:
360!   Use the value of deltatsol to determine
361!   a suitable value for deltals
362!   *********************************************
363    deltals=1./12.
364    if (0.6*deltatsol.ge.1/6.) deltals=1./6.
365    if (0.6*deltatsol.ge.0.25) deltals=0.25
366    if (0.6*deltatsol.ge.0.5) deltals=0.5
367    if (0.6*deltatsol.ge.0.75) deltals=0.75
368    if (0.6*deltatsol.ge.1) deltals=1.
369    if (0.6*deltatsol.ge.2) deltals=2.
370    if (0.6*deltatsol.ge.3) deltals=3.
371    if (0.6*deltatsol.ge.5) deltals=5.
372    ls0=lsgcm(1) ! First Ls date
373else
374    write(*,*) "New timestep in Ls (deg)"
375    read(*,*) deltals
376    write(*,*) "First Ls (deg)"
377    read(*,*) ls0
378    if (ls0.lt.lsgcm(1)) then
379       write(*,*) 'with this file, the earliest Ls is ',lsgcm(1),'let s use it'
380       ls0=lsgcm(1)
381    end if
382    write(*,*) "First Ls date (deg) = ", ls0  ! FF2013
383
384    do while((average.ne.'y').and.(average.ne.'n'))
385       write(*,*) "Average data within the Ls timestep (y/n) or interpolate ? "
386       read(*,*) average
387    end do
388endif
389  NLs=int(Int((lsgcm(timelen)-ls0)/deltals) +1)   ! FF2013
390!********************************************
391allocate(timels(Nls))
392
393! Build timels()
394timels(1) = 0.01*nint(100.*ls0) ! Ehouarn: what the !!!???!!!
395do k=2,Nls
396   timels(k) = timels(k-1) + deltals
397end do
398
399write(*,*) 'timestep in Ls (deg) ', deltals
400write(*,*) 'output data LS range : ', timels(1),timels(Nls)
401
402!==============================================================================
403! 2.6. Write timels() to output file
404!==============================================================================
405
406do k=1,Nls
407
408
409
410ierr= NF_PUT_VARA_REAL(nout,timevarout,k,1,timels(k))
411
412enddo
413
414
415!==============================================================================
416! 3. Read variables, interpolate them along the new time coordinate
417!    and send the result to output
418!==============================================================================
419
420do j=1,nbvar ! loop on the variables to read/interpolate/write
421
422!==============================================================================
423! 3.1 Check that variable exists, and get some of its attributes
424!==============================================================================
425   write(*,*) "variable ",var(j)
426   ! Get the variable's ID
427   ierr=NF_INQ_VARID(nid,var(j),varid)
428   if (ierr.NE.NF_NOERR) then
429      write(*,*) 'ERROR: Field <',var(j),'> not found'
430      write(*,*) NF_STRERROR(ierr)
431      stop ""
432   endif
433   ! Get the value of 'ndim' for this varriable
434   ierr=NF_INQ_VARNDIMS(nid,varid,ndim)
435   write(*,*) 'ndim', ndim
436
437!==============================================================================
438! 3.2 Prepare a few things in order to interpolate/write
439!==============================================================================
440
441   if (ndim==3) then
442      allocate(var3d(lonlen,latlen,timelen,1))
443      allocate(var3dls(lonlen,latlen,Nls,1))
444      allocate(var3dxy(timelen))
445
446      dim(1)=londimout
447      dim(2)=latdimout
448      dim(3)=timedimout
449
450      corner(1)=1
451      corner(2)=1
452      corner(3)=1
453      corner(4)=1
454
455      edges(1)=lonlen
456      edges(2)=latlen
457      edges(3)=Nls
458      edges(4)=1
459
460   else if (ndim==4) then
461      allocate(var3d(lonlen,latlen,altlen,timelen))
462      allocate(var3dls(lonlen,latlen,altlen,Nls))
463      allocate(var3dxy(timelen))
464
465      dim(1)=londimout
466      dim(2)=latdimout
467      dim(3)=altdimout
468      dim(4)=timedimout
469
470      corner(1)=1
471      corner(2)=1
472      corner(3)=1
473      corner(4)=1
474
475      edges(1)=lonlen
476      edges(2)=latlen
477      edges(3)=altlen
478      edges(4)=Nls
479   endif
480
481!==============================================================================
482! 3.3 Write this variable's definition and attributes to the output file
483!==============================================================================
484   units="                                                    "
485   title="                                                    "
486   ierr=NF_GET_ATT_TEXT(nid,varid,"title",title)
487   ierr=NF_GET_ATT_TEXT(nid,varid,"units",units)
488
489   call def_var(nout,var(j),title,units,ndim,dim,varidout,ierr)
490
491!==============================================================================
492! 3.4 Read variable
493!==============================================================================
494
495
496
497
498   ierr = NF_GET_VAR_REAL(nid,varid,var3d)
499   miss = NF_GET_ATT_REAL(nid,varid,"missing_value",missing)
500   miss = NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range)
501
502
503!==============================================================================
504! 3.6 interpolate or average
505!==============================================================================
506
507! 2D variable  (lon, lat, time)   
508! interpolation of var at timels
509   if (ndim==3) then
510      do x=1,lonlen
511       do y=1,latlen
512!        write(*,*) 'd: x, y', x, y
513        do l=1,timelen
514          var3dxy(l)=var3d(x,y,l,1)
515        enddo
516        do n=1,Nls
517          if(average.eq.'y') then
518            resultat=0.
519            Sls=0 ! (gcm data counter within each Ls timestep)
520            do l=1,timelen
521               if((lsgcm(l).ge.(timels(n)-deltals/2.)).and.   &
522               (lsgcm(l).lt.(timels(n)+deltals/2.))) then
523                  if(var3dxy(l) .ne. missing) then
524                    Sls= Sls +1
525                    resultat = resultat + var3dxy(l)
526                  end if
527               end if
528               if (Sls.ne.0) then
529                  var3dls(x,y,n,1)=resultat/float(Sls)
530               else
531                  var3dls(x,y,n,1)=missing
532               endif
533            enddo
534          else   ! average = n
535            call interpolf(timels(n),resultat,lsgcm,var3dxy,timelen)
536            var3dls(x,y,n,1)=resultat
537          endif
538        enddo
539       enddo
540      enddo
541! 3D variable (lon, lat, alt, time)     
542! interpolation of var at timels
543   else if (ndim==4) then
544      do x=1,lonlen
545       do y=1,latlen
546        do k=1,altlen
547        do l=1,timelen
548        var3dxy(l)=var3d(x,y,k,l)
549        enddo
550        do n=1,Nls
551          if(average.eq.'y') then
552            resultat=0.
553            Sls=0 ! (gcm data counter within each Ls timestep)
554            do l=1,timelen
555               if((lsgcm(l).ge.(timels(n)-deltals/2.)).and.   &
556               (lsgcm(l).lt.(timels(n)+deltals/2.))) then
557                  if(var3dxy(l) .ne. missing) then
558                    Sls= Sls +1
559                    resultat = resultat + var3dxy(l)
560                  end if
561               end if
562               if (Sls.ne.0) then
563                  var3dls(x,y,k,n)=resultat/float(Sls)
564               else
565                  var3dls(x,y,k,n)=missing
566               endif
567            enddo
568          else   ! average = n
569             call interpolf(timels(n),resultat,lsgcm,var3dxy,timelen)
570             var3dls(x,y,k,n)=resultat
571          end if
572         enddo
573        enddo
574       enddo
575      enddo
576   endif
577
578!==============================================================================
579! 3.7 Write variable to output file
580!==============================================================================
581
582
583
584
585   ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3dls)
586
587
588   if (ierr.ne.NF_NOERR) then
589      write(*,*) 'PUT_VAR ERROR: ',NF_STRERROR(ierr)
590      stop ""
591   endif
592
593! In case there is a "missing value" attribute and "valid range"
594   if (miss.eq.NF_NOERR) then
595      call missing_value(nout,varidout,valid_range,missing) 
596   end if
597
598   deallocate(var3d)
599   deallocate(var3dls)
600   deallocate(var3dxy)
601
602enddo ! of do j=1,nbvar loop
603
604deallocate(time)
605deallocate(lsgcm)
606
607! close input and output files
608ierr=nf_close(nid)
609ierr=NF_CLOSE(nout)
610
611
612!==============================================================================
613! 4. Build a companion file 'lslin.ctl', so that output file can be
614!    processed by Grads
615!==============================================================================
616
617!  ----------------------------------------------------
618!  Writing ctl file to directly read Ls coordinate in grads
619!  (because of bug in grads that refuse to read date like 0089 in .nc files)
620!open(33,file='lslin.ctl')
621 open(33,file=infile(1:len_trim(infile)-3)//"_Ls.ctl")
622 date= (timels(1)-int(timels(1)))*365.
623 mon='jan'
624 if(date.ge.32) mon='feb'
625 if(date.ge.60) mon='mar'
626 if(date.ge.91) mon='apr'
627 if(date.ge.121) mon='may'
628 if(date.ge.152) mon='jun'
629 if(date.ge.182) mon='jul'
630 if(date.ge.213) mon='aug'
631 if(date.ge.244) mon='sep'
632 if(date.ge.274) mon='oct'
633 if(date.ge.305) mon='nov'
634 if(date.ge.335) mon='dec'
635write(33,98) "^"//outfile
63698 format("DSET ",a)
637  write(33,99) Nls, 15,mon, int(timels(1)),nint(deltals*12),'mo'
63899 format("TDEF Time ",i5," LINEAR ", i2,a3,i4.4,1x,i2,a2)
639   
640   deallocate(timels)
641contains
642
643!************************************
644subroutine initiate (filename,lat,lon,alt,&
645                     nout,latdimout,londimout,altdimout,timedimout,timevarout,units_alt)
646!==============================================================================
647! Purpose:
648! Create and initialize a data file (NetCDF format)
649!==============================================================================
650! Remarks:
651! The NetCDF file (created in this subroutine) remains open
652!==============================================================================
653
654implicit none
655
656include "netcdf.inc" ! NetCDF definitions
657
658!==============================================================================
659! Arguments:
660!==============================================================================
661character (len=*), intent(in):: filename
662! filename(): the file's name
663real, dimension(:), intent(in):: lat
664! lat(): latitude
665real, dimension(:), intent(in):: lon
666! lon(): longitude
667real, dimension(:), intent(in):: alt
668! alt(): altitude
669integer, intent(out):: nout
670! nout: [netcdf] file ID
671integer, intent(out):: latdimout
672! latdimout: [netcdf] lat() (i.e.: latitude)  ID
673integer, intent(out):: londimout
674! londimout: [netcdf] lon()  ID
675integer, intent(out):: altdimout
676! altdimout: [netcdf] alt()  ID
677integer, intent(out):: timedimout
678! timedimout: [netcdf] "Time"  ID
679integer, intent(out):: timevarout
680! timevarout: [netcdf] "Time" (considered as a variable) ID
681character (len=50) :: units_alt
682! var(): units of altitude which can be m or Pa (after zrecast)
683
684
685!==============================================================================
686! Local variables:
687!==============================================================================
688!integer :: latdim,londim,altdim,timedim
689integer :: nvarid,ierr
690! nvarid: [netcdf] ID of a variable
691! ierr: [netcdf]  return error code (from called subroutines)
692
693!==============================================================================
694! 1. Create (and open) output file
695!==============================================================================
696write(*,*) "creating "//trim(adjustl(filename))//'...'
697ierr = NF_CREATE(filename,IOR(NF_CLOBBER,NF_64BIT_OFFSET),nout)
698! NB: setting NF_CLOBBER mode means that it's OK to overwrite an existing file
699if (ierr.NE.NF_NOERR) then
700   WRITE(*,*)'ERROR: Impossible to create the file.'
701   stop ""
702endif
703
704!==============================================================================
705! 2. Define/write "dimensions" and get their IDs
706!==============================================================================
707
708ierr = NF_DEF_DIM(nout, "latitude", size(lat), latdimout)
709ierr = NF_DEF_DIM(nout, "longitude", size(lon), londimout)
710ierr = NF_DEF_DIM(nout, "altitude", size(alt), altdimout)
711ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout)
712
713! End netcdf define mode
714ierr = NF_ENDDEF(nout)
715
716!==============================================================================
717! 3. Write "Time" and its attributes
718!==============================================================================
719
720call def_var(nout,"Time","Time","years since 0000-00-0 00:00:00",1,&
721             (/timedimout/),timevarout,ierr)
722
723!==============================================================================
724! 4. Write "latitude" (data and attributes)
725!==============================================================================
726
727call def_var(nout,"latitude","latitude","degrees_north",1,&
728             (/latdimout/),nvarid,ierr)
729
730ierr = NF_PUT_VAR_REAL (nout,nvarid,lat)
731
732!==============================================================================
733! 4. Write "longitude" (data and attributes)
734!==============================================================================
735
736call def_var(nout,"longitude","East longitude","degrees_east",1,&
737             (/londimout/),nvarid,ierr)
738
739ierr = NF_PUT_VAR_REAL (nout,nvarid,lon)
740
741!==============================================================================
742! 4. Write "altitude" (data and attributes)
743!==============================================================================
744
745! Switch to netcdf define mode
746ierr = NF_REDEF (nout)
747
748ierr = NF_DEF_VAR (nout,"altitude",NF_FLOAT,1,altdimout,nvarid)
749
750ierr = NF_PUT_ATT_TEXT (nout,nvarid,"long_name",8,"altitude")
751!********************* temlporaire ****************
752!ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',2,units_alt)
753!if((units_alt(1:2).eq.'Pa').or.(units_alt(1:2).eq.'pa')) then
754!   ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',2,"down")
755!else
756!   ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',2,"up")
757!end if
758!********************* temlporaire ****************
759ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',2,'Pa')
760   ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',4,'down')
761!********************* temlporaire ****************
762
763! End netcdf define mode
764ierr = NF_ENDDEF(nout)
765
766ierr = NF_PUT_VAR_REAL (nout,nvarid,alt)
767
768end Subroutine initiate
769!************************************
770subroutine def_var(nid,name,title,units,nbdim,dim,nvarid,ierr)
771!==============================================================================
772! Purpose: Write a variable (i.e: add a variable to a dataset)
773! called "name"; along with its attributes "title", "units"...
774! to a file (following the NetCDF format)
775!==============================================================================
776! Remarks:
777! The NetCDF file must be open
778!==============================================================================
779
780implicit none
781
782include "netcdf.inc" ! NetCDF definitions
783
784!==============================================================================
785! Arguments:
786!==============================================================================
787integer, intent(in) :: nid
788! nid: [netcdf] file ID #
789character (len=*), intent(in) :: name
790! name(): [netcdf] variable's name
791character (len=*), intent(in) :: title
792! title(): [netcdf] variable's "title" attribute
793character (len=*), intent(in) :: units
794! unit(): [netcdf] variable's "units" attribute
795integer, intent(in) :: nbdim
796! nbdim: number of dimensions of the variable
797integer, dimension(nbdim), intent(in) :: dim
798! dim(nbdim): [netcdf] dimension(s) ID(s)
799integer, intent(out) :: nvarid
800! nvarid: [netcdf] ID # of the variable
801integer, intent(out) :: ierr
802! ierr: [netcdf] subroutines returned error code
803
804! Switch to netcdf define mode
805ierr=NF_REDEF(nid)
806
807! Insert the definition of the variable
808ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid)
809
810! Write the attributes
811ierr=NF_PUT_ATT_TEXT(nid,nvarid,"title",len_trim(adjustl(title)),adjustl(title))
812ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units))
813
814! End netcdf define mode
815ierr=NF_ENDDEF(nid)
816
817end subroutine def_var
818!************************************
819subroutine sol2ls(sol,Ls)
820!==============================================================================
821! Purpose:
822! Convert a date/time, given in sol (martian day),
823! into solar longitude date/time, in Ls (in degrees),
824! where sol=0 is (by definition) the northern hemisphere
825!  spring equinox (where Ls=0).
826!==============================================================================
827! Notes:
828! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year,
829! "Ls" will be increased by N*360
830! Won't work as expected if sol is negative (then again,
831! why would that ever happen?)
832!==============================================================================
833
834implicit none
835
836!==============================================================================
837! Arguments:
838!==============================================================================
839real,intent(in) :: sol
840real,intent(out) :: Ls
841
842!==============================================================================
843! Local variables:
844!==============================================================================
845real year_day,peri_day,timeperi,e_elips,twopi,degrad
846data year_day /669./            ! # of sols in a martian year
847data peri_day /485.0/           
848data timeperi /1.9082314/! True anomaly at vernal equinox = 2*PI-Lsperi
849data e_elips  /0.093358/
850data twopi       /6.2831853/    ! 2.*pi
851data degrad   /57.2957795/      ! pi/180
852
853real zanom,xref,zx0,zdx,zteta,zz
854
855integer count_years
856integer iter
857
858!==============================================================================
859! 1. Compute Ls
860!==============================================================================
861
862zz=(sol-peri_day)/year_day
863zanom=twopi*(zz-nint(zz))
864xref=abs(zanom)
865
866!  The equation zx0 - e * sin (zx0) = xref, solved by Newton
867zx0=xref+e_elips*sin(xref)
868do iter=1,20 ! typically, 2 or 3 iterations are enough
869   zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0))
870   zx0=zx0+zdx
871   if(abs(zdx).le.(1.e-7)) then
872!      write(*,*)'iter:',iter,'     |zdx|:',abs(zdx)
873      exit
874   endif
875enddo
876
877if(zanom.lt.0.) zx0=-zx0
878
879zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
880Ls=zteta-timeperi
881
882if(Ls.lt.0.) then
883   Ls=Ls+twopi
884else
885   if(Ls.gt.twopi) then
886      Ls=Ls-twopi
887   endif
888endif
889
890Ls=degrad*Ls
891! Ls is now in degrees
892
893!==============================================================================
894! 1. Account for (eventual) years included in input date/time sol
895!==============================================================================
896
897count_years=0 ! initialize
898zz=sol  ! use "zz" to store (and work on) the value of sol
899do while (zz.ge.year_day)
900   count_years=count_years+1
901   zz=zz-year_day
902enddo
903
904! Add 360 degrees to Ls for every year
905if (count_years.ne.0) then
906   Ls=Ls+360.*count_years
907endif
908
909
910end subroutine sol2ls
911!************************************
912
913Subroutine interpolf(x,y,xd,yd,nd)
914!interpolation, give y = f(x) with array xd,yd known, size nd
915! Version with CONSTANT values oustide limits
916! Variable declaration
917! --------------------
918!  Arguments :
919real x,y
920real xd(*),yd(*)
921integer nd
922!  internal
923integer i,j
924real y_undefined
925                                                                               
926! run
927! ---
928      y_undefined=1.e20
929                                                                               
930      y=0.
931   if ((x.le.xd(1)).and.(x.le.xd(nd))) then
932    if (xd(1).lt.xd(nd)) y = yd(1) !y_undefined
933    if (xd(1).ge.xd(nd)) y = y_undefined ! yd(nd)
934   else if ((x.ge.xd(1)).and.(x.ge.xd(nd))) then
935    if (xd(1).lt.xd(nd)) y =  y_undefined ! yd(nd)
936    if (xd(1).ge.xd(nd)) y = y_undefined ! yd(1)
937    y = yd(nd)
938   else
939    do i=1,nd-1
940     if ( ( (x.ge.xd(i)).and.(x.lt.xd(i+1)) ).or.((x.le.xd(i)).and.(x.gt.xd(i+1))) ) then
941     y=yd(i)+(x-xd(i))*(yd(i+1)-yd(i))/(xd(i+1)-xd(i))
942     goto 99
943     end if
944    end do
945   end if
946
947!     write (*,*)' x , y' , x,y
948!     do i=1,nd
949!       write (*,*)' i, xd , yd' ,i, xd(i),yd(i)
950!     end do
951!     stop
952                                                                               
95399   continue
954                                                                               
955end subroutine interpolf
956     
957end
958
959
960
961!******************************************************************************
962!******************************************************************************
963subroutine  missing_value(nout,nvarid,valid_range,missing)
964!==============================================================================
965! Purpose:
966! Write "valid_range" and "missing_value" attributes (of a given
967! variable) to a netcdf file
968!==============================================================================
969! Remarks:
970! NetCDF file must be open
971! Variable (nvarid) ID must be set
972!==============================================================================
973
974implicit none
975
976include "netcdf.inc"  ! NetCDF definitions
977
978!==============================================================================
979! Arguments:
980!==============================================================================
981integer, intent(in) :: nout
982! nout: [netcdf] file ID #
983integer, intent(in) :: nvarid
984! varid: [netcdf] variable ID #
985real, dimension(2), intent(in) :: valid_range
986! valid_range(2): [netcdf] "valid_range" attribute (min and max)
987real, intent(in) :: missing
988! missing: [netcdf] "missing_value" attribute
989
990!==============================================================================
991! Local variables:
992!==============================================================================
993integer :: ierr
994! ierr: [netcdf] subroutine returned error code
995!      INTEGER nout,nvarid,ierr
996!      REAL missing
997!      REAL valid_range(2)
998
999! Switch to netcdf dataset define mode
1000ierr = NF_REDEF (nout)
1001if (ierr.ne.NF_NOERR) then
1002   print*,'missing_value: '
1003   print*, NF_STRERROR(ierr)
1004endif
1005
1006! Write valid_range() attribute
1007#ifdef NC_DOUBLE
1008ierr=NF_PUT_ATT_DOUBLE(nout,nvarid,'valid_range',NF_DOUBLE,2,valid_range)
1009#else
1010ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range)
1011#endif
1012
1013if (ierr.ne.NF_NOERR) then
1014   print*,'missing_value: valid_range attribution failed'
1015   print*, NF_STRERROR(ierr)
1016   !write(*,*) 'NF_NOERR', NF_NOERR
1017   !CALL abort
1018   stop ""
1019endif
1020
1021! Write "missing_value" attribute
1022#ifdef NC_DOUBLE
1023ierr= NF_PUT_ATT_DOUBLE(nout,nvarid,'missing_value',NF_DOUBLE,1,missing)
1024#else
1025ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,missing)
1026#endif
1027
1028if (ierr.NE.NF_NOERR) then
1029   print*, 'missing_value: missing value attribution failed'
1030   print*, NF_STRERROR(ierr)
1031!    WRITE(*,*) 'NF_NOERR', NF_NOERR
1032!    CALL abort
1033   stop ""
1034endif
1035
1036! End netcdf dataset define mode
1037ierr = NF_ENDDEF(nout)
1038
1039end subroutine  missing_value
1040!******************************************************************************
Note: See TracBrowser for help on using the repository browser.