program localtime ! ---------------------------------------------------------------------------- ! Program to redistribute and interpolate the variable a the same ! local times everywhere ! input : diagfi.nc / concat.nc / stats.nc kind of files ! author: F. Forget ! ---------------------------------------------------------------------------- implicit none include "netcdf.inc" ! NetCDF definitions character (len=100) file ! file(): input file(s) names(s) character (len=30), dimension(16) :: notprocessed ! notprocessed(): names of the (16) variables that won't be processed character (len=100), dimension(:), allocatable :: var ! var(): name(s) of variable(s) that will be processed character (len=100) :: tmpvar,long_name,units ! tmpvar(): used to temporarily store a variable name ! long_name(): [netcdf] long_name 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=50) :: altlong_name,altunits,altpositive ! altlong_name(): [netcdf] altitude "long_name" attribute ! altunits(): [netcdf] altitude "units" attribute ! altpositive(): [netcdf] altitude "positive" attribute integer :: nid,ierr,miss,validr ! 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 # real, dimension(:), allocatable:: lat,lon,alt,ctl,time ! lat(): array, stores latitude coordinates ! lon(): array, stores longitude coordinates ! alt(): array, stores altitude coordinates ! ctl(): array, stores controle array ! time(): array, stores time coordinates integer :: nbvar,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,ctldim,timedim ! latdim: [netcdf] "latitude" dim ID ! londim: [netcdf] "longitude" dim ID ! altdim: [netcdf] "altdim" dim ID ! ctldim: [netcdf] "controle" dim ID ! timedim: [netcdf] "timedim" dim ID integer :: gcmlayerdim ! NetCDF dimension ID for # of layers in GCM integer :: latvar,lonvar,altvar,ctlvar,timevar ! latvar: [netcdf] ID of "latitude" variable ! lonvar: [netcdf] ID of "longitude" variable ! altvar: [netcdf] ID of "altitude" variable ! ctlvar: [netcdf] ID of "controle" variable ! timevar: [netcdf] ID of "Time" variable integer :: latlen,lonlen,altlen,ctllen,timelen,timelen_lt,timelen_tot integer :: ilat,ilon,ialt,it ! latlen: # of elements of lat() array ! lonlen: # of elements of lon() array ! altlen: # of elements of alt() array ! ctllen: # of elements of ctl() array ! timelen: # of elemnets of time() array ! timelen_tot:# =timelen or timelen+1 (if 1 more time to interpolate needed) ! timelen_lt: # of elemnets of time() array in output integer :: nhour,nsol integer :: GCM_layers ! number of GCM atmospheric layers (may not be ! same as altlen if file is an output of zrecast) integer :: nout,latdimout,londimout,altdimout,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 :: layerdimout ! NetCDF dimension ID for # of layers in GCM integer :: varidout ! varidout: [netcdf] variable ID # (of a variable to write to the output file) integer :: Nnotprocessed,var_ok ! Nnotprocessed: # 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, dimension(:,:,:,:), allocatable :: var3d_lt ! var3D_lt(,,,): 4D array to store a field in local time coordinate real, dimension(:), allocatable :: lt_gcm,lt_hour real, dimension(:), allocatable :: lt_out,lt_outc real, dimension(:), allocatable :: var_gcm 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 logical :: stats ! stats=T when reading a "stats" kind of ile !============================================================================== ! 1.1. Get input file name(s) !============================================================================== write(*,*) write(*,*) "which file do you want to use? (diagfi... stats... concat...)" read(*,'(a50)') file if (len_trim(file).eq.0) then write(*,*) "no file... game over" stop endif !============================================================================== ! 1.3. Open the input file !============================================================================== ierr = NF_OPEN(file,NF_NOWRITE,nid) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Pb opening file '//trim(file) write(*,*) NF_STRERROR(ierr) stop endif ierr=NF_INQ_NVARS(nid,nbvarfile) ! nbvarfile now set to be the (total) number of variables in file if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Pb with NF_INQ_NVARS' write(*,*) NF_STRERROR(ierr) stop endif !============================================================================== ! 1.4. List of variables that should not be processed !============================================================================== notprocessed(1)='Time' notprocessed(2)='controle' notprocessed(3)='rlonu' notprocessed(4)='latitude' notprocessed(5)='longitude' notprocessed(6)='altitude' notprocessed(7)='rlatv' notprocessed(8)='aps' notprocessed(9)='bps' notprocessed(10)='ap' notprocessed(11)='bp' notprocessed(12)='soildepth' notprocessed(13)='cu' notprocessed(14)='cv' notprocessed(15)='aire' notprocessed(16)='phisinit' !============================================================================== ! 1.5. Get (and check) list of variables to process !============================================================================== write(*,*) Nnotprocessed=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,16 if (vartmp.eq.notprocessed(inter)) then var_ok=1 Nnotprocessed=Nnotprocessed+1 endif enddo if (var_ok.eq.0) write(*,*) trim(vartmp) enddo ! Nnotprocessed: # of variables that won't be processed ! nbvarfile: total # of variables in file allocate(var(nbvarfile-Nnotprocessed),stat=ierr) if (ierr.ne.0) then write(*,*) "Error: failed allocation of var(nbvarfile-Nnotprocessed)" write(*,*) " nbvarfile=",nbvarfile write(*,*) " Nnotprocessed=",Nnotprocessed stop endif write(*,*) write(*,*) "which variables do you want to redistribute ?" write(*,*) "all / list of (separated by s)" write(*,*) "(an empty line , i.e: just , implies end of list)" nbvar=0 read(*,'(a100)') tmpvar do while ((tmpvar/=' ').AND.(trim(tmpvar)/='all')) nbvar=nbvar+1 var(nbvar)=tmpvar read(*,'(a100)') tmpvar enddo if (tmpvar=="all") then nbvar=nbvarfile-Nnotprocessed do j=Nnotprocessed+1,nbvarfile ierr=nf_inq_varname(nid,j,var(j-Nnotprocessed)) enddo ! Variables names from the file are stored in var() nbvar=nbvarfile-Nnotprocessed do i=1,nbvar ierr=nf_inq_varname(nid,i+Nnotprocessed,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") !============================================================================== ! 1.6. Get output file name !============================================================================== filename=file(1:len_trim(file)-3)//"_LT.nc" write(*,*) "Output file :"//trim(filename) !============================================================================== ! 2.1. Open input file !============================================================================== write(*,*) write(*,*) "opening "//trim(file)//"..." ierr = NF_OPEN(file,NF_NOWRITE,nid) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Pb opening file '//trim(file) write(*,*) NF_STRERROR(ierr) stop endif !============================================================================== ! 2.2. Read (and check) dimensions of variables from input file !============================================================================== ierr=NF_INQ_DIMID(nid,"latitude",latdim) if (ierr.NE.NF_NOERR) then write(*,*) 'Dimension is missing, looking for lat' ierr=NF_INQ_DIMID(nid,"lat",latdim) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Dimension and is missing in file '//trim(file) stop endif endif ierr=NF_INQ_VARID(nid,"latitude",latvar) if (ierr.NE.NF_NOERR) then write(*,*) 'Field is missing, looking for lat' ierr=NF_INQ_VARID(nid,"lat",latvar) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Field and is missing in file '//trim(file) stop endif endif ierr=NF_INQ_DIMLEN(nid,latdim,latlen) ! write(*,*) "latlen: ",latlen ierr=NF_INQ_DIMID(nid,"longitude",londim) if (ierr.NE.NF_NOERR) then write(*,*) 'Field is missing, looking for lon' ierr=NF_INQ_DIMID(nid,"lon",londim) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Dimension and is missing in file '//trim(file) stop endif endif ierr=NF_INQ_VARID(nid,"longitude",lonvar) if (ierr.NE.NF_NOERR) then write(*,*) 'Field is missing, looking for lon' ierr=NF_INQ_VARID(nid,"lon",lonvar) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Field and is missing in file '//trim(file) stop endif endif ierr=NF_INQ_DIMLEN(nid,londim,lonlen) ! write(*,*) "lonlen: ",lonlen ierr=NF_INQ_DIMID(nid,"altitude",altdim) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Dimension is missing in file '//trim(file) stop endif ierr=NF_INQ_VARID(nid,"altitude",altvar) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Field is missing in file '//trim(file) stop endif ierr=NF_INQ_DIMLEN(nid,altdim,altlen) ! write(*,*) "altlen: ",altlen ! load size of aps() or sigma() (in case it is not altlen) ! default is that GCM_layers=altlen ! but for outputs of zrecast, it may be a different value ierr=NF_INQ_DIMID(nid,"GCM_layers",gcmlayerdim) if (ierr.ne.NF_NOERR) then ! didn't find a GCM_layers dimension; therefore we have: GCM_layers=altlen else ! load value of GCM_layers ierr=NF_INQ_DIMLEN(nid,gcmlayerdim,GCM_layers) endif ! write(*,*) "GCM_layers=",GCM_layers ierr=NF_INQ_DIMID(nid,"index",ctldim) if (ierr.NE.NF_NOERR) then write(*,*) 'Dimension is missing in file '//trim(file) write(*,*) 'Looking for controle_axe' ierr=NF_INQ_DIMID(nid,"controle_axe",ctldim) if (ierr.NE.NF_NOERR) then ctllen=0 endif !stop endif ! ierr=NF_INQ_DIMID(nid,"controle_axe",ctldim) ! if (ierr.NE.NF_NOERR) then ! !stop ! endif ierr=NF_INQ_VARID(nid,"controle",ctlvar) if (ierr.NE.NF_NOERR) then write(*,*) 'Field is missing in file '//trim(file) ctllen=0 !stop else ierr=NF_INQ_DIMLEN(nid,ctldim,ctllen) endif ! write(*,*) "controle: ",controle !============================================================================== ! 2.3. Read (and check compatibility of) dimensions of ! variables from input file !============================================================================== allocate(lat(latlen),stat=ierr) if (ierr.ne.0) then write(*,*) "Error: failed to allocate lat(latlen)" stop endif allocate(lon(lonlen),stat=ierr) if (ierr.ne.0) then write(*,*) "Error: failed to allocate lon(lonlen)" stop endif allocate(alt(altlen),stat=ierr) if (ierr.ne.0) then write(*,*) "Error: failed to allocate alt(altlen)" stop endif allocate(ctl(ctllen),stat=ierr) if (ierr.ne.0) then write(*,*) "Error: failed to allocate ctl(ctllen)" stop endif ierr = NF_GET_VAR_REAL(nid,latvar,lat) if (ierr.ne.0) then write(*,*) "Error: failed to load latitude" write(*,*) NF_STRERROR(ierr) stop endif ierr = NF_GET_VAR_REAL(nid,lonvar,lon) if (ierr.ne.0) then write(*,*) "Error: failed to load longitude" write(*,*) NF_STRERROR(ierr) stop endif ierr = NF_GET_VAR_REAL(nid,altvar,alt) if (ierr.ne.0) then write(*,*) "Error: failed to load altitude" write(*,*) NF_STRERROR(ierr) stop endif ! Get altitude attributes to handle files with any altitude type ierr=nf_get_att_text(nid,altvar,'long_name',altlong_name) if (ierr.ne.nf_noerr) then ! if no attribute "long_name", try "title" ierr=nf_get_att_text(nid,timevar,"title",long_name) endif ierr=nf_get_att_text(nid,altvar,'units',altunits) ierr=nf_get_att_text(nid,altvar,'positive',altpositive) if (ctllen .ne. 0) then ierr = NF_GET_VAR_REAL(nid,ctlvar,ctl) if (ierr.ne.0) then write(*,*) "Error: failed to load controle" write(*,*) NF_STRERROR(ierr) stop endif endif ! of if (ctllen .ne. 0) !============================================================================== ! 2.4. Handle "Time" dimension from input file !============================================================================== !============================================================================== ! 2.4.0 Read "Time" dimension from input file !============================================================================== ierr=NF_INQ_DIMID(nid,"Time",timedim) if (ierr.NE.NF_NOERR) then write(*,*) 'ERROR: Dimension