CFIO_SPutVar Subroutine

public subroutine CFIO_SPutVar(fid, vname, yyyymmdd, hhmmss, im, jm, kbeg, kount, grid, rc)

CFIO_SPutVar – Write a variable to the file

This routine is used to write a variable to an open CFIO stream. Multiple vertical levels can be written at one time provided they are contiguous in memory. Date and time must be consistent with the time increment and the starting date/time as defined in CFIO_Create. Times must fall on minute boundaries to allow GrADS to work. Error checking is done for dimensions that are out of bounds.

History

  • 1997.10.13 da Silva/Lucchesi Initial interface design.
  • 1998.02.10 Lucchesi Added support for applications running with 64-bit reals.
  • 1998.03.30 Lucchesi Documentation expanded. Clean-up of code.
  • 1998.07.02 Lucchesi Replaced vid with vname in argument list & made related mods to code.
  • 1998.09.24 Lucchesi Changed err codes, removed DIM_CHECK if-def
  • 1998.10.27 Lucchesi Added support for packing and range checks
  • 1998.12.15 Lucchesi Added support for skipping times (allTimes)
  • 1999.01.04 Lucchesi Fixed bug in skipping times (allTimes)/also changed variable initialization.
  • 1999.07.13 Lucchesi Changes for REAL or INT time dimension

Arguments

Type IntentOptional Attributes Name
integer :: fid

File handle

character(len=*) :: vname

Variable name

integer :: yyyymmdd

Year-month-day, e.g., 19971003

integer :: hhmmss

Hour-minute-second, e.g., 120000

integer :: im

size of longitudinal dimension

integer :: jm

size of latitudinal dimension

integer :: kbeg

first level to write; if 2-D grid use kbeg = 0.

integer :: kount

number of levels to write

real :: grid(im,kount)

Gridded data to write at this time

integer :: rc

Error return code: rc = 0 all is well rc = -2 time is inconsistent with increment rc = -3 number of levels is incompatible with file rc = -4 im is incompatible with file rc = -5 jm is incompatible with file rc = -6 time must fall on a minute boundary rc = -7 error in diffdate rc = -12 error determining default precision rc = -13 error determining variable type rc = -15 data outside of valid range rc = -16 data outside of packing range rc = -17 data outside of pack and valid range

NetCDF Errors


rc = -38 error from NF90_PUT_VAR (dimension variable) rc = -40 error from NF90_INQ_VARID rc = -41 error from NF90_INQ_DIMID or NF90_INQUIRE_DIMENSION (lat or lon) rc = -42 error from NF90_INQ_DIMID or NF90_INQUIRE_DIMENSION (lev) rc = -43 error from NF90_INQ_VARID (time variable) rc = -44 error from NF90_GET_ATT (time attribute) rc = -45 error from NF90_PUT_VAR rc = -46 error from NF90_GET_VAR rc = -52 error from NF90_INQUIRE_VARIABLE rc = -53 error from NF90_GET_ATT


Calls

proc~~cfio_sputvar~~CallsGraph proc~cfio_sputvar CFIO_SPutVar nf90_get_att nf90_get_att proc~cfio_sputvar->nf90_get_att nf90_get_var nf90_get_var proc~cfio_sputvar->nf90_get_var nf90_inq_dimid nf90_inq_dimid proc~cfio_sputvar->nf90_inq_dimid nf90_inq_varid nf90_inq_varid proc~cfio_sputvar->nf90_inq_varid nf90_inquire_dimension nf90_inquire_dimension proc~cfio_sputvar->nf90_inquire_dimension nf90_inquire_variable nf90_inquire_variable proc~cfio_sputvar->nf90_inquire_variable nf90_put_var nf90_put_var proc~cfio_sputvar->nf90_put_var proc~cfio_parseinttime CFIO_parseIntTime proc~cfio_sputvar->proc~cfio_parseinttime proc~diffdate DiffDate proc~cfio_sputvar->proc~diffdate proc~err err proc~cfio_sputvar->proc~err proc~diffdate->proc~cfio_parseinttime proc~julday julday proc~diffdate->proc~julday

Source Code

      subroutine CFIO_SPutVar ( fid, vname, yyyymmdd, hhmmss, &
                               im, jm, kbeg, kount, grid,     &
                               rc )
!
! !USES:

      Implicit NONE
!
! !INPUT PARAMETERS:
!
      integer        fid                 !! File handle
      character(len=*)  vname               !! Variable name
      integer        yyyymmdd            !! Year-month-day, e.g., 19971003
      integer        hhmmss              !! Hour-minute-second, e.g., 120000

      integer         im                 !! size of longitudinal dimension
      integer         jm                 !! size of latitudinal  dimension
      integer         kbeg               !! first level to write; if 2-D grid
                                         !!   use kbeg = 0.
      integer         kount              !! number of levels to write
      real            grid(im,kount)  !! Gridded data to write at this time


! !OUTPUT PARAMETERS:

      integer        rc  !! Error return code:
                         !!  rc =  0  all is well
                         !!  rc = -2  time is inconsistent with increment
                         !! rc = -3  number of levels is incompatible with file
                         !!  rc = -4  im is incompatible with file
                         !!  rc = -5  jm is incompatible with file
                         !!  rc = -6  time must fall on a minute boundary
                         !!  rc = -7  error in diffdate
                         !!  rc = -12  error determining default precision
                         !!  rc = -13  error determining variable type
                         !!  rc = -15  data outside of valid range
                         !!  rc = -16  data outside of packing range
                         !!  rc = -17  data outside of pack and valid range
                         !!
                         !!  NetCDF Errors
                         !!  -------------
                         !!  rc = -38  error from NF90_PUT_VAR (dimension variable)
                         !!  rc = -40  error from NF90_INQ_VARID
                         !!  rc = -41  error from NF90_INQ_DIMID or NF90_INQUIRE_DIMENSION (lat or lon)
                         !!  rc = -42  error from NF90_INQ_DIMID or NF90_INQUIRE_DIMENSION (lev)
                         !!  rc = -43  error from NF90_INQ_VARID (time variable)
                         !!  rc = -44  error from NF90_GET_ATT (time attribute)
                         !!  rc = -45  error from NF90_PUT_VAR
                         !!  rc = -46  error from NF90_GET_VAR
                         !!  rc = -52  error from NF90_INQUIRE_VARIABLE
                         !!  rc = -53  error from NF90_GET_ATT

!
!-------------------------------------------------------------------------

      integer timeid, timeDimId, dimSize, timeType
      character(len=MAXCHR) dimName
      integer corner(3), edges(3)
      integer vid
      integer(INT64) seconds
      integer timeIndex
      integer minutes                       ! added as a work-around
      integer i, k
      integer begDate, begTime, timInc
      integer hour,min,sec,incSecs
      integer, allocatable ::  allTimes(:)
      integer fillTime

! Variables for dealing with precision

      real(kind=REAL32), allocatable :: grid_32(:,:)
      real(kind=REAL64), allocatable :: grid_64(:,:)
      real(kind=REAL32) dummy32
      real(kind=REAL64) dummy64
      real   dummy

! Variables for NF90_INQUIRE_VARIABLE

      character(len=MAXCHR) varName
      integer type, nvDims, vdims(MAXVDIMS), nvAtts

! Variables for packing and range checking

      integer(kind=INT16), allocatable :: grid_16(:,:)
      real(kind=REAL32), allocatable :: fminutes_32(:)
      real(kind=REAL32) high_32, low_32, amiss_32
      real(kind=REAL32) scale_32, offset_32
      logical outRange
      logical outPRange

      _UNUSED_DUMMY(jm)

! Variable initialization

      outRange = .FALSE.
      outPRange = .FALSE.

! Check to make sure max string lengths are large enough.  NetCDF defines
! MAXNCNAM, but it can't be used in a character(len=MAXNCNAM) statement.

      if (MAXCHR .LT. MAXNCNAM) then
        print *, 'CFIO_PutVar warning: MAXNCNAM is larger than ', &
                'dimName array size.'
      endif

! Determine NetCDF variable ID.

      rc = NF90_INQ_VARID (fid, vname, vid)
      if (err("PutVar: variable not defined",rc,-40) .NE. 0) return

! Determine number of seconds since starting date/time.

      rc = NF90_INQ_VARID (fid, 'time', timeId)
      if (err("PutVar: time not defined",rc,-43) .NE. 0) return
      rc = NF90_GET_ATT(fid,timeId,'begin_date',begDate)
      if (err("PutVar: missing begin_date",rc,-44) .NE. 0) return
      rc = NF90_GET_ATT(fid,timeId,'begin_time',begTime)
      if (err("PutVar: missing begin_time",rc,-44) .NE. 0) return

      seconds = DiffDate (begDate, begTime, yyyymmdd, hhmmss)

      if (seconds .lt. 0) then
        print *, 'CFIO_PutVar: Error code from diffdate.  Problem with', &
                ' date/time.'
        rc = -7
        return
      endif
      if ( MOD (int(seconds),60) .eq. 0 ) then
        minutes = seconds / 60
      else
        print *, 'CFIO_PutVar: Currently, times must fall on minute ',&
                'boundaries.'
        rc = -6
        return
      endif

! Confirm that this time is consistent with the starting time coupled with
! the time increment.

      rc = NF90_GET_ATT(fid,timeId,'time_increment',timInc)
      if (err("PutVar: missing time increment",rc,-44) .NE. 0) return

! Convert time increment to seconds.

!ams      write (strBuf,203) timinc
!ams 203   format (I6)
!ams       read (strBuf,204) hour, min, sec
!ams 204   format (3I2)
      call CFIO_parseIntTime ( timinc, hour, min, sec )
      incSecs = hour*3600 + min*60 + sec

      if ( MOD (int(seconds), incSecs) .ne. 0 ) then
        print *, 'CFIO_putvar: Absolute time of ',seconds,' not ',&
                'possible with an interval of ',incSecs
        rc = -2
        return
      else
        timeIndex = seconds/incSecs + 1
      endif

! Load starting indicies.

      if ( kbeg .eq. 0 ) then
        corner(1)=1
        corner(2)=timeIndex
        edges(1)=im
        edges(2)=1
      else
        corner(1)=1
        corner(2)=kbeg
        corner(3)=timeIndex
        edges(1)=im
        edges(2)=kount
        edges(3)=1
      endif

! Check variable against valid range.

      rc = NF90_GET_ATT(fid,vid,'vmin',low_32)
      if (err("PutVar: can't get vmin",rc,-53) .NE. 0) return
      rc = NF90_GET_ATT(fid,vid,'vmax',high_32)
      if (err("PutVar: can't get vmax",rc,-53) .NE. 0) return
      rc = NF90_GET_ATT(fid,vid,'fmissing_value',amiss_32)
      if (err("PutVar: can't get fmissing_value",rc,-53) .NE. 0) return
      if (abs(low_32) .NE. amiss_32 .OR. high_32 .NE. amiss_32) then
        do k=1,kount
          do i=1,im
            if (grid(i,k) .GT. high_32 .OR. grid(i,k) .LT. &
                low_32) then
              outRange = .TRUE.
              goto 100
            endif
          enddo
        enddo
100     continue
      endif

! Determine if we are writing single- or double-precision.

      rc = NF90_INQUIRE_VARIABLE (fid, vid, varName, type, nvDims, vDims, nvAtts)
      if (err("PutVar: error in variable inquire",rc,-52) .NE. 0) return

! Write variable in the appropriate precision.

      if (HUGE(dummy) .EQ. HUGE(dummy32)) then        ! -r4
        if (type .EQ. NF90_FLOAT) then                     ! 32-bit
          rc = NF90_PUT_VAR(fid,vid,grid,corner,edges)
        else if (type .EQ. NF90_DOUBLE) then               ! 64-bit
          allocate (grid_64(im,kount))
          do k=1,kount
              do i=1,im
                grid_64(i,k) = grid(i,k)
              enddo
          enddo
          rc = NF90_PUT_VAR(fid,vid,grid_64,corner,edges)
          deallocate (grid_64)
        else if (type .EQ. NF90_SHORT) then
          rc = NF90_GET_ATT(fid,vid,'packmax',high_32)
          if (err("PutVar: error getting packmax",rc,-53) .NE. 0) return
          rc = NF90_GET_ATT(fid,vid,'packmin',low_32)
          if (err("PutVar: error getting packmin",rc,-53) .NE. 0) return
          rc = NF90_GET_ATT(fid,vid,'scale_factor',scale_32)
          if (err("PutVar: error getting scale",rc,-53) .NE. 0) return
          rc = NF90_GET_ATT(fid,vid,'add_offset',offset_32)
          if (err("PutVar: error getting offset",rc,-53) .NE. 0) return
          allocate (grid_16(im,kount))
          do k=1,kount
              do i=1,im
                if ( grid(i,k) .LT. low_32 .OR. grid(i,k) .GT. &
                high_32) then
                  grid_16(i,k) = PACK_FILL
                  outPRange = .TRUE.
                else
                  grid_16(i,k) = (grid(i,k) - offset_32)/scale_32
                endif
              enddo
          enddo
          rc = NF90_PUT_VAR(fid,vid,grid_16,corner,edges)
          deallocate (grid_16)
        else
          rc = -13
          return
        endif
      else if (HUGE(dummy) .EQ. HUGE(dummy64)) then   ! -r8
        if (type .EQ. NF90_FLOAT) then                     ! 32-bit
          allocate (grid_32(im,kount))
          do k=1,kount
              do i=1,im
                grid_32(i,k) = grid(i,k)
              enddo
          enddo
          rc = NF90_PUT_VAR(fid, vid, grid_32, corner, edges)
          deallocate (grid_32)
        else if (type .EQ. NF90_DOUBLE) then                ! 64-bit
          rc = NF90_PUT_VAR(fid,vid,grid,corner,edges)
        else if (type .EQ. NF90_SHORT) then
          rc = NF90_GET_ATT(fid,vid,'packmax',high_32)
          if (err("PutVar: error getting packmax",rc,-53) .NE. 0) return
          rc = NF90_GET_ATT(fid,vid,'packmin',low_32)
          if (err("PutVar: error getting packmin",rc,-53) .NE. 0) return
          rc = NF90_GET_ATT(fid,vid,'scale_factor',scale_32)
          if (err("PutVar: error getting scale",rc,-53) .NE. 0) return
          rc = NF90_GET_ATT(fid,vid,'add_offset',offset_32)
          if (err("PutVar: error getting offset",rc,-53) .NE. 0) return
          allocate (grid_16(im,kount))
          do k=1,kount
              do i=1,im
                if ( grid(i,k) .LT. low_32 .OR. grid(i,k) .GT. &
                high_32) then
                  grid_16(i,k) = PACK_FILL
                  outPRange = .TRUE.
                else
                  grid_16(i,k) = (grid(i,k) - offset_32)/scale_32
                endif
              enddo
          enddo
          rc = NF90_PUT_VAR(fid,vid,grid_16,corner,edges)
          deallocate (grid_16)
        else
          rc = -13
          return
        endif
      else
        rc = -12
        return
      endif
      if (err("PutVar: error writing variable",rc,-45) .NE. 0) return

! Read time dimension scale and fill all values up to the current time.
! This will insure missing times are defined with the proper time value.

      rc = NF90_INQ_DIMID(fid, 'time', timeDimId)
      rc = NF90_INQUIRE_DIMENSION (fid, timeDimId, dimName, dimSize)
      dimSize = dimSize - 1           ! We've already written the
                                      ! the new time.
      allocate ( allTimes (MAX(timeIndex,dimSize)) )
      allocate ( fminutes_32 (MAX(timeIndex,dimSize)) )

      if (dimSize .GT. 0) then
        ! Depending on the version of CFIO used to write the file, the Time
        ! dimension variable can either be floating point or integer.

        corner(1)=1
        edges(1)=dimSize

        rc = NF90_INQUIRE_VARIABLE (fid,timeId,dimName,timeType,nvDims,vDims,nvAtts)
        if (timeType .EQ. NF90_FLOAT) then
          rc = NF90_GET_VAR(fid,timeId,fminutes_32,corner,edges)
          do i=1,dimSize
            allTimes(i) = INT(fminutes_32(i))
          enddo
        else if (timeType .EQ. NF90_INT) then
          rc = NF90_GET_VAR(fid,timeId,allTimes,corner,edges)
        endif
        if (err("SPutVar: error reading times from file",rc,-46) .NE. 0) &
            return
      endif

      ! This loop fills the time dimension scale based on the time increment
      ! specified in CFIO_Create.  If CFIO ever changes to support variable
      ! time increments, this code MUST be changed.

      do i=1,timeIndex-1
        fillTime = (i-1) * incSecs/60
        allTimes(i) = fillTime
      enddo
      allTimes(timeIndex) = minutes

! Write filled time array to file.

      corner(1)=1
      edges(1)=timeIndex

      if (timeType .EQ. NF90_FLOAT) then
        do i=1,timeIndex
          fminutes_32(i) = INT(allTimes(i))
        enddo
        rc = NF90_PUT_VAR(fid,timeId,fminutes_32,corner,edges)
      else if (timeType .EQ. NF90_INT) then
        rc = NF90_PUT_VAR(fid,timeId,allTimes,corner,edges)
      endif
      if (err("PutVar: error writing time",rc,-38) .NE. 0) return

      if (outRange .AND. outPRange) then
        rc = -17
      else if (outPRange) then
        rc = -16
      else if (outRange) then
        rc = -15
      else
        rc = 0
      endif

      deallocate ( allTimes )
      deallocate ( fminutes_32 )

      return
      end subroutine CFIO_SPutVar