program concatnc ! ******************************************************** ! Program to concatenate data from Netcdf "diagfi" files, ! outputs from the martian GCM ! input : diagfi.nc / concat.nc / stats.nc kind of files ! author: Y. Wanherdrick ! + aps(), bps() and phisinit() are now also written ! to output file (E. Millour, september 2006) ! if available (F. Forget, october 2006) ! + ap(), bp() (F. Forget, February 2008) ! ******************************************************** implicit none include "netcdf.inc" ! NetCDF definitions character (len=50), dimension(200) :: file ! file(): input file(s) names(s) character (len=30), dimension(15) :: notconcat ! notconcat(): names of the (15) variables that won't be concatenated character (len=50), dimension(:), allocatable :: var ! var(): name(s) of variable(s) that will be concatenated character (len=50) :: tmpvar,tmpfile,title,units ! tmpvar(): used to temporarily store a variable name ! tmpfile(): used to temporarily store a file name ! title(): [netcdf] title attribute ! units(): [netcdf] units attribute character (len=100) :: filename,vartmp ! filename(): output file name ! vartmp(): temporary variable name (read from netcdf input file) !character (len=1) :: ccopy ! ccpy: 'y' or 'n' answer character (len=4) :: axis ! axis: "ls" or "sol" integer :: nid,ierr,miss ! nid: [netcdf] file ID # ! ierr: [netcdf] subroutine returned error code ! miss: [netcdf] subroutine returned error code integer :: i,j,k,inter ! for various loops integer :: varid ! varid: [netcdf] variable ID # integer :: memolatlen,memolonlen,memoaltlen ! memolatlen: # of elements of lat(), read from the first input file ! memolonlen: # of elements of lon(), read from the first input file ! memoaltlen: # of elements of alt(), read from the first input file real, dimension(:), allocatable:: lat,lon,alt,time ! lat(): array, stores latitude coordinates ! lon(): array, stores longitude coordinates ! alt(): array, stores altitude coordinates ! time(): array, stores time coordinates integer :: nbvar,nbfile,nbvarfile,ndim ! nbvar: # of variables to concatenate ! nbfile: # number of input file(s) ! nbvarfile: total # of variables in an input file ! ndim: [netcdf] # (3 or 4) of dimensions (for variables) integer :: latdim,londim,altdim,timedim ! latdim: [netcdf] "latitude" dim ID ! londim: [netcdf] "longitude" dim ID ! altdim: [netcdf] "altdim" dim ID ! timedim: [netcdf] "timedim" dim ID integer :: latvar,lonvar,altvar,timevar ! latvar: [netcdf] ID of "latitude" variable ! lonvar: [netcdf] ID of "longitude" variable ! altvar: [netcdf] ID of "altitude" variable ! timevar: [netcdf] ID of "Time" variable integer :: latlen,lonlen,altlen,timelen ! latlen: # of elements of lat() array ! lonlen: # of elements of lon() array ! altvar: # of elements of alt() array ! timelen: # of elemnets of time() array integer :: nout,latdimout,londimout,altdimout,apdimout,timedimout,timevarout ! nout: [netcdf] output file ID ! latdimout: [netcdf] output latitude (dimension) ID ! londimout: [netcdf] output longitude (dimension) ID ! altdimout: [netcdf] output altitude (dimension) ID ! timedimout: [netcdf] output time (dimension) ID ! timevarout: [netcdf] ID of output "Time" variable integer :: reptime,rep,varidout ! reptime: total length of concatenated time() arrays ! rep: # number of elements of a time() array to write to the output file ! varidout: [netcdf] variable ID # (of a variable to write to the output file) integer :: Nnotconcat,var_ok ! Nnotconcat: # of (leading)variables that won't be concatenated ! var_ok: flag (0 or 1) integer, dimension(4) :: corner,edges,dim ! corner: [netcdf] ! edges: [netcdf] ! dim: [netcdf] real, dimension(:,:,:,:), allocatable :: var3d ! var3D(,,,): 4D array to store a field real :: memotime ! memotime: (cumulative) time value, in martian days (sols) real :: missing !PARAMETER(missing=1E+20) ! missing: [netcdf] to handle "missing" values when reading/writing files real, dimension(2) :: valid_range ! valid_range(2): [netcdf] interval in which a value is considered valid !============================================================================== ! 1.1. Get input file name(s) !============================================================================== write(*,*) write(*,*) "-- Program written specifically for NetCDF 'diagfi' files from the martian GCM --" write(*,*) write(*,*) "which files do you want to use?" write(*,*) " when list ok" nbfile=0 read(*,'(a50)') tmpfile do While (len_trim(tmpfile).ne.0) nbfile=nbfile+1 file(nbfile)=tmpfile read(*,'(a50)') tmpfile enddo if(nbfile==0) then write(*,*) "no file... game over" stop "" endif !============================================================================== ! 1.2. Ask for starting day value (memotime) !============================================================================== write(*,*) !write(*,*) "Beginning day of the first specified file?" write(*,*) "Starting day of the run stored in the first input file?" write(*,*) "(e.g.: 100 if that run started at time=100 sols)" read(*,*) memotime !============================================================================== ! 1.3. Open the first input file !============================================================================== ierr = NF_OPEN(file(1),NF_NOWRITE,nid) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Pb opening file '//file(1) stop "" endif ierr=NF_INQ_NVARS(nid,nbvarfile) ! nbvarfile now set to be the (total) number of variables in file !============================================================================== ! 1.4. Ask for (output) "Time" axis type !============================================================================== write(*,*) "Warning: to read the result with grads, choose sol" write(*,*) "Warning: use ferret to read the non linear scale ls" write(*,*) "Which time axis should be given in the output? (sol/ls)" read(*,*) axis ! loop as long as axis is neither "sol" nor "ls" do while ((axis/="sol").AND.(axis/="ls")) read(*,*) axis enddo ! The following variables don't need to be concatenated notconcat(1)='Time' notconcat(2)='controle' notconcat(3)='rlonu' notconcat(4)='latitude' notconcat(5)='longitude' notconcat(6)='altitude' notconcat(7)='rlatv' notconcat(8)='aps' notconcat(9)='bps' notconcat(10)='ap' notconcat(11)='bp' notconcat(12)='cu' notconcat(13)='cv' notconcat(14)='aire' notconcat(15)='phisinit' !============================================================================== ! 1.5. Get (and check) list of variables to concatenate !============================================================================== write(*,*) Nnotconcat=0 do i=1,nbvarfile ierr=NF_INQ_VARNAME(nid,i,vartmp) ! vartmp now contains the "name" of variable of ID # i var_ok=0 do inter=1,15 if (vartmp.eq.notconcat(inter)) then var_ok=1 Nnotconcat=Nnotconcat+1 endif enddo if (var_ok.eq.0) write(*,*) trim(vartmp) enddo ! Nnotconcat: # of variables that won't be concatenated ! nbvarfile: total # of variables in file allocate(var(nbvarfile-Nnotconcat)) write(*,*) write(*,*) "which variables do you want to concatenate?" write(*,*) "all / list of (separated by s)" write(*,*) "(an empty line , i.e: just , implies end of list)" nbvar=0 read(*,'(a50)') tmpvar do while ((tmpvar/=' ').AND.(trim(tmpvar)/='all')) nbvar=nbvar+1 var(nbvar)=tmpvar read(*,'(a50)') tmpvar enddo if (tmpvar=="all") then if (axis=="ls") then ! write(*,*) "Do you want to keep the original file? (y/n)" ! read(*,*) ccopy ! if ((ccopy=="n").or.(ccopy=="N")) then ! do i=1,nbfile ! ierr=NF_CLOSE(nid) ! ierr = NF_OPEN(file(1),NF_WRITE,nid) ! call change_time_axis(nid,ierr) ! ierr=NF_CLOSE(nid) ! STOP "" ! enddo ! else nbvar=nbvarfile-Nnotconcat do j=Nnotconcat+1,nbvarfile ierr=nf_inq_varname(nid,j,var(j-Nnotconcat)) enddo ! endif ! of if ((ccopy=="n").or.(ccopy=="N")) endif ! of if (axis=="ls") ! Variables names from the file are catched nbvar=nbvarfile-Nnotconcat do i=1,nbvar ierr=nf_inq_varname(nid,i+Nnotconcat,var(i)) write(*,'(a9,1x,i2,1x,a1,1x,a50)') "variable ",i,":",var(i) enddo else if(nbvar==0) then write(*,*) "no variable... game over" stop "" endif ! of if (tmpvar=="all") ! Name of the new file !========================================================== !filename=var(1) !do i=2, nbvar ! filename=trim(adjustl(filename))//"_"//var(i) !enddo !filename=trim(adjustl(filename))//".nc" !============================================================================== ! 1.6. Get output file name !============================================================================== filename="concat.nc" !============================================================================== ! 2. Concatenate input file(s) into output file !============================================================================== reptime=0 do i=1,nbfile !============================================================================== ! 2.1. Open input file !============================================================================== if (i/=1) then write(*,*) write(*,*) "opening "//trim(file(i))//"..." ierr = NF_OPEN(file(i),NF_NOWRITE,nid) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Pb opening file '//file(i) write(*,*) NF_STRERROR(ierr) stop "" endif endif !============================================================================== ! 2.2. Read (and check) dimensions of variables from input file !============================================================================== ierr=NF_INQ_DIMID(nid,"latitude",latdim) ierr=NF_INQ_VARID(nid,"latitude",latvar) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Field is missing in file'//file(i) stop "" endif ierr=NF_INQ_DIMLEN(nid,latdim,latlen) ! write(*,*) "latlen: ",latlen ierr=NF_INQ_DIMID(nid,"longitude",londim) ierr=NF_INQ_VARID(nid,"longitude",lonvar) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Field is missing in file'//file(i) stop "" endif ierr=NF_INQ_DIMLEN(nid,londim,lonlen) ! write(*,*) "lonlen: ",lonlen ierr=NF_INQ_DIMID(nid,"altitude",altdim) ierr=NF_INQ_VARID(nid,"altitude",altvar) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Field is missing in file'//file(i) stop "" endif ierr=NF_INQ_DIMLEN(nid,altdim,altlen) ! write(*,*) "altlen: ",altlen !============================================================================== ! 2.3. Read (and check compatibility of) dimensions of ! variables from input file !============================================================================== if (i==1) then ! First call; initialize/allocate memolatlen=latlen memolonlen=lonlen memoaltlen=altlen allocate(lat(latlen)) allocate(lon(lonlen)) allocate(alt(altlen)) #ifdef NC_DOUBLE ierr = NF_GET_VAR_DOUBLE(nid,latvar,lat) ierr = NF_GET_VAR_DOUBLE(nid,lonvar,lon) ierr = NF_GET_VAR_DOUBLE(nid,altvar,alt) #else ierr = NF_GET_VAR_REAL(nid,latvar,lat) ierr = NF_GET_VAR_REAL(nid,lonvar,lon) ierr = NF_GET_VAR_REAL(nid,altvar,alt) #endif ! Initialize output file's lat,lon,alt and time dimensions call initiate (filename,lat,lon,alt,nout,& latdimout,londimout,altdimout,apdimout,timedimout,timevarout) ! Initialize output file's aps,bps,ap,bp and phisinit variables call init2(nid,lonlen,latlen,altlen,& nout,londimout,latdimout,altdimout,apdimout) else ! Not a first call, ! Check that latitude,longitude and altitude of current input file ! are identical to those of the output file if (memolatlen/=latlen) then write(*,*) "ERROR: Not the same latitude axis" stop "" else if (memolonlen/=lonlen) then write(*,*) "ERROR: Not the same longitude axis" stop "" else if (memoaltlen/=altlen) then write(*,*) "ERROR: Not the same altitude axis" stop "" endif endif ! of if (i==1) !============================================================================== ! 2.4. Handle "Time" dimension from input file !============================================================================== !============================================================================== ! 2.4.1 Read "Time" dimension from input file !============================================================================== ierr=NF_INQ_DIMID(nid,"Time",timedim) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Dimension