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

Last change on this file since 137 was 137, checked in by aslmd, 14 years ago

ATTENTION COMMIT MAJEUR NON-CONSERVATIF

  • CHANGEMENT ARBORESCENCE POUR SEPARATION CLAIRE DES COMPOSANTES et MODELES
  • UTILISATION DE LIENS SYMBOLIQUES POUR GARDER UNE BASE COMMUNE LMDZ.COMMON
  • VOIR LE FICHIER DOC/000-MODELS POUR EN SAVOIR PLUS !

EN CAS DE PROBLEME AVEC svn update IL EST CONSEILLE DE REVENIR A UN svn co

DERNIERE REVISION AVEC L'ANCIENNE ARBORESCENCE : 132

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