ESMF_XGridConstruct Subroutine

private subroutine ESMF_XGridConstruct(xgtype, sideA, sideB, sparseMatA2X, sparseMatX2A, sparseMatB2X, sparseMatX2B, offline, mesh, internal_alloc, rc)

Arguments

Type IntentOptional Attributes Name
type(ESMF_XGridType), intent(inout) :: xgtype
type(ESMF_XGridGeomBase), intent(in) :: sideA(:)
type(ESMF_XGridGeomBase), intent(in) :: sideB(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatA2X(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatX2A(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatB2X(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatX2B(:)
logical, intent(in), optional :: offline
type(ESMF_Mesh), intent(inout), optional :: mesh
logical, intent(in), optional :: internal_alloc
integer, intent(out), optional :: rc

Source Code

subroutine ESMF_XGridConstruct(xgtype, sideA, sideB, &
    sparseMatA2X, sparseMatX2A, sparseMatB2X, sparseMatX2B, offline, &
    mesh, internal_alloc, rc)
!
! !ARGUMENTS:
type(ESMF_XGridType), intent(inout)        :: xgtype
type(ESMF_XGridGeomBase), intent(in)       :: sideA(:), sideB(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatA2X(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatX2A(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatB2X(:)
type(ESMF_XGridSpec), intent(in), optional :: sparseMatX2B(:)
logical, intent(in), optional              :: offline
type(ESMF_Mesh), intent(inout), optional   :: mesh
logical, intent(in), optional              :: internal_alloc
integer, intent(out), optional             :: rc 

!
! !DESCRIPTION:
!     Construct internals of xgtype from input
!
!     The arguments are:
!     \begin{description}
!     \item [xgtype]
!           the {ESMF\_XGridType} object.
!     \item [sideA]
!           2D Grids on side A.
!     \item [sideB]
!           2D Grids on side B.
!     \item [{[sparseMatA2X]}]
!           indexlist from a Grid index space on side A to xgrid index space;
!           indexFactorlist from a Grid index space on side A to xgrid index space.
!     \item [{[sparseMatX2A]}]
!           indexlist from xgrid index space to a Grid index space on side A;
!           indexFactorlist from xgrid index space to a Grid index space on side A.
!     \item [{[sparseMatB2X]}]
!           indexlist from a Grid index space on side B to xgrid index space;
!           indexFactorlist from a Grid index space on side B to xgrid index space.
!     \item [{[sparseMatX2B]}]
!           indexlist from xgrid index space to a Grid index space on side B;
!           indexFactorlist from xgrid index space to a Grid index space on side B.
!     \item [{[offline]}]
!           online generation optimization turned on/off (default off)
!     \item [{[mesh]}]
!           online generation with mesh
!     \item [{[rc]}]
!           Return code; equals {\tt ESMF\_SUCCESS} only if successful.
!     \end{description}
!
!EOPI

  integer :: localrc, ngrid_a, ngrid_b
  integer :: i
  logical :: l_offline
  real(ESMF_KIND_R8), pointer :: xgrid_frac(:)

  localrc = ESMF_SUCCESS

  ! Initialize return code   
  if(present(rc)) rc = ESMF_RC_NOT_IMPL
  l_offline = .true.
  if(present(offline)) l_offline = offline

  ngrid_a = size(sideA, 1)
  ngrid_b = size(sideB, 1)

  ! check and copy all the sparse matrix spec structures
  if(present(sparseMatA2X) .and. l_offline) then
      call ESMF_SparseMatca(sparseMatA2X, xgtype%sparseMatA2X, ngrid_a, &
        'sparseMatA2X', rc=localrc)
      if (ESMF_LogFoundAllocError(localrc, &
          msg="- Initializing xgtype%sparseMatX2A ", &
          ESMF_CONTEXT, rcToReturn=rc)) return
  endif

  if(present(sparseMatX2A) .and. l_offline) then
      call ESMF_SparseMatca(sparseMatX2A, xgtype%sparseMatX2A, ngrid_a, &
        'sparseMatX2A', rc=localrc)
      if (ESMF_LogFoundAllocError(localrc, &
          msg="- Initializing xgtype%sparseMatX2A ", &
          ESMF_CONTEXT, rcToReturn=rc)) return
  endif

  ! TODO:
  ! if both A2X and X2A are present, check the sequence index list of X are identical
  ! this checking will be collective since the indices needs to be gathered
  ! if(present(sparseMatA2X) .and. present(sparseMatX2A)) then
  ! endif
  ! Another approach is to create 2 distgrids and use distgridMatch to compare
  ! the result Distgrid as discussed.

  if(present(sparseMatB2X) .and. l_offline) then
      call ESMF_SparseMatca(sparseMatB2X, xgtype%sparseMatB2X, ngrid_b, &
        'sparseMatB2X', rc=localrc)
      if (ESMF_LogFoundAllocError(localrc, &
          msg="- Initializing xgtype%sparseMatB2X ", &
          ESMF_CONTEXT, rcToReturn=rc)) return
  endif

  if(present(sparseMatX2B) .and. l_offline) then
      call ESMF_SparseMatca(sparseMatX2B, xgtype%sparseMatX2B, ngrid_b, &
        'sparseMatX2B', rc=localrc)
      if (ESMF_LogFoundAllocError(localrc, &
          msg="- Initializing xgtype%sparseMatX2B ", &
          ESMF_CONTEXT, rcToReturn=rc)) return
  endif

  ! TODO:
  ! if both B2X and X2B are present, check the sequence index list of X are identical
  ! this checking will be collective since the indices needs to be gathered
  ! if(present(sparseMatA2X) .and. present(sparseMatX2A)) then
  ! endif

  ! create the distgrids
  if((.not. l_offline) .and. present(mesh)) then
    call ESMF_XGridDistGridsOnline(xgtype, mesh, rc=localrc)
    if (ESMF_LogFoundError(localrc, &
        ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
  else
    call ESMF_XGridDistGrids(xgtype, rc=localrc)
    if (ESMF_LogFoundError(localrc, &
        ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
  endif

  ! Create the fracX here because we have a non-vanishing distgridM at this point and
  ! we know its entries should always be 1.0. This could be left for the user but it's 
  ! provided here so regridstore call retrieve this Field directly either as src or dst Frac.
  if(.not. l_offline) then
    xgtype%fracX = ESMF_ArrayCreate(xgtype%distgridM, typekind=ESMF_TYPEKIND_R8, rc=localrc)
    if (ESMF_LogFoundError(localrc, &
        ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
    ! mesh distgrid should always 1 de/pet
    call ESMF_ArrayGet(xgtype%fracX, localDe=0, farrayPtr=xgrid_frac, rc=localrc)
    if (ESMF_LogFoundError(localrc, &
        ESMF_ERR_PASSTHRU, &
        ESMF_CONTEXT, rcToReturn=rc)) return
    xgrid_frac = 1.0
  endif

  if(present(rc)) rc = ESMF_SUCCESS

end subroutine ESMF_XGridConstruct