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

Last change on this file since 2394 was 2394, checked in by cmathe, 4 years ago

improvedco2clouds_mod: clean

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