ESMF_GridspecReadMosaic Subroutine

public subroutine ESMF_GridspecReadMosaic(filename, mosaic, tileFilePath, rc)

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: filename
type(ESMF_Mosaic), intent(inout) :: mosaic
character(len=*), intent(in), optional :: tileFilePath
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_GridspecReadMosaic(filename, mosaic, tileFilePath, rc)

! !ARGUMENTS:
 
    character(len=*), intent(in)           :: filename
    type(ESMF_Mosaic), intent(inout)       :: mosaic
    character(len=*), intent(in), optional :: tileFilePath
    integer, optional, intent(out)         :: rc

! !DESCRIPTIONS:
!   find a variable with the standard_name attribute set to "grid_mosaic_spec"
!   use the "children" attribute to find the variable that defines tile names
!   use the "contact_regions" attribute to find the variable that defines the contacts
!   The value of the variable is the mosaic name used in the contact list
!   
!   Variable "gridfiles" defines the path where the tile files are stored. If the 
!   optional argument tileFilePath is given, use it instead of the path defined in the
!   mosaic file.  Variable "gridfiles" defines the name of the tile files.  It is a 1D 
!   array with the same size as the "children" variable.
!
!   the "children" varible is a 1D array containing a list of tile names.  The dimension
!   of the variable defines the number of tiles
!
!   the "contact" variable defines the contacts of the tiles, it should have the following 
!   attributes:
!    standard_name = "grid_contact_spec"
!    contact_type = "boundary"
!    contact_index = variable name that defines how each contact is connected
!   This variable is a 1D array and its size defines the number of contacts in this mosaic
!   Each entry is a character string with the format:
!     mosaic_name:tilex::mosaic_name::tiley
!   Assuming mosaic_name should match with the value defined in the grid_mosaic variable, 
!    and the tilename is one of the tiles defined in the children variable.
!
!   the "contact_index" variable defines the two edges of the neighboring tiles at the boundary
!   It is a 1D array with the same size as the "contact" variable
!   Each entry has the following format:
!    i1,i2:j1,j2::i3,i4,j3,j4
!   where (i1,j1) and (i2, j2) are the two end points of edge in the first tile and 
!   (i3,j3) and (i4,j4) are the two end points of the edge in the second tile.  The index are
!   both local indices within the tile with 1-based index.
!   (i1,j1) is the same point as (i3,j3) and (i2,j2) is the same point as (i4,j4)
!
!EOPI
#ifdef ESMF_NETCDF

    integer   :: ncid, varid
    integer   :: i, j, nvars, attlen
    integer   :: ntiles, ncontacts
    character(len=128) :: attstr
    character(len=128) :: mosaicname
    character(len=ESMF_MAXPATHLEN) :: tempname
    integer            :: strlen
    character(len=1), allocatable :: temptilenames(:,:)
    character(len=1),  allocatable :: tilefilenames(:,:)
    integer, pointer :: contact(:,:)
    integer, pointer :: connindex(:,:,:)
    integer   :: dimids(2), ndims, dims(2)
    integer   :: ncStatus, localrc
    integer   :: totallen
    integer   :: k

    !initialize mosaic values
    mosaic%ntiles = 0

    if (present(rc)) rc=ESMF_FAILURE

    ncStatus = nf90_open(path=filename, mode=nf90_nowrite, ncid=ncid)
    if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, &
        trim(filename), &
        rc)) return

    ncStatus = nf90_inquire(ncid, nVariables=nvars)
    if (CDFCheckError (ncStatus, &
        ESMF_METHOD,  &
        ESMF_SRCLINE, &
        trim(filename), &
        rc)) return
     do i=1,nvars
       ! Check its standard_name attribute
       ncStatus = nf90_inquire_attribute(ncid, i, 'standard_name', len=attlen)
       if (ncStatus == nf90_noerr) then
          ncStatus = nf90_get_att(ncid, i, 'standard_name', values=attstr)
          if (attstr(1:attlen) .eq. 'grid_mosaic_spec') then
            ! get mosaic name
            ncStatus = nf90_get_var(ncid, i, mosaicname)
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               trim(filename), &
               rc)) return
            mosaic%name = trim(mosaicname)
            ! Get attributes children and contact_regions
            ncStatus = nf90_get_att(ncid, i, 'children', values=attstr)
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "children attribute does not exit", &
               rc)) return
            ! Get the tilenames from the variable defined in "children" variable
            ncStatus = nf90_inq_varid(ncid, trim(attstr), varid)
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "children variable does not exit", &
               rc)) return
            ! Find the dimension of the tile variable
            ncStatus = nf90_inquire_variable(ncid, varid, ndims=ndims, dimids=dimids)
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "fail to inquire the variable defined by children attribute", &
               rc)) return
            if (ndims /= 2) then
                call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg="- tile variable dimension greater than 1", & 
                 ESMF_CONTEXT, rcToReturn=rc) 
               return
            endif
            ! query it's dimension to find out number of tiles
            ncStatus = nf90_inquire_dimension(ncid, dimids(2), len=ntiles)
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "tile dimension inquire", &
               rc)) return
            ncStatus = nf90_inquire_dimension(ncid, dimids(1), len=strlen)
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "tile dimension inquire", &
               rc)) return
             mosaic%ntiles = ntiles
             ! get the tile names
             allocate(mosaic%tilenames(ntiles), temptilenames(strlen, ntiles))
             ncStatus = nf90_get_var(ncid, varid, temptilenames, start=(/1,1/), count=(/strlen, ntiles/))
             if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "get tile names", &
               rc)) return
              ! replace null character by blank
             do j=1,ntiles
              mosaic%tilenames(j) = ESMF_UtilArray2String(temptilenames(:,j))
              call trim_null(mosaic%tilenames(j))
             enddo
             
            if (ntiles > 1) then
              ! Find contact regions
              ncStatus = nf90_get_att(ncid, i, 'contact_regions', values=attstr)
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD,  &
                 ESMF_SRCLINE, &
                 "contact_regions attribute does not exit", &
                 rc)) return
              ! Get the contact info from the variable defined in "contact_regions" variable
              ncStatus = nf90_inq_varid(ncid, trim(attstr), varid)
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD,  &
                 ESMF_SRCLINE, &
                 "contact_regions variable does not exit", &
                 rc)) return
              ! Find the dimension of this variable
              ncStatus = nf90_inquire_variable(ncid, varid, ndims=ndims, dimids=dimids)
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD,  &
                 ESMF_SRCLINE, &
                 "fail to inquire the variable defined by contact_regions", &
                 rc)) return
              if (ndims /= 2) then
                 call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                     msg="- contact region variable dimension greater than 1", & 
                     ESMF_CONTEXT, rcToReturn=rc) 
                 return
              endif
              ! query it's dimension to find out the string length
              ncStatus = nf90_inquire_dimension(ncid, dimids(1), len=dims(1))
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD,  &
                 ESMF_SRCLINE, &
                 "contact dimension inquire", &
                 rc)) return
              ! query it's dimension to find out number of contacts
              ncStatus = nf90_inquire_dimension(ncid, dimids(2), len=dims(2))
              if (CDFCheckError (ncStatus, &
                 ESMF_METHOD,  &
                 ESMF_SRCLINE, &
                 "contact dimension inquire", &
                 rc)) return
              allocate(mosaic%contact(2,dims(2)))
              allocate(mosaic%connindex(2,4,dims(2)))
              call readContacts(ncid, varid, dims, mosaicname, mosaic%tilenames, &
                   mosaic%contact, mosaic%connindex, localrc)
              if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
                 ESMF_CONTEXT, rcToReturn=rc)) return
              mosaic%ncontacts=dims(2)
            endif
         endif
       endif
     enddo

     ! find the tile file path
     if (present(tileFilePath)) then
       mosaic%tileDirectory = tileFilePath
     else
       ncStatus = nf90_inq_varid(ncid, 'gridlocation', varid)
       if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "inquire variable gridlocation", &
               rc)) return
       ! Get the directory name where the tile files are stored
       ncStatus = nf90_inquire_variable(ncid, varid, ndims=ndims, dimids=dimids)
       if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "inquire variable gridlocation", &
               rc)) return
       ncStatus = nf90_inquire_dimension(ncid, dimids(1), len=dims(1))
            if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "contact dimension inquire", &
               rc)) return
       ! initialize the string to null characters
       do k=1,ESMF_MAXPATHLEN
           mosaic%tileDirectory(k:k)=char(0)
       enddo
       ncStatus = nf90_get_var(ncid, varid, mosaic%tileDirectory, start=(/1/), count=(/dims(1)/))
       if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "fail to get gridlocation", &
               rc)) return
       call trim_null(mosaic%tileDirectory)
     endif

     ! get the tile file names
     ncStatus = nf90_inq_varid(ncid, 'gridfiles', varid)
     if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "inquire variable gridfiles", &
               rc)) return
       
     ! Find the dimension of this variable
     ncStatus = nf90_inquire_variable(ncid, varid, ndims=ndims, dimids=dimids)
     if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "inquire variable gridfiles error", &
               rc)) return
     if (ndims /= 2) then
        call ESMF_LogSetError(rcToCheck=ESMF_FAILURE, &
                 msg="- gridfiles variable dimension /= 2", & 
                 ESMF_CONTEXT, rcToReturn=rc) 
        return
     endif    
     ! query it's dimension to find out the string length
     ncStatus = nf90_inquire_dimension(ncid, dimids(1), len=dims(1))
     if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "gridfiles dimension inquire", &
               rc)) return
     allocate(tilefilenames(dims(1), ntiles))
     allocate(mosaic%filenames(ntiles))
     ! initialize the string to null characters
     do j=1,ntiles
        do k=1,ESMF_MAXPATHLEN
            mosaic%filenames(ntiles)(k:k)=char(0)
        enddo
     enddo
     ncStatus = nf90_get_var(ncid, varid, tilefilenames, start=(/1,1/), count=(/dims(1), ntiles/))
     if (CDFCheckError (ncStatus, &
               ESMF_METHOD,  &
               ESMF_SRCLINE, &
               "fail to get gridfiles", &
               rc)) return
     ! replace null character by blank
     do j=1,ntiles
          mosaic%filenames(j)=ESMF_UtilArray2String(tilefilenames(:,j))
          call trim_null(mosaic%filenames(j))
     enddo

      ncStatus = nf90_close(ncid)
      if (CDFCheckError (ncStatus, &
           ESMF_METHOD,  &
           ESMF_SRCLINE, &
           "close mosaic file", &
           rc)) return

      ! check if we found the mosaic variables
      if (mosaic%ntiles == 0) then
         call ESMF_LogSetError(ESMF_FAILURE, & 
                 msg="- Failed to parse GridSpec Mosaic file", & 
                 ESMF_CONTEXT, rcToReturn=rc) 
         return    
       else
         ! find the dimension of the tile by reading one of the tilefiles
         ! assuming all the tilefiles have the same size
         tempname = trim(mosaic%TileDirectory)//trim(mosaic%filenames(1))
         call ESMF_GridSpecQueryTileSize(tempname, mosaic%nx, mosaic%ny, rc=localrc)
         if (ESMF_LogFoundError(localrc, ESMF_ERR_PASSTHRU, &
            ESMF_CONTEXT, rcToReturn=rc)) return
         if (present(rc)) rc=ESMF_SUCCESS
         return
       endif

#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 ESMF_GridspecReadMosaic