CreateDstVar Subroutine

public subroutine CreateDstVar(srcFile, dstFile, fileType, srcVarName, dstVarName, varDims, varRank, varStr, locStr, varDimids, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: srcFile
character(len=*), intent(in) :: dstFile
type(ESMF_FileFormat_Flag), intent(in) :: fileType
character(len=*), intent(in) :: srcVarName
character(len=*), intent(in) :: dstVarName
integer, intent(in) :: varDims(:)
integer, intent(in) :: varRank
character(len=*), intent(in) :: varStr
character(len=*), intent(in) :: locStr
integer, intent(in) :: varDimids(:)
integer, intent(out) :: rc

Source Code

subroutine CreateDstVar(srcFile, dstFile, fileType, srcVarName, dstVarName, &
           varDims, varRank, varStr, locStr, varDimids, rc)
    character(len=*), intent(in) :: srcFile
    character(len=*), intent(in) :: dstFile
    type(ESMF_FILEFORMAT_FLAG), intent(in) :: fileType
    character(len=*), intent(in) :: srcVarName
    character(len=*), intent(in) :: dstVarName
    integer, intent(in) :: varDims(:)
    integer, intent(in) :: varRank
    character(len=*), intent(in) :: varStr
    character(len=*), intent(in) :: locStr
    integer, intent(in) :: varDimids(:)
    integer, intent(out) :: rc

    integer:: gridid, gridid1, varid, varid1, dimid
    integer:: extradim, dimval
    integer:: i, j, vartype, ndims, nattrs, rank
    integer, pointer :: dimids(:), srcdimids(:)
    character(len=40) ::dimnames(MAX_VARDIMS), dimname, attname
    integer :: ncStatus
    character(len=128) :: errmsg
    integer, parameter :: nf90_noerror = 0

#ifdef ESMF_NETCDF
    ! First check how many extra dimensions to be created
    if (fileType == ESMF_FILEFORMAT_UGRID) then
      extradim = varRank - 1
      rank = 1
    else
      extradim = varRank - 2
      rank = 2
    endif

    allocate(dimids(varRank))
    dimids(1:rank)=varDimids(1:rank)
    ! Open the destination file
    ncStatus = nf90_open (path=trim(dstFile), mode=nf90_write, ncid=gridid1)
    if (CDFCheckError (ncStatus, &
           ESMF_METHOD, &
           ESMF_SRCLINE,&
           trim(dstFile), &
           rc)) return

    ncStatus = nf90_open (path=trim(srcFile), mode=nf90_nowrite, ncid=gridid)
    if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        trim(srcFile), &
        rc)) return
    ncStatus = nf90_inq_varid( gridid, srcVarName, varid)
    errmsg = "variable "//trim(srcVarName)// " in "//trim(srcFile)
    if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg, &
        rc)) return
    ncStatus = nf90_inquire_variable(gridid, varid, xtype=vartype, ndims=ndims, &
             nAtts=nattrs)
    errmsg = "Variable "//trim(srcVarName)//" in "//trim(srcFile)
    if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg, &
        rc)) return
    if (extradim > 0) then
      ! get the dimension id and name for all the extra dimensions
      ! Open the src grid file
      allocate(srcdimids(ndims))
      ncStatus = nf90_inquire_variable(gridid, varid, &
               dimids=srcdimids)
      errmsg = "Variable "//trim(srcVarName)//" in "//trim(srcFile)
      if (CDFCheckError (ncStatus, &
        ESMF_METHOD, &
        ESMF_SRCLINE,&
        errmsg, &
        rc)) return
      do i=1,extradim
        ncStatus = nf90_inquire_dimension(gridid,srcdimids(ndims-(extradim-i)),&
                        name=dimnames(i))
        if (CDFCheckError (ncStatus, &
          ESMF_METHOD, &
          ESMF_SRCLINE,&
          errmsg, &
          rc)) return
      enddo
      ! check if the last dimension is "time"
      if (dimnames(extradim) .eq. "time") then
        ! the last dimension of the srcVariable is time, check if there is one in the dstFile
        ncStatus = nf90_inq_dimid(gridid1, 'time', dimid)
        if (ncStatus == nf90_noerror) then
           ncStatus = nf90_inquire_dimension(gridid1, dimid, len=dimval)
           errmsg = 'inquire time dimension in'//trim(dstFile)
           if (CDFCheckError (ncStatus, &
             ESMF_METHOD, &
             ESMF_SRCLINE,&
             errmsg, &
             rc)) return
           if (dimval /= varDims(varRank)) then
             call ESMF_LogSetError(ESMF_FAILURE, &
                 msg="- the time dimension values do not match in the src and dst files", &
                 ESMF_CONTEXT, rcToReturn=rc)
             return
           endif
           extradim = extradim - 1
           dimids(varRank) = dimid
        else
           ! create time dimension and other extra dimensions in the dstFile
           ncStatus = nf90_redef(gridid1)
           ncStatus = nf90_def_dim(gridid1, 'time', varDims(varRank) , dimids(varRank))
           errmsg = 'define time dimension in '//trim(dstFile)
           if (CDFCheckError (ncStatus, &
             ESMF_METHOD, &
             ESMF_SRCLINE,&
             errmsg, &
             rc)) return
           extradim = extradim - 1
           ncStatus = nf90_enddef(gridid1)
         endif
       endif ! there is not time dimension in the srcVar
       !if there are more extra dimension, create them
       if (extradim > 0) then
         j = 1
         do i=1,extradim
           write(dimname, '("extradim", I1)') j
           ! check if the dimension name already exist, if so and the value
           ! matches, use it, otherwise use another name
           ncStatus = nf90_inq_dimid(gridid1, dimname, dimid)
           do while (ncStatus == nf90_noerror)
              j=j+1
              write(dimname, '("extradim", I1)') j
              ncStatus = nf90_inq_dimid(gridid1, dimname, dimid)
           enddo        
           ncStatus = nf90_redef(gridid1)
           errmsg = 'enter redefine mode '//trim(dstFile)
           if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
           ncStatus = nf90_def_dim(gridid1, dimname, varDims(rank+i), &
                    dimids(rank+i))
           errmsg = 'define extra dimension in '//trim(dstFile)
           if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
         enddo
         ncStatus = nf90_enddef(gridid1)
         errmsg = 'end redefine mode '//trim(dstFile)
         if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
       endif
     endif ! End extradim > 0
     ! create the destination variable now
     ncStatus = nf90_redef(gridid1)
     errmsg = 'enter redefine mode '//trim(dstFile)
     if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
     ncStatus = nf90_def_var(gridid1, dstVarName, vartype, dimids, varid1)      
     errmsg = 'define destination variable '//trim(dstVarName)
     if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
     ! define the attributes for the variable, copy the attributes from the
     ! the srcVar and create new coordinates attribute
     do i=1,nattrs
       ncStatus = nf90_inq_attname(gridid, varid, i, attname)
       errmsg = 'attribute for '//trim(srcVarName)
       if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
       if ((attname .ne. 'coordinates') .and. (attname .ne. 'mesh')) then
         ncStatus = nf90_copy_att(gridid, varid, attname, gridid1, varid1)
         errmsg = 'adding attribute '//trim(attname)//' to variable '//trim(dstVarName)
         if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
       endif
     enddo
     if (fileType == ESMF_FILEFORMAT_GRIDSPEC .or. &
         fileType == ESMF_FILEFORMAT_MOSAIC ) then
       ! Add coordinates attribute
       ncStatus = nf90_put_att(gridid1, varid1, 'coordinates', varStr)
       errmsg = 'attribute for '//trim(dstVarName)
       if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
     else
       ! If the file is a UGRID, need to add mesh attribute too
       ncStatus = nf90_put_att(gridid1, varid1, 'mesh', varStr)
       errmsg = 'define attribute mesh for '//trim(dstVarName)
       if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
       ! Need to know if the variable is on face or node
       ncStatus = nf90_put_att(gridid1, varid1, 'location', locStr)
       errmsg = 'define attribute location for '//trim(dstVarName)
       if (CDFCheckError (ncStatus, &
               ESMF_METHOD, &
               ESMF_SRCLINE,&
               errmsg, &
               rc)) return
       ! Need to know if the variable is on face or node
     endif
     ncStatus = nf90_close(gridid)
     ncStatus = nf90_close(gridid1)

    rc = ESMF_SUCCESS
    return
#else
    call ESMF_LogSetError(ESMF_RC_LIB_NOT_PRESENT, &
                 msg="- ESMF_NETCDF not defined when lib was compiled", &
                 ESMF_CONTEXT, rcToReturn=rc)
    return
#endif

end subroutine CreateDstVar