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

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

Mars GCM utility concatnc : minor upgrade. Handle cases when ap and bp are not available but aps and bps are.
FF

  • Property svn:executable set to *
File size: 43.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=50), 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,NF_CLOBBER,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() and phisinit() from input file to outpout file
716!==============================================================================
717! Remarks:
718! The NetCDF files must be open
719!==============================================================================
720
721implicit none
722
723include "netcdf.inc" ! NetCDF definitions
724
725!==============================================================================
726! Arguments:
727!==============================================================================
728integer, intent(in) :: infid  ! NetCDF output file ID
729integer, intent(in) :: lonlen ! # of grid points along longitude
730integer, intent(in) :: latlen ! # of grid points along latitude
731integer, intent(in) :: altlen ! # of grid points along latitude
732integer, intent(in) :: GCM_layers ! # of GCM atmospheric layers
733integer, intent(in) :: outfid ! NetCDF output file ID
734integer, intent(in) :: londimout ! longitude dimension ID
735integer, intent(in) :: latdimout ! latitude dimension ID
736integer, intent(in) :: altdimout ! altitude dimension ID
737integer, intent(in) :: layerdimout ! GCM_layers dimension ID
738integer, intent(in) :: interlayerdimout ! GCM_layers+1 dimension ID
739!==============================================================================
740! Local variables:
741!==============================================================================
742real,dimension(:),allocatable :: aps,bps ! hybrid vertical coordinates
743real,dimension(:),allocatable :: ap,bp ! hybrid vertical coordinates
744real,dimension(:),allocatable :: sigma ! sigma levels
745real,dimension(:,:),allocatable :: phisinit ! Ground geopotential
746integer :: apsid,bpsid,sigmaid,phisinitid
747integer :: apid,bpid
748integer :: ierr
749integer :: tmpvarid ! temporary variable ID
750logical :: phis ! is "phisinit" available ?
751logical :: hybrid ! are "aps" and "bps" available ?
752logical :: apbp ! are "ap" and "bp" available ?
753
754!==============================================================================
755! 1. Read data from input file
756!==============================================================================
757
758! hybrid coordinate aps
759!write(*,*) "aps: altlen=",altlen," GCM_layers=",GCM_layers
760allocate(aps(GCM_layers),stat=ierr)
761if (ierr.ne.0) then
762  write(*,*) "init2: failed to allocate aps!"
763  stop
764endif
765ierr=NF_INQ_VARID(infid,"aps",tmpvarid)
766if (ierr.ne.NF_NOERR) then
767  write(*,*) "Ooops. Failed to get aps ID. OK, will look for sigma coord."
768  hybrid=.false.
769else
770  ierr=NF_GET_VAR_REAL(infid,tmpvarid,aps)
771  hybrid=.true.
772  if (ierr.ne.NF_NOERR) then
773    stop "init2 Error: Failed reading aps"
774  endif
775
776  ! hybrid coordinate bps
777!  write(*,*) "bps: altlen=",altlen," GCM_layers=",GCM_layers
778  allocate(bps(GCM_layers),stat=ierr)
779  if (ierr.ne.0) then
780    write(*,*) "init2: failed to allocate bps!"
781    stop
782  endif
783  ierr=NF_INQ_VARID(infid,"bps",tmpvarid)
784  if (ierr.ne.NF_NOERR) then
785    stop "init2 Error: Failed to get bps ID."
786  endif
787  ierr=NF_GET_VAR_REAL(infid,tmpvarid,bps)
788  if (ierr.ne.NF_NOERR) then
789    stop "init2 Error: Failed reading bps"
790  endif
791endif
792
793! hybrid coordinate ap
794allocate(ap(GCM_layers+1),stat=ierr)
795if (ierr.ne.0) then
796  write(*,*) "init2: failed to allocate ap!"
797  stop
798else
799  ierr=NF_INQ_VARID(infid,"ap",tmpvarid)
800  if (ierr.ne.NF_NOERR) then
801    write(*,*) "Ooops. Failed to get ap ID. OK."
802    apbp=.false.
803  else
804    ierr=NF_GET_VAR_REAL(infid,tmpvarid,ap)
805    apbp=.true.
806    if (ierr.ne.NF_NOERR) then
807      stop "Error: Failed reading ap"
808    endif
809  endif
810endif
811
812! hybrid coordinate bp
813allocate(bp(GCM_layers+1),stat=ierr)
814if (ierr.ne.0) then
815  write(*,*) "init2: failed to allocate bp!"
816  stop
817else
818  ierr=NF_INQ_VARID(infid,"bp",tmpvarid)
819  if (ierr.ne.NF_NOERR) then
820    write(*,*) "Ooops. Failed to get bp ID. OK."
821    apbp=.false.
822  else
823    ierr=NF_GET_VAR_REAL(infid,tmpvarid,bp)
824    apbp=.true.
825    if (ierr.ne.NF_NOERR) then
826      stop "Error: Failed reading bp"
827    endif
828  endif
829endif
830
831! sigma levels (if any)
832if (.not.hybrid) then
833  allocate(sigma(GCM_layers),stat=ierr)
834  if (ierr.ne.0) then
835    write(*,*) "init2: failed to allocate sigma"
836    stop
837  endif
838  ierr=NF_INQ_VARID(infid,"sigma",tmpvarid)
839  ierr=NF_GET_VAR_REAL(infid,tmpvarid,sigma)
840  if (ierr.ne.NF_NOERR) then
841    stop "init2 Error: Failed reading sigma"
842  endif
843endif ! of if (.not.hybrid)
844
845! ground geopotential phisinit
846allocate(phisinit(lonlen,latlen),stat=ierr)
847if (ierr.ne.0) then
848  write(*,*) "init2: failed to allocate phisinit!"
849  stop
850endif
851ierr=NF_INQ_VARID(infid,"phisinit",tmpvarid)
852if (ierr.ne.NF_NOERR) then
853  write(*,*)"init2 warning: Failed to get phisinit ID."
854  phis = .false.
855else
856  ierr=NF_GET_VAR_REAL(infid,tmpvarid,phisinit)
857  if (ierr.ne.NF_NOERR) then
858    stop "init2 Error: Failed reading phisinit"
859  endif
860  phis = .true.
861endif
862
863!==============================================================================
864! 2. Write
865!==============================================================================
866
867!==============================================================================
868! 2.2. Hybrid coordinates ap() , bp(), aps() and bps()
869!==============================================================================
870if(hybrid) then
871! define aps
872  call def_var(nout,"aps","hybrid pressure at midlayers"," ",1,&
873             (/layerdimout/),apsid,ierr)
874  if (ierr.ne.NF_NOERR) then
875    stop "init2 Error: Failed to def_var aps"
876  endif
877
878! write aps
879#ifdef NC_DOUBLE
880  ierr=NF_PUT_VAR_DOUBLE(outfid,apsid,aps)
881#else
882  ierr=NF_PUT_VAR_REAL(outfid,apsid,aps)
883#endif
884  if (ierr.ne.NF_NOERR) then
885    stop "init2 Error: Failed to write aps"
886  endif
887
888! define bps
889  call def_var(nout,"bps","hybrid sigma at midlayers"," ",1,&
890             (/layerdimout/),bpsid,ierr)
891  if (ierr.ne.NF_NOERR) then
892    stop "init2 Error: Failed to def_var bps"
893  endif
894
895! write bps
896#ifdef NC_DOUBLE
897  ierr=NF_PUT_VAR_DOUBLE(outfid,bpsid,bps)
898#else
899  ierr=NF_PUT_VAR_REAL(outfid,bpsid,bps)
900#endif
901  if (ierr.ne.NF_NOERR) then
902    stop "init2 Error: Failed to write bps"
903  endif
904
905  if (apbp) then
906!   define ap
907    call def_var(nout,"ap","hybrid sigma at interlayers"," ",1,&
908             (/interlayerdimout/),apid,ierr)
909    if (ierr.ne.NF_NOERR) then
910      stop "Error: Failed to def_var ap"
911    endif
912
913! write ap
914#ifdef NC_DOUBLE
915    ierr=NF_PUT_VAR_DOUBLE(outfid,apid,ap)
916#else
917    ierr=NF_PUT_VAR_REAL(outfid,apid,ap)
918#endif
919    if (ierr.ne.NF_NOERR) then
920      stop "Error: Failed to write ap"
921    endif
922
923! define bp
924    call def_var(nout,"bp","hybrid sigma at interlayers"," ",1,&
925             (/interlayerdimout/),bpid,ierr)
926    if (ierr.ne.NF_NOERR) then
927      stop "Error: Failed to def_var bp"
928    endif
929
930!   write bp
931#ifdef NC_DOUBLE
932    ierr=NF_PUT_VAR_DOUBLE(outfid,bpid,bp)
933#else
934    ierr=NF_PUT_VAR_REAL(outfid,bpid,bp)
935#endif
936    if (ierr.ne.NF_NOERR) then
937      stop "Error: Failed to write bp"
938    endif
939  endif ! of if (apbp)
940
941else
942! define sigma
943  call def_var(nout,"sigma","sigma at midlayers"," ",1,&
944             (/layerdimout/),sigmaid,ierr)
945  if (ierr.ne.NF_NOERR) then
946    stop "init2 Error: Failed to def_var sigma"
947  endif
948! write sigma
949#ifdef NC_DOUBLE
950  ierr=NF_PUT_VAR_DOUBLE(outfid,sigmaid,sigma)
951#else
952  ierr=NF_PUT_VAR_REAL(outfid,sigmaid,sigma)
953#endif
954  if (ierr.ne.NF_NOERR) then
955    stop "init2 Error: Failed to write sigma"
956  endif
957endif ! of if (hybrid)
958
959!==============================================================================
960! 2.2. phisinit()
961!==============================================================================
962
963IF (phis) THEN
964
965  !define phisinit
966   call def_var(nout,"phisinit","Ground level geopotential"," ",2,&
967            (/londimout,latdimout/),phisinitid,ierr)
968   if (ierr.ne.NF_NOERR) then
969     stop "init2 Error: Failed to def_var phisinit"
970   endif
971
972  ! write phisinit
973#ifdef NC_DOUBLE
974  ierr=NF_PUT_VAR_DOUBLE(outfid,phisinitid,phisinit)
975#else
976  ierr=NF_PUT_VAR_REAL(outfid,phisinitid,phisinit)
977#endif
978  if (ierr.ne.NF_NOERR) then
979    stop "init2 Error: Failed to write phisinit"
980  endif
981
982ENDIF ! of IF (phis)
983
984
985! Cleanup
986if (allocated(aps)) deallocate(aps)
987if (allocated(bps)) deallocate(bps)
988if (allocated(ap)) deallocate(ap)
989if (allocated(bp)) deallocate(bp)
990if (allocated(sigma)) deallocate(sigma)
991if (allocated(phisinit)) deallocate(phisinit)
992
993end subroutine init2
994!******************************************************************************
995subroutine def_var(nid,name,title,units,nbdim,dim,nvarid,ierr)
996!==============================================================================
997! Purpose: Write a variable (i.e: add a variable to a dataset)
998! called "name"; along with its attributes "title", "units"...
999! to a file (following the NetCDF format)
1000!==============================================================================
1001! Remarks:
1002! The NetCDF file must be open
1003!==============================================================================
1004
1005implicit none
1006
1007include "netcdf.inc" ! NetCDF definitions
1008
1009!==============================================================================
1010! Arguments:
1011!==============================================================================
1012integer, intent(in) :: nid
1013! nid: [netcdf] file ID #
1014character (len=*), intent(in) :: name
1015! name(): [netcdf] variable's name
1016character (len=*), intent(in) :: title
1017! title(): [netcdf] variable's "title" attribute
1018character (len=*), intent(in) :: units
1019! unit(): [netcdf] variable's "units" attribute
1020integer, intent(in) :: nbdim
1021! nbdim: number of dimensions of the variable
1022integer, dimension(nbdim), intent(in) :: dim
1023! dim(nbdim): [netcdf] dimension(s) ID(s)
1024integer, intent(out) :: nvarid
1025! nvarid: [netcdf] ID # of the variable
1026integer, intent(out) :: ierr
1027! ierr: [netcdf] subroutines returned error code
1028
1029! Switch to netcdf define mode
1030ierr=NF_REDEF(nid)
1031
1032! Insert the definition of the variable
1033#ifdef NC_DOUBLE
1034ierr=NF_DEF_VAR(nid,adjustl(name),NF_DOUBLE,nbdim,dim,nvarid)
1035#else
1036ierr=NF_DEF_VAR(nid,adjustl(name),NF_FLOAT,nbdim,dim,nvarid)
1037#endif
1038
1039! Write the attributes
1040ierr=NF_PUT_ATT_TEXT(nid,nvarid,"title",len_trim(adjustl(title)),adjustl(title))
1041ierr=NF_PUT_ATT_TEXT(nid,nvarid,"units",len_trim(adjustl(units)),adjustl(units))
1042
1043! End netcdf define mode
1044ierr=NF_ENDDEF(nid)
1045
1046end subroutine def_var
1047!******************************************************************************
1048subroutine change_time_axis(nid,ierr)
1049!==============================================================================
1050! Purpose:
1051! Read "time" variable from a dataset, convert it from "sol" (martian
1052! days) to "Ls" (solar longitude) and write the result back in the file
1053!==============================================================================
1054! Remarks:
1055! The NetCDF file must be opened before this subroutine is called
1056!==============================================================================
1057
1058implicit none
1059
1060include "netcdf.inc" ! NetCDF definitions
1061
1062!==============================================================================
1063! Arguments:
1064!==============================================================================
1065integer, intent(in) :: nid
1066! nid: [netcdf] file ID
1067integer, intent(out) :: ierr
1068! ierr: [netcdf]  return error code
1069
1070!==============================================================================
1071! Local variables:
1072!==============================================================================
1073integer :: nvarid
1074! nvarid: ID of the "Time" variable
1075integer :: timelen
1076! timelen: size of the arrays
1077integer :: timedim
1078! timedim: ID of the "Time" dimension
1079integer i
1080
1081real, dimension(:), allocatable :: time,ls
1082! time(): time, given in sols
1083! ls(): time, given in Ls (solar longitude)
1084
1085!==============================================================================
1086! 1. Read
1087!==============================================================================
1088
1089ierr=NF_INQ_DIMID(nid,"Time",timedim)
1090ierr=NF_INQ_VARID(nid,"Time",nvarid)
1091if (ierr.NE.NF_NOERR) then
1092   write(*,*) 'ERROR in change_time_axis: Field <Time> not found'
1093   print*, NF_STRERROR(ierr)
1094   stop ""
1095endif
1096
1097ierr=NF_INQ_DIMLEN(nid,timedim,timelen)
1098allocate(time(timelen),ls(timelen))
1099#ifdef NC_DOUBLE
1100ierr = NF_GET_VAR_DOUBLE(nid,nvarid,time)
1101#else
1102ierr = NF_GET_VAR_REAL(nid,nvarid,time)
1103#endif
1104if (ierr.ne.NF_NOERR) then
1105  write(*,*) "ERROR in change_time_axis: Failed to load Time"
1106  stop
1107endif
1108
1109!==============================================================================
1110! 2. Convert sols to Ls
1111!==============================================================================
1112
1113do i=1,timelen
1114   call sol2ls(time(i),ls(i))
1115enddo
1116
1117!==============================================================================
1118! 2.1 Check if there are not jumps in Ls (rounding problems in sol2ls)
1119!==============================================================================
1120
1121do i=1,timelen-1
1122   if ((ls(i+1)-ls(i)) > 350) then
1123       write(*,*) "+ 360° Ls jump solved:", ls(i), ls(i+1), "at timestep", i
1124      ls(i+1) = ls(i+1) - 360
1125       write(*,*) " corrected to now be:", ls(i), ls(i+1)
1126   else if ((ls(i)-ls(i+1)) > 350) then
1127       write(*,*) "- 360° Ls jump solved:", ls(i), ls(i+1), "at timestep", i
1128      ls(i+1) = ls(i+1) + 360
1129       write(*,*) " corrected to now be:", ls(i), ls(i+1)
1130   endif
1131enddo
1132
1133!==============================================================================
1134! 3. Write
1135!==============================================================================
1136
1137#ifdef NC_DOUBLE
1138ierr = NF_PUT_VAR_DOUBLE(nid,nvarid,ls)
1139#else
1140ierr = NF_PUT_VAR_REAL(nid,nvarid,ls)
1141#endif
1142if (ierr.ne.NF_NOERR) then
1143  write(*,*) "ERROR in change_time_axis: Failed to write Ls"
1144  stop
1145endif
1146
1147end subroutine change_time_axis
1148!******************************************************************************
1149subroutine sol2ls(sol,Ls)
1150!==============================================================================
1151! Purpose:
1152! Convert a date/time, given in sol (martian day),
1153! into solar longitude date/time, in Ls (in degrees),
1154! where sol=0 is (by definition) the northern hemisphere
1155!  spring equinox (where Ls=0).
1156!==============================================================================
1157! Notes:
1158! Even though "Ls" is cyclic, if "sol" is greater than N (martian) year,
1159! "Ls" will be increased by N*360
1160! Won't work as expected if sol is negative (then again,
1161! why would that ever happen?)
1162!==============================================================================
1163
1164implicit none
1165
1166!==============================================================================
1167! Arguments:
1168!==============================================================================
1169real,intent(in) :: sol
1170real,intent(out) :: Ls
1171
1172!==============================================================================
1173! Local variables:
1174!==============================================================================
1175real year_day,peri_day,timeperi,e_elips,twopi,degrad
1176data year_day /669./            ! # of sols in a martian year
1177data peri_day /485.0/           
1178data timeperi /1.9082314/
1179data e_elips  /0.093358/
1180data twopi       /6.2831853/    ! 2.*pi
1181data degrad   /57.2957795/      ! pi/180
1182
1183real zanom,xref,zx0,zdx,zteta,zz
1184
1185integer count_years
1186integer iter
1187
1188!==============================================================================
1189! 1. Compute Ls
1190!==============================================================================
1191
1192zz=(sol-peri_day)/year_day
1193zanom=twopi*(zz-nint(zz))
1194xref=abs(zanom)
1195
1196!  The equation zx0 - e * sin (zx0) = xref, solved by Newton
1197zx0=xref+e_elips*sin(xref)
1198do iter=1,20 ! typically, 2 or 3 iterations are enough
1199   zdx=-(zx0-e_elips*sin(zx0)-xref)/(1.-e_elips*cos(zx0))
1200   zx0=zx0+zdx
1201   if(abs(zdx).le.(1.e-7)) then
1202!      write(*,*)'iter:',iter,'     |zdx|:',abs(zdx)
1203      exit
1204   endif
1205enddo
1206
1207if(zanom.lt.0.) zx0=-zx0
1208
1209zteta=2.*atan(sqrt((1.+e_elips)/(1.-e_elips))*tan(zx0/2.))
1210Ls=zteta-timeperi
1211
1212if(Ls.lt.0.) then
1213   Ls=Ls+twopi
1214else
1215   if(Ls.gt.twopi) then
1216      Ls=Ls-twopi
1217   endif
1218endif
1219
1220Ls=degrad*Ls
1221! Ls is now in degrees
1222
1223!==============================================================================
1224! 1. Account for (eventual) years included in input date/time sol
1225!==============================================================================
1226
1227count_years=0 ! initialize
1228zz=sol  ! use "zz" to store (and work on) the value of sol
1229do while (zz.ge.year_day)
1230   count_years=count_years+1
1231   zz=zz-year_day
1232enddo
1233
1234! Add 360 degrees to Ls for every year
1235if (count_years.ne.0) then
1236   Ls=Ls+360.*count_years
1237endif
1238
1239
1240end subroutine sol2ls
1241!******************************************************************************
1242subroutine  missing_value(nout,nvarid,valid_range,missing)
1243!==============================================================================
1244! Purpose:
1245! Write "valid_range" and "missing_value" attributes (of a given
1246! variable) to a netcdf file
1247!==============================================================================
1248! Remarks:
1249! NetCDF file must be open
1250! Variable (nvarid) ID must be set
1251!==============================================================================
1252
1253implicit none
1254
1255include "netcdf.inc"  ! NetCDF definitions
1256
1257!==============================================================================
1258! Arguments:
1259!==============================================================================
1260integer, intent(in) :: nout
1261! nout: [netcdf] file ID #
1262integer, intent(in) :: nvarid
1263! varid: [netcdf] variable ID #
1264real, dimension(2), intent(in) :: valid_range
1265! valid_range(2): [netcdf] "valid_range" attribute (min and max)
1266real, intent(in) :: missing
1267! missing: [netcdf] "missing_value" attribute
1268
1269!==============================================================================
1270! Local variables:
1271!==============================================================================
1272integer :: ierr
1273! ierr: [netcdf] subroutine returned error code
1274!      INTEGER nout,nvarid,ierr
1275!      REAL missing
1276!      REAL valid_range(2)
1277
1278! Switch to netcdf dataset define mode
1279ierr = NF_REDEF (nout)
1280if (ierr.ne.NF_NOERR) then
1281   print*,'missing_value: '
1282   print*, NF_STRERROR(ierr)
1283endif
1284
1285! Write valid_range() attribute
1286#ifdef NC_DOUBLE
1287ierr=NF_PUT_ATT_DOUBLE(nout,nvarid,'valid_range',NF_DOUBLE,2,valid_range)
1288#else
1289ierr=NF_PUT_ATT_REAL(nout,nvarid,'valid_range',NF_FLOAT,2,valid_range)
1290#endif
1291
1292if (ierr.ne.NF_NOERR) then
1293   print*,'missing_value: valid_range attribution failed'
1294   print*, NF_STRERROR(ierr)
1295   !write(*,*) 'NF_NOERR', NF_NOERR
1296   !CALL abort
1297   stop ""
1298endif
1299
1300! Write "missing_value" attribute
1301#ifdef NC_DOUBLE
1302ierr= NF_PUT_ATT_DOUBLE(nout,nvarid,'missing_value',NF_DOUBLE,1,missing)
1303#else
1304ierr= NF_PUT_ATT_REAL(nout,nvarid,'missing_value',NF_FLOAT,1,missing)
1305#endif
1306
1307if (ierr.NE.NF_NOERR) then
1308   print*, 'missing_value: missing value attribution failed'
1309   print*, NF_STRERROR(ierr)
1310!    WRITE(*,*) 'NF_NOERR', NF_NOERR
1311!    CALL abort
1312   stop ""
1313endif
1314
1315! End netcdf dataset define mode
1316ierr = NF_ENDDEF(nout)
1317
1318end subroutine  missing_value
1319!******************************************************************************
1320
1321end program concatnc
Note: See TracBrowser for help on using the repository browser.