source: trunk/LMDZ.MARS/util/concatnc.F90 @ 2434

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

Mars GCM:
Update utilities :
concat.F90

  • rewrite and simplify the handling of time and offset so that any file can be concatenated, including files from different years or stats file
  • use modulo to shift starting sols in concat.nc to values between 0 and 669
  • add a "adls" option in addition to "sol" or "ls" to add a Ls 1D variable

while keeping "sol" as the Time axis

  • add conservation of altitude attributes (long_name,units,positive)
  • enable absence of both hybrid (aps,bps) and sigma levels

solzenangle.F90

  • improve calculation for 1 sol file (stats) to use all local time data
  • read the starting sol in the file instead of asking it to the user (except for stats file) ; keep the possibility for user to change it
  • ask mean Ls value for stats file instead of sol number
  • fix crash when using "all variable" option (ticket #66)
  • fix bug on aps,bps by using GCM_layers dimension instead of altitude

localtime.F90

  • fix crash when using "all variable" option (ticket #66)
  • fix bug on aps,bps by using GCM_layers dimension instead of altitude

For all 3 utilities :

  • handle all ndim cases for the variables (ticket #52)
  • change "title" attribute into "long_name" by default (ticket #48)
  • extend size of long_name character string (ticket #57)
  • remove the use of #ifdef NC_DOUBLE (ticket #67)

AB

  • Property svn:executable set to *
File size: 49.7 KB
Line 
1program concatnc
2
3
4! ********************************************************
5! Program to concatenate data from Netcdf "diagfi" files,
6! outputs from the martian GCM
7! input : diagfi.nc  / concat.nc / stats.nc kind of files
8! author: Y. Wanherdrick
9! + aps(), bps() and phisinit() are now also written
10!   to output file (E. Millour, september 2006)
11!   if available (F. Forget, october 2006)
12! + handle 1D data (EM, January 2007)
13! + ap(), bp()  (F. Forget, February 2008)
14! + handle the possibility that number of GCM layers (aps, bps
15!   or sigma) may be different from number of vertical levels
16!   of data (which is the case for outputs from 'zrecast')
17!   (EM, April 2010)
18! + handle absence of ap() and bp() if aps and bps are available
19!    (case of stats file) FF, November 2011
20! + read and write controle field, if available. TN, October 2013
21! + possibility to concatenate concat files with a coherent
22!   time axis in the output, and to redistribute the time axis
23!   from an offset given by the user. AB, April 2020
24! ********************************************************
25
26implicit none
27
28include "netcdf.inc" ! NetCDF definitions
29
30character (len=80), dimension(1000) :: file
31! file(): input file(s) names(s)
32character (len=30), dimension(16) :: notconcat
33! notconcat(): names of the (16) variables that won't be concatenated
34character (len=50), dimension(:), allocatable :: var
35! var(): name(s) of variable(s) that will be concatenated
36character (len=100) :: tmpvar,tmpfile,long_name,units
37! tmpvar(): used to temporarily store a variable name
38! tmpfile(): used to temporarily store a file name
39! long_name(): [netcdf] long_name attribute
40! units(): [netcdf] units attribute
41character (len=100) :: filename,vartmp
42! filename(): output file name
43! vartmp(): temporary variable name (read from netcdf input file)
44character (len=50) :: altlong_name,altunits,altpositive
45! altlong_name(): [netcdf] altitude "long_name" attribute
46! altunits(): [netcdf] altitude "units" attribute
47! altpositive(): [netcdf] altitude "positive" attribute
48character (len=4) :: axis
49! axis:  "sol" or "ls" or "adls"
50integer :: nid,ierr,miss
51! nid: [netcdf] file ID #
52! ierr: [netcdf] subroutine returned error code
53! miss: [netcdf] subroutine returned error code
54integer :: i,j,k,inter
55! for various loops
56integer :: varid
57! varid: [netcdf] variable ID #
58integer :: memolatlen=0,memolonlen=0,memoaltlen=0,memoctllen=0
59! memolatlen: # of elements of lat(), read from the first input file
60! memolonlen: # of elements of lon(), read from the first input file
61! memoaltlen: # of elements of alt(), read from the first input file
62! memoctllen: # of elements of ctl(), read from the first input file
63real, dimension(:), allocatable:: lat,lon,alt,ctl,time
64! lat(): array, stores latitude coordinates
65! lon(): array, stores longitude coordinates
66! alt(): array, stores altitude coordinates
67! ctl(): array, stores controle coordinates
68! time(): array, stores time coordinates
69integer :: nbvar,nbfile,nbvarfile,ndim
70! nbvar: # of variables to concatenate
71! nbfile: # number of input file(s)
72! nbvarfile: total # of variables in an input file
73! ndim: [netcdf] # (3 or 4) of dimensions (for variables)
74integer :: latdim,londim,altdim,ctldim,timedim
75! latdim: [netcdf] "latitude" dim ID
76! londim: [netcdf] "longitude" dim ID
77! altdim: [netcdf] "altdim" dim ID
78! ctldim: [netcdf] "ctldim" dim ID
79! timedim: [netcdf] "timedim" dim ID
80integer :: gcmlayerdim ! NetCDF dimension ID for # of layers in GCM
81integer :: latvar,lonvar,altvar,ctlvar,timevar
82! latvar: [netcdf] ID of "latitude" variable
83! lonvar: [netcdf] ID of "longitude" variable
84! altvar: [netcdf] ID of "altitude" variable
85! ctlvar: [netcdf] ID of "controle" variable
86! timevar: [netcdf] ID of "Time" variable
87integer :: latlen,lonlen,altlen,ctllen,timelen
88! latlen: # of elements of lat() array
89! lonlen: # of elements of lon() array
90! altlen: # of elements of alt() array
91! ctllen: # of elements of ctl() array
92! timelen: # of elemnets of time() array
93integer :: GCM_layers ! number of GCM atmospheric layers (may not be
94! same as altlen if file is an output of zrecast)
95integer :: nout,latdimout,londimout,altdimout,timedimout,timevarout
96! nout: [netcdf] output file ID
97! latdimout: [netcdf] output latitude (dimension) ID
98! londimout: [netcdf] output longitude (dimension) ID
99! altdimout: [netcdf] output altitude (dimension) ID
100! timedimout: [netcdf] output time (dimension) ID
101! timevarout: [netcdf] ID of output "Time" variable
102integer :: layerdimout ! NetCDF dimension ID for # of layers in GCM
103integer :: interlayerdimout ! dimension ID for # of interlayers in GCM
104integer :: reptime,rep,varidout
105! reptime: total length of concatenated time() arrays
106! rep: # number of elements of a time() array to write to the output file
107! varidout: [netcdf] variable ID # (of a variable to write to the output file)
108integer :: Nnotconcat,var_ok
109! Nnotconcat: # of (leading)variables that won't be concatenated
110! var_ok: flag (0 or 1)
111integer, dimension(4) :: corner,edges,dim
112! corner: [netcdf]
113! edges: [netcdf]
114! dim: [netcdf]
115real, dimension(:,:,:,:), allocatable :: var3d
116! var3D(,,,): 4D array to store a field
117real :: ctlsol
118! ctlsol: starting sol value of the file, stored in the variable ctl()
119real :: previous_last_output_time, output_time
120! previous_last_output_time: last time stored in concat.nc after treating a file
121logical :: firstsol
122! firstsol:  true if the first file's starting day must be used in the output
123real :: startsol
124! startsol: sol of the firstday as requested by the users (used if firstsol="y")
125real :: starttimeoffset
126! starttimeoffset: offset given by the user for the output time axis (0 if ! firstsol="n")
127real :: missing
128!PARAMETER(missing=1E+20)
129! missing: [netcdf] to handle "missing" values when reading/writing files
130real, dimension(2) :: valid_range
131! valid_range(2): [netcdf] interval in which a value is considered valid
132real :: year_day = 669.
133
134!==============================================================================
135! 1.1. Get input file name(s)
136!==============================================================================
137write(*,*)
138write(*,*) "-- Program written specifically for NetCDF 'diagfi'/'concat' files"
139write(*,*) "                     from the martian GCM --"
140write(*,*)
141write(*,*) "which files do you want to use?"
142write(*,*) "<Enter> when list ok"
143
144nbfile=0
145read(*,'(a50)') tmpfile
146do While (len_trim(tmpfile).ne.0)
147   nbfile=nbfile+1
148   file(nbfile)=tmpfile
149   read(*,'(a50)') tmpfile
150enddo
151
152if(nbfile==0) then
153   write(*,*) "no file... game over"
154   stop ""
155endif
156
157!==============================================================================
158! 1.2. Ask for starting day value (starttimeoffset)
159!==============================================================================
160
161 write(*,*) "Starting time in the concat file ?"
162 write(*,*) "   If you really wish to change it, write it now (not recommanded!)"
163 write(*,*) "   If you do not want to change it,  just write 'n' (recommanded)"
164
165read(*,*,iostat=ierr) startsol
166
167firstsol= (ierr.eq.0)
168if (firstsol) then
169   ! if nothing or a character that is not a float is read
170   write(*,*) "1st input file's starting day will serve as starting day for the outfile."
171else ! if the user gave a number
172   write(*,*) "Your input day will serve as starting day for the outfile."
173endif
174
175!==============================================================================
176! 1.3. Open the first input file
177!==============================================================================
178
179write(*,*)
180write(*,*) "opening "//trim(file(1))//"..."
181ierr = NF_OPEN(file(1),NF_NOWRITE,nid)
182if (ierr.NE.NF_NOERR) then
183   write(*,*) 'ERROR: Pb opening file '//file(1)
184   stop ""
185endif
186
187ierr=NF_INQ_NVARS(nid,nbvarfile)
188! nbvarfile now set to be the (total) number of variables in file
189
190!==============================================================================
191! 1.4. Ask for (output) "Time" axis type
192!==============================================================================
193
194write(*,*)
195write(*,*) "Which time axis should be given in the output? (sol/ls/adls)"
196write(*,*) "        (Warning: to read the result with grads, choose sol)"
197write(*,*) "  To keep sol as the time axis, and just add Ls as a variable, type adls  "
198read(*,*) axis
199! loop as long as axis is neither "sol" nor "ls"
200do while ((axis/="sol").AND.(axis/="ls").AND.(axis/="adls"))
201   read(*,*) axis
202enddo
203
204! The following variables don't need to be concatenated
205notconcat(1)='Time'
206notconcat(2)='controle'
207notconcat(3)='rlonu'
208notconcat(4)='latitude'
209notconcat(5)='longitude'
210notconcat(6)='altitude'
211notconcat(7)='rlatv'
212notconcat(8)='aps'
213notconcat(9)='bps'
214notconcat(10)='ap'
215notconcat(11)='bp'
216notconcat(12)='cu'
217notconcat(13)='cv'
218notconcat(14)='aire'
219notconcat(15)='phisinit'
220notconcat(16)='soildepth'
221
222
223!==============================================================================
224! 1.5. Get (and check) list of variables to concatenate
225!==============================================================================
226write(*,*)
227   Nnotconcat=0
228do i=1,nbvarfile
229   ierr=NF_INQ_VARNAME(nid,i,vartmp)
230   ! vartmp now contains the "name" of variable of ID # i
231   var_ok=0
232   do inter=1,size(notconcat)
233      if (vartmp.eq.notconcat(inter)) then
234         var_ok=1
235         Nnotconcat=Nnotconcat+1
236      endif
237   enddo     
238   if (var_ok.eq.0)  write(*,*) trim(vartmp)
239enddo
240
241! Nnotconcat: # of variables that won't be concatenated
242! nbvarfile: total # of variables in file
243allocate(var(nbvarfile-Nnotconcat))
244
245
246write(*,*)
247write(*,*) "which variables do you want to concatenate?"
248write(*,*) "all / list of <variables> (separated by <Enter>s)"
249write(*,*) "(an empty line , i.e: just <Enter>, implies end of list)"
250nbvar=0
251read(*,'(a50)') tmpvar
252do while ((tmpvar/=' ').AND.(trim(tmpvar)/='all'))
253   nbvar=nbvar+1
254   var(nbvar)=tmpvar
255   read(*,'(a50)') tmpvar
256enddo
257
258if (tmpvar=="all") then
259   if (axis=="ls") then
260         nbvar=nbvarfile-Nnotconcat
261         do j=Nnotconcat+1,nbvarfile
262            ierr=nf_inq_varname(nid,j,var(j-Nnotconcat))
263         enddo
264   endif ! of if (axis=="ls")
265! Variables names from the file are catched
266   nbvar=nbvarfile-Nnotconcat
267   j=1
268   do i=1,nbvarfile
269      ierr=nf_inq_varname(nid,i,vartmp)
270      var_ok=0
271      do inter=1,size(notconcat)
272        if (vartmp.eq.notconcat(inter)) then
273           var_ok=1
274        endif
275      enddo
276      if (var_ok.eq.0) then
277         if (j .gt. nbvar) then
278           write(*,*) "PROBLEM HERE !", var
279           stop
280         endif
281         var(j) = vartmp
282         write(*,'(a9,1x,i2,1x,a1,1x,a50)') "variable ",j,":",var(j)
283         j=j+1
284      endif
285   enddo
286else if(nbvar==0) then
287   write(*,*) "no variable... game over"
288   stop ""
289endif ! of if (tmpvar=="all")
290
291
292!==============================================================================
293! 1.6. Get output file name
294!==============================================================================
295filename="concat.nc"
296
297
298!==============================================================================
299! 2. Concatenate input file(s) into output file
300!==============================================================================
301
302reptime=0
303
304do i=1,nbfile
305
306!==============================================================================
307! 2.1. Open input file
308!==============================================================================
309
310   if (i/=1) then
311      write(*,*)
312      write(*,*) "opening "//trim(file(i))//"..."
313      ierr = NF_OPEN(file(i),NF_NOWRITE,nid)
314      if (ierr.NE.NF_NOERR) then
315         write(*,*) 'ERROR: Pb opening file '//file(i)
316         write(*,*) NF_STRERROR(ierr)
317         stop ""
318      endif
319   endif
320
321!==============================================================================
322! 2.2. Read (and check) dimensions of variables from input file
323!==============================================================================
324!==============================================================================
325! 2.2.1 Lat, lon, alt, GCM_layers
326!==============================================================================
327   ierr=NF_INQ_DIMID(nid,"latitude",latdim)
328   ierr=NF_INQ_VARID(nid,"latitude",latvar)
329   if (ierr.NE.NF_NOERR) then
330      write(*,*) 'ERROR: Field <latitude> is missing in file '//file(i)
331      stop "" 
332   endif
333   ierr=NF_INQ_DIMLEN(nid,latdim,latlen)
334!  write(*,*) "latlen: ",latlen
335
336   ierr=NF_INQ_DIMID(nid,"longitude",londim)
337   ierr=NF_INQ_VARID(nid,"longitude",lonvar)
338   if (ierr.NE.NF_NOERR) then
339      write(*,*) 'ERROR: Field <longitude> is missing in file '//file(i)
340      stop ""
341   endif
342   ierr=NF_INQ_DIMLEN(nid,londim,lonlen)
343!  write(*,*) "lonlen: ",lonlen
344
345   ierr=NF_INQ_DIMID(nid,"altitude",altdim)
346   ierr=NF_INQ_VARID(nid,"altitude",altvar)
347   if (ierr.NE.NF_NOERR) then
348      write(*,*) 'ERROR: Field <altitude> is missing in file '//file(i)
349      stop ""
350   endif
351   ierr=NF_INQ_DIMLEN(nid,altdim,altlen)
352!  write(*,*) "altlen: ",altlen
353
354! load size of aps() or sigma() (in case it is not altlen)
355   ! default is that GCM_layers=altlen
356   ! but for outputs of zrecast, it may be a different value
357   ierr=NF_INQ_DIMID(nid,"GCM_layers",gcmlayerdim)
358   if (ierr.ne.NF_NOERR) then
359     ! didn't find a GCM_layers dimension; therefore we have:
360     GCM_layers=altlen
361   else
362     ! load value of GCM_layers
363     ierr=NF_INQ_DIMLEN(nid,gcmlayerdim,GCM_layers)
364   endif
365!   write(*,*) "GCM_layers=",GCM_layers
366
367!==============================================================================
368! 2.2.1 Check presence and read "controle" dimension & variable from input file
369!==============================================================================
370   ierr=NF_INQ_DIMID(nid,"index",ctldim)
371   if (ierr.NE.NF_NOERR) then
372      write(*,*) 'Dimension <index> is missing in file '//file(i)
373      write(*,*) "The program continues..."
374      ctllen=0
375      !stop ""
376   else
377      ierr=NF_INQ_VARID(nid,"controle",ctlvar)
378      if (ierr.NE.NF_NOERR) then
379         write(*,*) 'Field <controle> is missing in file '//file(i)
380         write(*,*) "The program continues..."
381         ctllen=0
382         !stop ""
383      else
384         ierr=NF_INQ_DIMLEN(nid,ctldim,ctllen)
385      endif
386   endif
387
388   if (ctllen .ne. 0) then ! variable controle
389      allocate(ctl(ctllen))
390      ierr = NF_GET_VAR_REAL(nid,ctlvar,ctl)
391      ctlsol = ctl(4)
392   endif
393
394!==============================================================================
395! 2.3 Read "Time" dimension & variable from input file
396!==============================================================================
397   ierr=NF_INQ_DIMID(nid,"Time",timedim)
398   if (ierr.NE.NF_NOERR) then
399      write(*,*) 'ERROR: Dimension <Time> is missing in file '//file(i)
400      stop ""
401   endif
402   ierr=NF_INQ_VARID(nid,"Time",timevar)
403   if (ierr.NE.NF_NOERR) then
404      write(*,*) 'ERROR: Field <Time> is missing in file '//file(i)
405      stop ""
406   endif
407   ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
408!  write(*,*) "timelen: ",timelen
409
410   ! allocate time() array and fill it with values from input file
411   allocate(time(timelen))
412
413   ierr = NF_GET_VAR_REAL(nid,timevar,time)
414
415!==============================================================================
416! 2.3. Read (and check compatibility of) dimensions of
417!       variables from input file
418!==============================================================================
419
420   if (i==1) then ! First File ; initialize/allocate
421      memolatlen=latlen
422      memolonlen=lonlen
423      memoaltlen=altlen
424      memoctllen=ctllen
425      allocate(lat(latlen))
426      allocate(lon(lonlen))
427      allocate(alt(altlen))
428      ierr = NF_GET_VAR_REAL(nid,latvar,lat)
429      ierr = NF_GET_VAR_REAL(nid,lonvar,lon)
430      ierr = NF_GET_VAR_REAL(nid,altvar,alt)
431!   Get altitude attributes to handle files with any altitude type
432      ierr=nf_get_att_text(nid,altvar,'long_name',altlong_name)
433      if (ierr.ne.nf_noerr) then
434      ! if no attribute "long_name", try "title"
435        ierr=nf_get_att_text(nid,altvar,"title",altlong_name)
436      endif
437      ierr=nf_get_att_text(nid,altvar,'units',altunits)
438      ierr=nf_get_att_text(nid,altvar,'positive',altpositive)
439      write(*,*)
440     
441      if (ctllen .ne. 0) then
442       
443         if (firstsol) then
444           starttimeoffset = float(int(ctlsol)+int(time(1) + ctl(27))) - startsol
445         else ! if firstsol.eq..false.
446!          if there is no user chosen first sol, nevertheless use
447!          starttimeoffset to lower the sol number (same exact season)
448           startsol = modulo(int(int(ctlsol)+int(time(1) + ctl(27))),int(year_day))
449           starttimeoffset = float(int(ctlsol)+int(time(1) + ctl(27))) - startsol
450         endif
451
452!        artificially estimating "previous_last_output_time" for first file read:
453         previous_last_output_time= startsol
454
455         ctl(4) = 0.  ! values written in the output
456         ctl(27) = 0. ! file controle variable (resp. day_ini and hour_ini)
457         
458      else ! if no variable "controle"
459         if (firstsol) then
460           starttimeoffset = int(time(1)) - startsol
461         else ! if firstsol.eq.false
462           starttimeoffset = 0
463         endif
464
465!        artificially estimating "previous_last_output_time" for first file read:
466         previous_last_output_time=time(1) -(time(2) - time(1)) - starttimeoffset
467      endif
468      write(*,*) "Original starting day of 1st input file", previous_last_output_time + starttimeoffset
469      write(*,*) "Starting day of the concat file : ", previous_last_output_time
470
471       
472   ! Initialize output file's lat,lon,alt and time dimensions
473      write(*,*)
474      call initiate(filename,lat,lon,alt,ctl,GCM_layers,&
475       altlong_name,altunits,altpositive,&
476       nout,latdimout,londimout,altdimout,timedimout,&
477       layerdimout,interlayerdimout,timevarout,axis)
478   ! Initialize output file's aps,bps,ap,bp and phisinit variables
479      call init2(nid,lonlen,latlen,altlen,GCM_layers,&
480                nout,londimout,latdimout,altdimout,&
481                layerdimout,interlayerdimout)
482   
483   else ! Not first file to concatenate
484   ! Check that latitude,longitude and altitude of current input file
485   ! are identical to those of the output file
486   ! and that either all the files or no file contain the controle variable
487      if (memolatlen/=latlen) then
488           write(*,*) "ERROR: Not the same latitude axis"
489           stop ""
490       else if (memolonlen/=lonlen) then
491           write(*,*) "ERROR: Not the same longitude axis"
492           stop ""
493       else if (memoaltlen/=altlen) then
494           write(*,*) "ERROR: Not the same altitude axis"
495           stop ""
496       else if (memoctllen/=ctllen) then
497           write(*,*) "ERROR: Not the same controle axis"
498           stop ""
499       endif
500   endif ! of if (i==1)
501
502
503!==============================================================================
504! 2.4 Write/extend "Time" dimension/values in output file
505!==============================================================================
506
507   rep=0
508
509   do k=reptime+1,reptime+timelen
510       rep=rep+1
511       if (ctllen.ne.0) then
512          output_time=time(rep)+ctlsol-starttimeoffset
513!         correction if the file to concatenate seems "before" previous file
514!         (for instance to concatenate diagfi from the previous year at the enf of a year)
515          do while (output_time.lt.previous_last_output_time) 
516                   output_time = output_time + year_day
517          end do
518       else ! if no control, just add timestep after timestep
519          output_time= previous_last_output_time + (time(2)-time(1)) &
520         + (time(rep)-time(1))
521       end if
522       if (rep.eq.1) write(*,*) "Sol", int(output_time) 
523
524       ierr= NF_PUT_VARA_REAL(nout,timevarout,k,1,output_time)
525   end do
526!  use the last output_time value to update memotime   
527   previous_last_output_time = output_time
528
529
530!==============================================================================
531! 2.5 Read/write variables
532!==============================================================================
533
534   do j=1,nbvar ! loop on variables to read/write
535
536!==============================================================================
537! 2.5.1 Check that variable to be read existe in input file
538!==============================================================================
539
540      write(*,*) "variable ",var(j)
541      ierr=nf_inq_varid(nid,var(j),varid)
542      if (ierr.NE.NF_NOERR) then
543         write(*,*) 'ERROR: Field <',var(j),'> not found in file '//file(i)
544         stop ""
545      endif
546      ierr=nf_inq_varndims(nid,varid,ndim)
547
548      ! Check that it is a 1D, 3D or 4D variable
549      if ((ndim.ne.1).and.(ndim.ne.3).and.(ndim.ne.4)) then
550        write(*,*) "Error:",trim(var(j))," is not a time-dependent variable,"
551        write(*,*) "so it can't be concatenated."
552        write(*,*) "We'll skip that variable..."
553        CYCLE ! go directly to the next variable
554      endif
555     
556!==============================================================================
557! 2.5.2 Prepare things in order to read/write the variable
558!==============================================================================
559
560      ! build dim(),corner() and edges() arrays
561      ! and allocate var3d() array
562      if (ndim==1) then
563         allocate(var3d(timelen,1,1,1))
564         dim(1)=timedimout
565
566         ! start indexes (where data values will be written)
567         corner(1)=reptime+1
568         corner(2)=1
569         corner(3)=1
570         corner(4)=1
571
572         ! length (along dimensions) of block of data to be written
573         edges(1)=timelen
574         edges(2)=1
575         edges(3)=1
576         edges(4)=1
577
578      else if (ndim==3) then
579         allocate(var3d(lonlen,latlen,timelen,1))
580         dim(1)=londimout
581         dim(2)=latdimout
582         dim(3)=timedimout
583
584         ! start indexes (where data values will be written)
585         corner(1)=1
586         corner(2)=1
587         corner(3)=reptime+1
588         corner(4)=1
589
590         ! length (along dimensions) of block of data to be written
591         edges(1)=lonlen
592         edges(2)=latlen
593         edges(3)=timelen
594         edges(4)=1
595   
596      else if (ndim==4) then
597         allocate(var3d(lonlen,latlen,altlen,timelen))
598         dim(1)=londimout
599         dim(2)=latdimout
600         dim(3)=altdimout
601         dim(4)=timedimout
602
603         ! start indexes (where data values will be written)
604         corner(1)=1
605         corner(2)=1
606         corner(3)=1
607         corner(4)=reptime+1
608
609         ! length (along dimensions) of block of data to be written
610         edges(1)=lonlen
611         edges(2)=latlen
612         edges(3)=altlen
613         edges(4)=timelen
614      endif
615
616      if (i==1) then ! First call: write some definitions to output file
617         units="                                                    "
618         long_name="                                                    "
619         ierr=nf_get_att_text(nid,varid,"long_name",long_name)
620         if (ierr.ne.nf_noerr) then
621         ! if no attribute "long_name", try "title"
622           ierr=nf_get_att_text(nid,varid,"title",long_name)
623         endif
624         ierr=nf_get_att_text(nid,varid,"units",units)
625         call def_var(nout,var(j),long_name,units,ndim,dim,varidout,ierr)
626      else
627         ierr=NF_INQ_VARID(nout,var(j),varidout)
628      endif
629
630!==============================================================================
631! 2.5.3 Read from input file and write (append) to the output file
632!==============================================================================
633
634      ierr = NF_GET_VAR_REAL(nid,varid,var3d)
635      ierr= NF_PUT_VARA_REAL(nout,varidout,corner,edges,var3d)
636      miss=NF_GET_ATT_REAL(nid,varid,"missing_value",missing)
637      miss=NF_GET_ATT_REAL(nid,varid,"valid_range",valid_range)
638
639      if (ierr.ne.NF_NOERR) then
640         write(*,*) 'PUT_VAR ERROR: ',NF_STRERROR(ierr)
641         stop ""
642      endif
643
644! In case there is a "valid_range" attribute
645      ! Write "valid_range" and "missing_value" attributes in output file
646      if (miss.eq.NF_NOERR) then
647         call missing_value(nout,varidout,valid_range,missing)
648      endif
649
650      ! free var3d() array
651      deallocate(var3d)
652
653   enddo ! of do j=1,nbvar
654
655   ! Free time() and compute/store array length (for further concatenations)
656   deallocate(time)
657   reptime=reptime+timelen
658
659   ! Deallocate controle if needed
660   if (ctllen.ne.0) deallocate(ctl)
661
662   ! Close input file
663   ierr=nf_close(nid)
664
665enddo ! of i=1,nbfile
666
667!==============================================================================
668! 3. If required, change time axis (from sols to Ls)
669!==============================================================================
670
671if ((axis=="ls").or.(axis=="adls")) then
672   call change_time_axis(nout,ierr,axis)
673endif
674
675! Close output file
676ierr=NF_CLOSE(nout)
677
678write(*,*);write(*,*) "Program concatnc completed !"
679contains
680
681!******************************************************************************
682subroutine initiate (filename,lat,lon,alt,ctl,GCM_layers,&
683                     altlong_name,altunits,altpositive,&
684                     nout,latdimout,londimout,altdimout,timedimout,&
685                     layerdimout,interlayerdimout,timevarout,axis)
686!==============================================================================
687! Purpose:
688! Create and initialize a data file (NetCDF format)
689!==============================================================================
690! Remarks:
691! The NetCDF file (created in this subroutine) remains open
692!==============================================================================
693
694implicit none
695
696include "netcdf.inc" ! NetCDF definitions
697
698!==============================================================================
699! Arguments:
700!==============================================================================
701character (len=*), intent(in):: filename
702! filename(): the file's name
703real, dimension(:), intent(in):: lat
704! lat(): latitude
705real, dimension(:), intent(in):: lon
706! lon(): longitude
707real, dimension(:), intent(in):: alt
708! alt(): altitude
709real, dimension(:), intent(in):: ctl
710! ctl(): controle
711integer,intent(in) :: GCM_layers ! number of GCM layers
712character (len=*), intent(in) :: altlong_name
713! altlong_name(): [netcdf] altitude "long_name" attribute
714character (len=*), intent(in) :: altunits
715! altunits(): [netcdf] altitude "units" attribute
716character (len=*), intent(in) :: altpositive
717! altpositive(): [netcdf] altitude "positive" attribute
718integer, intent(out):: nout
719! nout: [netcdf] file ID
720integer, intent(out):: latdimout
721! latdimout: [netcdf] lat() (i.e.: latitude)  ID
722integer, intent(out):: londimout
723! londimout: [netcdf] lon()  ID
724integer, intent(out):: altdimout
725! altdimout: [netcdf] alt()  ID
726integer, intent(out):: timedimout
727! timedimout: [netcdf] "Time"  ID
728integer,intent(out) :: layerdimout
729! layerdimout: [netcdf] "GCM_layers" ID
730integer,intent(out) :: interlayerdimout
731! layerdimout: [netcdf] "GCM_layers+1" ID
732integer, intent(out):: timevarout
733! timevarout: [netcdf] "Time" (considered as a variable) ID
734character (len=4),intent(in) :: axis
735! axis: "ls" or "sol"
736
737!==============================================================================
738! Local variables:
739!==============================================================================
740!integer :: latdim,londim,altdim,timedim
741integer :: nvarid,ierr
742integer :: ctldimout
743! nvarid: [netcdf] ID of a variable
744! ierr: [netcdf]  return error code (from called subroutines)
745
746!==============================================================================
747! 1. Create (and open) output file
748!==============================================================================
749write(*,*) "creating "//trim(adjustl(filename))//'...'
750ierr = NF_CREATE(filename,IOR(NF_CLOBBER,NF_64BIT_OFFSET),nout)
751! NB: setting NF_CLOBBER mode means that it's OK to overwrite an existing file
752if (ierr.NE.NF_NOERR) then
753   WRITE(*,*)'ERROR: Impossible to create the file.'
754   stop ""
755endif
756
757!==============================================================================
758! 2. Define/write "dimensions" and get their IDs
759!==============================================================================
760
761ierr = NF_DEF_DIM(nout, "latitude", size(lat), latdimout)
762ierr = NF_DEF_DIM(nout, "longitude", size(lon), londimout)
763ierr = NF_DEF_DIM(nout, "altitude", size(alt), altdimout)
764if (ctllen.ne.0) ierr = NF_DEF_DIM(nout, "index", size(ctl), ctldimout)
765ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout)
766ierr = NF_DEF_DIM(nout, "GCM_layers", GCM_layers, layerdimout)
767ierr = NF_DEF_DIM(nout, "GCM_interlayers",GCM_layers+1,interlayerdimout)
768
769! End netcdf define mode
770ierr = NF_ENDDEF(nout)
771
772!==============================================================================
773! 3. Write "Time" and its attributes
774!==============================================================================
775if (axis=="ls") then
776  call def_var(nout,"Time","Solar longitude","degree",1,&
777             (/timedimout/),timevarout,ierr)
778else
779  call def_var(nout,"Time","Time","days since 0000-00-0 00:00:00",1,&
780             (/timedimout/),timevarout,ierr)
781endif
782!==============================================================================
783! 4. Write "latitude" (data and attributes)
784!==============================================================================
785
786call def_var(nout,"latitude","latitude","degrees_north",1,&
787             (/latdimout/),nvarid,ierr)
788
789ierr = NF_PUT_VAR_REAL (nout,nvarid,lat)
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)
799
800!==============================================================================
801! 5. Write "altitude" (data and attributes)
802!==============================================================================
803
804! Switch to netcdf define mode
805ierr = NF_REDEF (nout)
806
807ierr = NF_DEF_VAR (nout,"altitude",NF_FLOAT,1,altdimout,nvarid)
808
809ierr = NF_PUT_ATT_TEXT (nout,nvarid,'long_name',len_trim(adjustl(altlong_name)),adjustl(altlong_name))
810ierr = NF_PUT_ATT_TEXT (nout,nvarid,'units',len_trim(adjustl(altunits)),adjustl(altunits))
811ierr = NF_PUT_ATT_TEXT (nout,nvarid,'positive',len_trim(adjustl(altpositive)),adjustl(altpositive))
812
813! End netcdf define mode
814ierr = NF_ENDDEF(nout)
815
816ierr = NF_PUT_VAR_REAL (nout,nvarid,alt)
817
818!==============================================================================
819! 6. Write "controle" (data and attributes)
820!==============================================================================
821
822if (ctllen.ne.0) then
823   ! Switch to netcdf define mode
824   ierr = NF_REDEF (nout)
825
826   ierr = NF_DEF_VAR (nout,"controle",NF_FLOAT,1,ctldimout,nvarid)
827
828   ierr = NF_PUT_ATT_TEXT (nout,nvarid,"long_name",18,"Control parameters")
829
830   ! End netcdf define mode
831   ierr = NF_ENDDEF(nout)
832
833   ierr = NF_PUT_VAR_REAL (nout,nvarid,ctl)
834endif
835
836end Subroutine initiate
837!******************************************************************************
838subroutine init2(infid,lonlen,latlen,altlen,GCM_layers, &
839                 outfid,londimout,latdimout,altdimout, &
840                 layerdimout,interlayerdimout)
841!==============================================================================
842! Purpose:
843! Copy ap(), bp(), aps(), bps(), aire() and phisinit()
844! from input file to output file
845!==============================================================================
846! Remarks:
847! The NetCDF files must be open
848!==============================================================================
849
850implicit none
851
852include "netcdf.inc" ! NetCDF definitions
853
854!==============================================================================
855! Arguments:
856!==============================================================================
857integer, intent(in) :: infid  ! NetCDF output file ID
858integer, intent(in) :: lonlen ! # of grid points along longitude
859integer, intent(in) :: latlen ! # of grid points along latitude
860integer, intent(in) :: altlen ! # of grid points along altitude
861integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers
862integer, intent(in) :: outfid ! NetCDF output file ID
863integer, intent(in) :: londimout ! longitude dimension ID
864integer, intent(in) :: latdimout ! latitude dimension ID
865integer, intent(in) :: altdimout ! altitude dimension ID
866integer, intent(in) :: layerdimout ! GCM_layers dimension ID
867integer, intent(in) :: interlayerdimout ! GCM_layers+1 dimension ID
868!==============================================================================
869! Local variables:
870!==============================================================================
871real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates
872real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates
873real,dimension(:),allocatable :: sigma ! sigma levels
874real,dimension(:,:),allocatable :: aire ! mesh areas
875real,dimension(:,:),allocatable :: phisinit ! Ground geopotential
876integer :: apsid,bpsid,sigmaid,phisinitid
877integer :: apid,bpid
878integer :: ierr
879integer :: tmpvarid ! temporary variable ID
880logical :: area ! is "aire" available ?
881logical :: phis ! is "phisinit" available ?
882logical :: hybrid ! are "aps" and "bps" available ?
883logical :: apbp ! are "ap" and "bp" available ?
884logical :: sig ! is "sigma" available ?
885
886!==============================================================================
887! 1. Read data from input file
888!==============================================================================
889
890! hybrid coordinate aps
891!write(*,*) "aps: altlen=",altlen," GCM_layers=",GCM_layers
892allocate(aps(GCM_layers),stat=ierr)
893if (ierr.ne.0) then
894  write(*,*) "init2: failed to allocate aps!"
895  stop
896endif
897ierr=NF_INQ_VARID(infid,"aps",tmpvarid)
898if (ierr.ne.NF_NOERR) then
899  write(*,*) "Ooops. Failed to get aps ID. OK, will look for sigma coord."
900  hybrid=.false.
901else
902  ierr=NF_GET_VAR_REAL(infid,tmpvarid,aps)
903  hybrid=.true.
904  if (ierr.ne.NF_NOERR) then
905    stop "init2 Error: Failed reading aps"
906  endif
907
908  ! hybrid coordinate bps
909!  write(*,*) "bps: altlen=",altlen," GCM_layers=",GCM_layers
910  allocate(bps(GCM_layers),stat=ierr)
911  if (ierr.ne.0) then
912    write(*,*) "init2: failed to allocate bps!"
913    stop
914  endif
915  ierr=NF_INQ_VARID(infid,"bps",tmpvarid)
916  if (ierr.ne.NF_NOERR) then
917    stop "init2 Error: Failed to get bps ID."
918  endif
919  ierr=NF_GET_VAR_REAL(infid,tmpvarid,bps)
920  if (ierr.ne.NF_NOERR) then
921    stop "init2 Error: Failed reading bps"
922  endif
923endif
924
925! hybrid coordinate ap
926allocate(ap(GCM_layers+1),stat=ierr)
927if (ierr.ne.0) then
928  write(*,*) "init2: failed to allocate ap!"
929  stop
930endif
931ierr=NF_INQ_VARID(infid,"ap",tmpvarid)
932if (ierr.ne.NF_NOERR) then
933  write(*,*) "Ooops. Failed to get ap ID. OK."
934  apbp=.false.
935else
936  ierr=NF_GET_VAR_REAL(infid,tmpvarid,ap)
937  apbp=.true.
938  if (ierr.ne.NF_NOERR) then
939    stop "Error: Failed reading ap"
940  endif
941endif
942
943! hybrid coordinate bp
944allocate(bp(GCM_layers+1),stat=ierr)
945if (ierr.ne.0) then
946  write(*,*) "init2: failed to allocate bp!"
947  stop
948endif
949ierr=NF_INQ_VARID(infid,"bp",tmpvarid)
950if (ierr.ne.NF_NOERR) then
951  write(*,*) "Ooops. Failed to get bp ID. OK."
952  apbp=.false.
953else
954  ierr=NF_GET_VAR_REAL(infid,tmpvarid,bp)
955  apbp=.true.
956  if (ierr.ne.NF_NOERR) then
957    stop "Error: Failed reading bp"
958  endif
959endif
960
961! sigma levels (if any)
962if (.not.hybrid) then
963  allocate(sigma(GCM_layers),stat=ierr)
964  if (ierr.ne.0) then
965    write(*,*) "init2: failed to allocate sigma"
966    stop
967  endif
968  ierr=NF_INQ_VARID(infid,"sigma",tmpvarid)
969  if (ierr.ne.NF_NOERR) then
970    write(*,*) "Ooops. Failed to get sigma ID. OK."
971    sig=.false.
972  else
973    ierr=NF_GET_VAR_REAL(infid,tmpvarid,sigma)
974    sig=.true.
975    if (ierr.ne.NF_NOERR) then
976      stop "init2 Error: Failed reading sigma"
977    endif
978  endif
979endif ! of if (.not.hybrid)
980
981! mesh area
982allocate(aire(lonlen,latlen),stat=ierr)
983if (ierr.ne.0) then
984  write(*,*) "init2: failed to allocate aire!"
985  stop
986endif
987ierr=NF_INQ_VARID(infid,"aire",tmpvarid)
988if (ierr.ne.NF_NOERR) then
989  write(*,*)"init2 warning: Failed to get aire ID."
990  area = .false.
991else
992  ierr=NF_GET_VAR_REAL(infid,tmpvarid,aire)
993  if (ierr.ne.NF_NOERR) then
994    stop "init2 Error: Failed reading aire"
995  endif
996  area = .true.
997endif
998
999! ground geopotential phisinit
1000allocate(phisinit(lonlen,latlen),stat=ierr)
1001if (ierr.ne.0) then
1002  write(*,*) "init2: failed to allocate phisinit!"
1003  stop
1004endif
1005ierr=NF_INQ_VARID(infid,"phisinit",tmpvarid)
1006if (ierr.ne.NF_NOERR) then
1007  write(*,*)"init2 warning: Failed to get phisinit ID."
1008  phis = .false.
1009else
1010  ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
1011  if (ierr.ne.NF_NOERR) then
1012    stop "init2 Error: Failed reading phisinit"
1013  endif
1014  phis = .true.
1015endif
1016
1017!==============================================================================
1018! 2. Write
1019!==============================================================================
1020
1021!==============================================================================
1022! 2.2. Hybrid coordinates ap() , bp(), aps() and bps()
1023!==============================================================================
1024if(hybrid) then
1025! define aps
1026  call def_var(nout,"aps","hybrid pressure at midlayers"," ",1,&
1027             (/layerdimout/),apsid,ierr)
1028  if (ierr.ne.NF_NOERR) then
1029    stop "init2 Error: Failed to def_var aps"
1030  endif
1031
1032! write aps
1033  ierr=NF_PUT_VAR_REAL(outfid,apsid,aps)
1034  if (ierr.ne.NF_NOERR) then
1035    stop "init2 Error: Failed to write aps"
1036  endif
1037
1038! define bps
1039  call def_var(nout,"bps","hybrid sigma at midlayers"," ",1,&
1040             (/layerdimout/),bpsid,ierr)
1041  if (ierr.ne.NF_NOERR) then
1042    stop "init2 Error: Failed to def_var bps"
1043  endif
1044
1045! write bps
1046  ierr=NF_PUT_VAR_REAL(outfid,bpsid,bps)
1047  if (ierr.ne.NF_NOERR) then
1048    stop "init2 Error: Failed to write bps"
1049  endif
1050
1051  if (apbp) then
1052!   define ap
1053    call def_var(nout,"ap","hybrid sigma at interlayers"," ",1,&
1054             (/interlayerdimout/),apid,ierr)
1055    if (ierr.ne.NF_NOERR) then
1056      stop "Error: Failed to def_var ap"
1057    endif
1058
1059! write ap
1060    ierr=NF_PUT_VAR_REAL(outfid,apid,ap)
1061    if (ierr.ne.NF_NOERR) then
1062      stop "Error: Failed to write ap"
1063    endif
1064
1065! define bp
1066    call def_var(nout,"bp","hybrid sigma at interlayers"," ",1,&
1067             (/interlayerdimout/),bpid,ierr)
1068    if (ierr.ne.NF_NOERR) then
1069      stop "Error: Failed to def_var bp"
1070    endif
1071
1072!   write bp
1073    ierr=NF_PUT_VAR_REAL(outfid,bpid,bp)
1074    if (ierr.ne.NF_NOERR) then
1075      stop "Error: Failed to write bp"
1076    endif
1077  endif ! of if (apbp)
1078
1079else if (sig) then
1080! define sigma
1081  call def_var(nout,"sigma","sigma at midlayers"," ",1,&
1082             (/layerdimout/),sigmaid,ierr)
1083  if (ierr.ne.NF_NOERR) then
1084    stop "init2 Error: Failed to def_var sigma"
1085  endif
1086! write sigma
1087  ierr=NF_PUT_VAR_REAL(outfid,sigmaid,sigma)
1088  if (ierr.ne.NF_NOERR) then
1089    stop "init2 Error: Failed to write sigma"
1090  endif
1091endif ! of if (hybrid)
1092
1093!==============================================================================
1094! 2.2. aire() and phisinit()
1095!==============================================================================
1096
1097if (area) then
1098  ! define aire
1099  call def_var(nout,"aire","Mesh area","m2",2,&
1100           (/londimout,latdimout/),tmpvarid,ierr)
1101  if (ierr.ne.NF_NOERR) then
1102     stop "init2 Error: Failed to def_var aire"
1103  endif
1104 
1105  ! write aire
1106  ierr=NF_PUT_VAR_REAL(outfid,tmpvarid,aire)
1107  if (ierr.ne.NF_NOERR) then
1108    stop "init2 Error: Failed to write aire"
1109  endif
1110endif ! of if (area)
1111
1112IF (phis) THEN
1113
1114  !define phisinit
1115   call def_var(nout,"phisinit","Ground level geopotential"," ",2,&
1116            (/londimout,latdimout/),phisinitid,ierr)
1117   if (ierr.ne.NF_NOERR) then
1118     stop "init2 Error: Failed to def_var phisinit"
1119   endif
1120
1121  ! write phisinit
1122  ierr=NF_PUT_VAR_REAL(outfid,phisinitid,phisinit)
1123  if (ierr.ne.NF_NOERR) then
1124    stop "init2 Error: Failed to write phisinit"
1125  endif
1126
1127ENDIF ! of IF (phis)
1128
1129
1130! Cleanup
1131if (allocated(aps)) deallocate(aps)
1132if (allocated(bps)) deallocate(bps)
1133if (allocated(ap)) deallocate(ap)
1134if (allocated(bp)) deallocate(bp)
1135if (allocated(sigma)) deallocate(sigma)
1136if (allocated(phisinit)) deallocate(phisinit)
1137if (allocated(aire)) deallocate(aire)
1138
1139end subroutine init2
1140!******************************************************************************
1141subroutine def_var(nid,name,long_name,units,nbdim,dim,nvarid,ierr)
1142!==============================================================================
1143! Purpose: Write a variable (i.e: add a variable to a dataset)
1144! called "name"; along with its attributes "long_name", "units"...
1145! to a file (following the NetCDF format)
1146!==============================================================================
1147! Remarks:
1148! The NetCDF file must be open
1149!==============================================================================
1150
1151implicit none
1152
1153include "netcdf.inc" ! NetCDF definitions
1154
1155!==============================================================================
1156! Arguments:
1157!==============================================================================
1158integer, intent(in) :: nid
1159! nid: [netcdf] file ID #
1160character (len=*), intent(in) :: name
1161! name(): [netcdf] variable's name
1162character (len=*), intent(in) :: long_name
1163! long_name(): [netcdf] variable's "long_name" attribute
1164character (len=*), intent(in) :: units
1165! unit(): [netcdf] variable's "units" attribute
1166integer, intent(in) :: nbdim
1167! nbdim: number of dimensions of the variable
1168integer, dimension(nbdim), intent(in) :: dim
1169! dim(nbdim): [netcdf] dimension(s) ID(s)
1170integer, intent(out) :: nvarid
1171! nvarid: [netcdf] ID # of the variable
1172integer, intent(out) :: ierr
1173! ierr: [netcdf] subroutines returned error code
1174
1175! Switch to netcdf define mode
1176ierr=NF_REDEF(nid)
1177
1178! Insert the definition of the variable
1179ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid)
1180
1181! Write the attributes
1182ierr=NF_PUT_ATT_TEXT(nid,nvarid,"long_name",len_trim(adjustl(long_name)),adjustl(long_name))
1183ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units))
1184
1185! End netcdf define mode
1186ierr=NF_ENDDEF(nid)
1187
1188end subroutine def_var
1189!******************************************************************************
1190subroutine change_time_axis(nid,ierr,axis)
1191!==============================================================================
1192! Purpose:
1193! Read "time" variable from a dataset, convert it from "sol" (martian
1194! days) to "Ls" (solar longitude) and write the result back in the file
1195!==============================================================================
1196! Remarks:
1197! The NetCDF file must be opened before this subroutine is called
1198!==============================================================================
1199
1200implicit none
1201
1202include "netcdf.inc" ! NetCDF definitions
1203
1204!==============================================================================
1205! Arguments:
1206!==============================================================================
1207integer, intent(in) :: nid
1208! nid: [netcdf] file ID
1209integer, intent(out) :: ierr
1210! ierr: [netcdf]  return error code
1211character (len=4),intent(in) :: axis
1212! axis: "ls" or "sol"
1213
1214!==============================================================================
1215! Local variables:
1216!==============================================================================
1217integer :: nvarid
1218! nvarid: ID of the "Time" variable
1219integer :: timelen
1220! timelen: size of the arrays
1221integer :: timedim
1222! timedim: ID of the "Time" dimension
1223integer i
1224
1225real, dimension(:), allocatable :: time,ls
1226! time(): time, given in sols
1227! ls(): time, given in Ls (solar longitude)
1228
1229!==============================================================================
1230! 1. Read
1231!==============================================================================
1232
1233ierr=NF_INQ_DIMID(nid,"Time",timedim)
1234ierr=NF_INQ_VARID(nid,"Time",nvarid)
1235if (ierr.NE.NF_NOERR) then
1236   write(*,*) 'ERROR in change_time_axis: Field <Time> not found'
1237   print*, NF_STRERROR(ierr)
1238   stop ""
1239endif
1240
1241ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
1242allocate(time(timelen),ls(timelen))
1243ierr = NF_GET_VAR_REAL(nid,nvarid,time)
1244if (ierr.ne.NF_NOERR) then
1245  write(*,*) "ERROR in change_time_axis: Failed to load Time"
1246  stop
1247endif
1248
1249!==============================================================================
1250! 2. Convert sols to Ls
1251!==============================================================================
1252
1253do i=1,timelen
1254   call sol2ls(time(i),ls(i))
1255enddo
1256
1257!==============================================================================
1258! 2.1 Check if there are not jumps in Ls (rounding problems in sol2ls)
1259!==============================================================================
1260
1261do i=1,timelen-1
1262   if ((ls(i+1)-ls(i)) > 350) then
1263       write(*,*) "+ 360 deg Ls jump solved:", ls(i), ls(i+1), "at timestep", i
1264      ls(i+1) = ls(i+1) - 360
1265       write(*,*) " corrected to now be:", ls(i), ls(i+1)
1266   else if ((ls(i)-ls(i+1)) > 350) then
1267       write(*,*) "- 360 deg Ls jump solved:", ls(i), ls(i+1), "at timestep", i
1268      ls(i+1) = ls(i+1) + 360
1269       write(*,*) " corrected to now be:", ls(i), ls(i+1)
1270   endif
1271enddo
1272
1273!==============================================================================
1274! 3. Write
1275!==============================================================================
1276
1277if (axis.eq."ls") then
1278   write(*,*) "Converting the time axis from sol to Ls..."
1279else ! if (axis.eq."adls")
1280   write(*,*) "Adding Ls as a 1D time variable"
1281   call def_var(nout,"Ls","Solar Longitude","degree",1, (/timedimout/),nvarid,ierr)
1282end if
1283
1284ierr = NF_PUT_VAR_REAL(nid,nvarid,ls)
1285if (ierr.ne.NF_NOERR) then
1286   write(*,*) "ERROR in change_time_axis: Failed to write Ls"
1287   stop
1288endif
1289
1290end subroutine change_time_axis
1291!******************************************************************************
1292subroutine sol2ls(sol,Ls)
1293!==============================================================================
1294! Purpose:
1295! Convert a date/time, given in sol (martian day),
1296! into solar longitude date/time, in Ls (in degrees),
1297! where sol=0 is (by definition) the northern hemisphere
1298!  spring equinox (where Ls=0).
1299!==============================================================================
1300! Notes:
1301! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year,
1302! "Ls" will be increased by N*360
1303! Won't work as expected if sol is negative (then again,
1304! why would that ever happen?)
1305!==============================================================================
1306
1307implicit none
1308
1309!==============================================================================
1310! Arguments:
1311!==============================================================================
1312real,intent(in) :: sol
1313real,intent(out) :: Ls
1314
1315!==============================================================================
1316! Local variables:
1317!==============================================================================
1318real year_day,peri_day,timeperi,e_elips,twopi,degrad
1319data year_day /669./            ! # of sols in a martian year
1320data peri_day /485.0/           
1321data timeperi /1.9082314/
1322data e_elips  /0.093358/
1323data twopi       /6.2831853/    ! 2.*pi
1324data degrad   /57.2957795/      ! pi/180
1325
1326real zanom,xref,zx0,zdx,zteta,zz
1327
1328integer count_years
1329integer iter
1330
1331!==============================================================================
1332! 1. Compute Ls
1333!==============================================================================
1334
1335zz=(sol-peri_day)/year_day
1336zanom=twopi*(zz-nint(zz))
1337xref=abs(zanom)
1338
1339!  The equation zx0 - e * sin (zx0) = xref, solved by Newton
1340zx0=xref+e_elips*sin(xref)
1341do iter=1,20 ! typically, 2 or 3 iterations are enough
1342   zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0))
1343   zx0=zx0+zdx
1344   if(abs(zdx).le.(1.e-7)) then
1345!      write(*,*)'iter:',iter,'     |zdx|:',abs(zdx)
1346      exit
1347   endif
1348enddo
1349
1350if(zanom.lt.0.) zx0=-zx0
1351
1352zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
1353Ls=zteta-timeperi
1354
1355if(Ls.lt.0.) then
1356   Ls=Ls+twopi
1357else
1358   if(Ls.gt.twopi) then
1359      Ls=Ls-twopi
1360   endif
1361endif
1362
1363Ls=degrad*Ls
1364! Ls is now in degrees
1365
1366!==============================================================================
1367! 1. Account for (eventual) years included in input date/time sol
1368!==============================================================================
1369
1370count_years=0 ! initialize
1371zz=sol  ! use "zz" to store (and work on) the value of sol
1372do while (zz.ge.year_day)
1373   count_years=count_years+1
1374   zz=zz-year_day
1375enddo
1376
1377! Add 360 degrees to Ls for every year
1378if (count_years.ne.0) then
1379   Ls=Ls+360.*count_years
1380endif
1381
1382
1383end subroutine sol2ls
1384!******************************************************************************
1385subroutine  missing_value(nout,nvarid,valid_range,missing)
1386!==============================================================================
1387! Purpose:
1388! Write "valid_range" and "missing_value" attributes (of a given
1389! variable) to a netcdf file
1390!==============================================================================
1391! Remarks:
1392! NetCDF file must be open
1393! Variable (nvarid) ID must be set
1394!==============================================================================
1395
1396implicit none
1397
1398include "netcdf.inc"  ! NetCDF definitions
1399
1400!==============================================================================
1401! Arguments:
1402!==============================================================================
1403integer, intent(in) :: nout
1404! nout: [netcdf] file ID #
1405integer, intent(in) :: nvarid
1406! varid: [netcdf] variable ID #
1407real, dimension(2), intent(in) :: valid_range
1408! valid_range(2): [netcdf] "valid_range" attribute (min and max)
1409real, intent(in) :: missing
1410! missing: [netcdf] "missing_value" attribute
1411
1412!==============================================================================
1413! Local variables:
1414!==============================================================================
1415integer :: ierr
1416! ierr: [netcdf] subroutine returned error code
1417!      INTEGER nout,nvarid,ierr
1418!      REAL missing
1419!      REAL valid_range(2)
1420
1421! Switch to netcdf dataset define mode
1422ierr = NF_REDEF (nout)
1423if (ierr.ne.NF_NOERR) then
1424   print*,'missing_value: '
1425   print*, NF_STRERROR(ierr)
1426endif
1427
1428! Write valid_range() attribute
1429ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range)
1430
1431if (ierr.ne.NF_NOERR) then
1432   print*,'missing_value: valid_range attribution failed'
1433   print*, NF_STRERROR(ierr)
1434   !write(*,*) 'NF_NOERR', NF_NOERR
1435   !CALL abort
1436   stop ""
1437endif
1438
1439! Write "missing_value" attribute
1440ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,missing)
1441
1442if (ierr.NE.NF_NOERR) then
1443   print*, 'missing_value: missing value attribution failed'
1444   print*, NF_STRERROR(ierr)
1445!    WRITE(*,*) 'NF_NOERR', NF_NOERR
1446!    CALL abort
1447   stop ""
1448endif
1449
1450! End netcdf dataset define mode
1451ierr = NF_ENDDEF(nout)
1452
1453end subroutine  missing_value
1454!******************************************************************************
1455
1456end program concatnc
Note: See TracBrowser for help on using the repository browser.