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

Last change on this file since 2308 was 2308, checked in by abierjon, 5 years ago

Mars GCM:
Following r2303 (concatenate concat files) and truly resolving Ticket #46 :
1) little bugs fixed from the previous revision ; 2) add the possibility to change the
offset of the output file time axis by asking the user the new offset and the level of
priority to put on it over the starting day stored in the first input file
(NB: this last point also implies a change in concatnc.def, but retrocompatibility
with old concatnc.def files has been ensured)

AB

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