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

Last change on this file since 632 was 632, checked in by emillour, 13 years ago

Mars GCM utilities: minor improvements:

  • zrecast also looks for a 'diagfi1.nc' file when looking for 'phisinit' variable
  • concatnc now also includes 'aires' in its output file

EM

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