Ignore:
Timestamp:
Oct 17, 2013, 4:19:26 PM (11 years ago)
Author:
tnavarro
Message:

Better handling of the first date of the file : Read and write the controle field, if possible, for zrecast, localtime and concatnc + cosmetic change in lslin. This is retrocompatible with previous versions.

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/LMDZ.MARS/util/localtime.F90

    r410 r1073  
    3535integer :: varid
    3636! varid: [netcdf] variable ID #
    37 real, dimension(:), allocatable:: lat,lon,alt,time
     37real, dimension(:), allocatable:: lat,lon,alt,ctl,time
    3838! lat(): array, stores latitude coordinates
    3939! lon(): array, stores longitude coordinates
    4040! alt(): array, stores altitude coordinates
     41! ctl(): array, stores controle array
    4142! time(): array, stores time coordinates
    4243integer :: nbvar,nbvarfile,ndim
     
    4546! nbvarfile: total # of variables in an input file
    4647! ndim: [netcdf] # (3 or 4) of dimensions (for variables)
    47 integer :: latdim,londim,altdim,timedim
     48integer :: latdim,londim,altdim,ctldim,timedim
    4849! latdim: [netcdf] "latitude" dim ID
    4950! londim: [netcdf] "longitude" dim ID
    5051! altdim: [netcdf] "altdim" dim ID
     52! ctldim: [netcdf] "controle" dim ID
    5153! timedim: [netcdf] "timedim" dim ID
    52 integer :: latvar,lonvar,altvar,timevar
     54integer :: latvar,lonvar,altvar,ctlvar,timevar
    5355! latvar: [netcdf] ID of "latitude" variable
    5456! lonvar: [netcdf] ID of "longitude" variable
    5557! altvar: [netcdf] ID of "altitude" variable
     58! ctlvar: [netcdf] ID of "controle" variable
    5659! timevar: [netcdf] ID of "Time" variable
    57 integer :: latlen,lonlen,altlen,timelen,timelen_lt,timelen_tot
     60integer :: latlen,lonlen,altlen,ctllen,timelen,timelen_lt,timelen_tot
    5861integer :: ilat,ilon,ialt,it
    5962! latlen: # of elements of lat() array
    6063! lonlen: # of elements of lon() array
    6164! altvar: # of elements of alt() array
     65! ctlvar: # of elements of ctl() array
    6266! timelen: # of elemnets of time() array
    6367! timelen_tot:# =timelen or timelen+1 (if 1 more time to interpolate needed)
     
    246250!  write(*,*) "altlen: ",altlen
    247251
     252   ierr=NF_INQ_DIMID(nid,"index",ctldim)
     253   ierr=NF_INQ_VARID(nid,"controle",ctlvar)
     254   if (ierr.NE.NF_NOERR) then
     255      write(*,*) 'Field <controle> is missing in file'//file
     256      ctllen=0
     257      !stop ""
     258   else
     259      ierr=NF_INQ_DIMLEN(nid,ctldim,ctllen)
     260   endif
     261!  write(*,*) "controle: ",controle
     262
    248263!==============================================================================
    249264! 2.3. Read (and check compatibility of) dimensions of
     
    255270      allocate(lon(lonlen))
    256271      allocate(alt(altlen))
     272      allocate(ctl(ctllen))
    257273#ifdef NC_DOUBLE
    258274      ierr = NF_GET_VAR_DOUBLE(nid,latvar,lat)
    259275      ierr = NF_GET_VAR_DOUBLE(nid,lonvar,lon)
    260276      ierr = NF_GET_VAR_DOUBLE(nid,altvar,alt)
     277      if (ctllen .ne. -1) ierr = NF_GET_VAR_DOUBLE(nid,ctlvar,ctl)
    261278#else
    262279      ierr = NF_GET_VAR_REAL(nid,latvar,lat)
    263280      ierr = NF_GET_VAR_REAL(nid,lonvar,lon)
    264281      ierr = NF_GET_VAR_REAL(nid,altvar,alt)
     282      if (ctllen .ne. -1) ierr = NF_GET_VAR_REAL(nid,ctlvar,ctl)
    265283#endif
    266284!==============================================================================
     
    351369
    352370   ! Initialize output file's lat,lon,alt and time dimensions
    353       call initiate (filename,lat,lon,alt,nout,&
     371      call initiate (filename,lat,lon,alt,ctl,nout,&
    354372           latdimout,londimout,altdimout,timedimout,timevarout)
    355373   ! Initialize output file's aps,bps and phisinit variables
     
    561579
    562580!******************************************************************************
    563 Subroutine initiate (filename,lat,lon,alt,&
     581Subroutine initiate (filename,lat,lon,alt,ctl,&
    564582                     nout,latdimout,londimout,altdimout,timedimout,timevarout)
    565583!==============================================================================
     
    586604real, dimension(:), intent(in):: alt
    587605! alt(): altitude
     606real, dimension(:), intent(in):: ctl
     607! ctl(): controle
    588608integer, intent(out):: nout
    589609! nout: [netcdf] file ID
     
    604624!integer :: latdim,londim,altdim,timedim
    605625integer :: nvarid,ierr
     626integer :: ctldimout
    606627! nvarid: [netcdf] ID of a variable
    607628! ierr: [netcdf]  return error code (from called subroutines)
     
    625646ierr = NF_DEF_DIM(nout, "longitude", size(lon), londimout)
    626647ierr = NF_DEF_DIM(nout, "altitude", size(alt), altdimout)
     648if (size(ctl).ne.0) ierr = NF_DEF_DIM(nout, "index", size(ctl), ctldimout)
    627649ierr = NF_DEF_DIM(nout, "Time", NF_UNLIMITED, timedimout)
    628650
     
    664686
    665687!==============================================================================
    666 ! 4. Write "altitude" (data and attributes)
     688! 5. Write "altitude" (data and attributes)
    667689!==============================================================================
    668690
     
    688710ierr = NF_PUT_VAR_REAL (nout,nvarid,alt)
    689711#endif
     712
     713!==============================================================================
     714! 6. Write "controle" (data and attributes)
     715!==============================================================================
     716
     717if (size(ctl).ne.0) then
     718   ! Switch to netcdf define mode
     719   ierr = NF_REDEF (nout)
     720
     721   #ifdef NC_DOUBLE
     722   ierr = NF_DEF_VAR (nout,"controle",NF_DOUBLE,1,ctldimout,nvarid)
     723   #else
     724   ierr = NF_DEF_VAR (nout,"controle",NF_FLOAT,1,ctldimout,nvarid)
     725   #endif
     726
     727   ierr = NF_PUT_ATT_TEXT (nout,nvarid,"title",18,"Control parameters")
     728
     729   ! End netcdf define mode
     730   ierr = NF_ENDDEF(nout)
     731
     732   #ifdef NC_DOUBLE
     733   ierr = NF_PUT_VAR_DOUBLE (nout,nvarid,ctl)
     734   #else
     735   ierr = NF_PUT_VAR_REAL (nout,nvarid,ctl)
     736   #endif
     737endif
    690738
    691739end Subroutine initiate
Note: See TracChangeset for help on using the changeset viewer.