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

Last change on this file since 3325 was 3077, checked in by emillour, 14 months ago

Mars PCM utility:
Improve concatnc's handling of XIOS output files (handle change in some
variable names, e.g. "area" or "phisfi").
EM

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