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

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

Update utilitaire "concatnc" du GCM Mars.
EM

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