MOAB: Mesh Oriented datABase  (version 5.1.1)
Adjacencies
+ Collaboration diagram for Adjacencies:

Functions

void iMesh_getAdjTable (iMesh_Instance instance, int **adjacency_table, int *adjacency_table_allocated, int *adjacency_table_size, int *err)
 Get the adjacency table for this implementation.
void iMesh_setAdjTable (iMesh_Instance instance, int *adjacency_table, int adjacency_table_size, int *err)
 Set the adjacency table as requested by the application.
void iMesh_getEntArrAdj (iMesh_Instance instance, const iBase_EntityHandle *entity_handles, const int entity_handles_size, const int entity_type_requested, iBase_EntityHandle **adjacentEntityHandles, int *adjacentEntityHandles_allocated, int *adj_entity_handles_size, int **offset, int *offset_allocated, int *offset_size, int *err)
 Get entities of specified type adjacent to entities.
void iMesh_getEntArr2ndAdj (iMesh_Instance instance, iBase_EntityHandle const *entity_handles, int entity_handles_size, int bridge_entity_type, int requested_entity_type, iBase_EntityHandle **adj_entity_handles, int *adj_entity_handles_allocated, int *adj_entity_handles_size, int **offset, int *offset_allocated, int *offset_size, int *err)
 Get "2nd order" adjacencies to an array of entities.
void iMesh_getAdjEntIndices (iMesh_Instance instance, iBase_EntitySetHandle entity_set_handle, int entity_type_requester, int entity_topology_requester, int entity_type_requested, iBase_EntityHandle **entity_handles, int *entity_handles_allocated, int *entity_handles_size, iBase_EntityHandle **adj_entity_handles, int *adj_entity_handles_allocated, int *adj_entity_handles_size, int **adj_entity_indices, int *adj_entity_indices_allocated, int *adj_entity_indices_size, int **offset, int *offset_allocated, int *offset_size, int *err)
 Get indexed representation of mesh or subset of mesh.
void iMesh_getEntAdj (iMesh_Instance instance, const iBase_EntityHandle entity_handle, const int entity_type_requested, iBase_EntityHandle **adj_entity_handles, int *adj_entity_handles_allocated, int *adj_entity_handles_size, int *err)
 Get entities of specified type adjacent to an entity.
void iMesh_getEnt2ndAdj (iMesh_Instance instance, iBase_EntityHandle entity_handle, int bridge_entity_type, int requested_entity_type, iBase_EntityHandle **adjacent_entities, int *adjacent_entities_allocated, int *adjacent_entities_size, int *err)
 Get "2nd order" adjacencies to an entity.

Function Documentation

void iMesh_getAdjEntIndices ( iMesh_Instance  instance,
iBase_EntitySetHandle  entity_set_handle,
int  entity_type_requester,
int  entity_topology_requester,
int  entity_type_requested,
iBase_EntityHandle **  entity_handles,
int *  entity_handles_allocated,
int *  entity_handles_size,
iBase_EntityHandle **  adj_entity_handles,
int *  adj_entity_handles_allocated,
int *  adj_entity_handles_size,
int **  adj_entity_indices,
int *  adj_entity_indices_allocated,
int *  adj_entity_indices_size,
int **  offset,
int *  offset_allocated,
int *  offset_size,
int *  err 
)

Get indexed representation of mesh or subset of mesh.

Given an entity set and optionally a type or topology, return:

  • The entities in the set of the specified type or topology
  • The entities adjacent to those entities with a specified type, as a list of unique handles.
  • For each entity in the first list, the adjacent entities, specified as indices into the second list.

Note 1: Because 'adjacent' as defined by the iMesh data model refers to those entities that bound another, the entities being queried here (in entity_set_handle arg) are NEVER ALSO returned in adj_entity_handles even if the entity_type_requested matches the entity type(s) in entity_set_handle. Note 2: The entities adjacent to the ith entity in entity_handles are found in adj_entity_handles running from offset[i] to offset[i+1] - 1. This implies that the offset_size will be entity_handles_size + 1. Note 3: This function will fail and return an error if the caller passes a combination of entity_type and entity_topology that are not consistent.

Parameters:
[in]instanceiMesh instance handle
[in]entity_set_handleThe set of entities from which to query
[in]entity_type_requesterIf not iBase_ALL_TYPES, act only on the subset of entities contained in entity_set_handle of the specified type.
[in]entity_topology_requesterIf not iMesh_ALL_TOPOLOGIES, act only on the subset of entities contained in entity_set_handle of specified topology.
[in]entity_type_requestedThe type of the adjacent entities to return
[in,out]entity_handlesThe handles of the (non-struct) subset of the entities contained in entity_set_handle indicated by the optional type and topology filtering arguments. Array pointer, allocated and occupied sizes argument trio)
[in,out]entity_handles_allocatedAllocated size of entity_handles array
[out]entity_handles_sizeOccupied size of entity_handles array
[in,out]adj_entity_handlesThe union of the unique entities of type requested_entity_type adjacent to each entity in entity_handles. Note that the implicit INTERLEAVED storage order rule applies (see iBase_StorageOrder) Array pointer, allocated and occupied sizes argument trio)
[in,out]adj_entity_handles_allocatedAllocated size of adj_entity_handles array
[out]adj_entity_handles_sizeOccupied size of adj_entity_handles array
[in,out]adj_entity_indicesFor each entity in entity_handles, the adjacent entities of type entity_type_requested, specified as indices into adj_entity_handles. The values are concatenated into a single array in the order of the entity handles in entity_handles. Array pointer, allocated and occupied sizes argument trio)
[in,out]adj_entity_indices_allocatedAllocated size of adj_entity_indices array
[out]adj_entity_indices_sizeOccupied size of adj_entity_indices array
[in,out]offsetFor each entity in the corresponding position in entity_handles, the position in adj_entity_indices at which values for that entity are stored Array pointer, allocated and occupied sizes argument trio)
[in,out]offset_allocatedAllocated size of offset array
[out]offset_sizeOccipied size of offset array
[out]errReturned Error status (see iBase_ErrorType)

Definition at line 830 of file iMesh_MOAB.cpp.

References iBase_BAD_ARRAY_DIMENSION, iBase_MEMORY_ALLOCATION_FAILED, iBase_SUCCESS, iMesh_getEntArrAdj, iMesh_getEntities, and size.

  {
    const int allocated_entity_handles = (*entity_handles_allocated == 0);
    const int allocated_indices = (*adj_entity_indices_allocated == 0);
    const int allocated_offset = (*offset_allocated == 0);

    // get source entities
    iMesh_getEntities( instance,
                       entity_set_handle,
                       entity_type_requestor,
                       entity_topology_requestor,
                       entity_handles,
                       entity_handles_allocated,
                       entity_handles_size,
                       err );
    if (iBase_SUCCESS != *err)
      return;

    // get adjacencies
    iBase_EntityHandle* all_adj_handles = 0;
    int size = 0, alloc = 0;
    iMesh_getEntArrAdj( instance,
                        *entity_handles, *entity_handles_size,
                        entity_type_requested,
                        &all_adj_handles, &alloc, &size,
                        offset, offset_allocated, offset_size,
                        err );
    if (*err != iBase_SUCCESS) {
      if (allocated_entity_handles) {
        free( *entity_handles );
        *entity_handles = 0;
        *entity_handles_allocated = 0;
      }
      return;
    }

    // allocate or check size of adj_entity_indices
    *adj_entity_indices_size = size;
    if (allocated_indices) {
      *adj_entity_indices = (int*)malloc(sizeof(iBase_EntityHandle)*size);
      if (!*adj_entity_indices)
        *err = iBase_MEMORY_ALLOCATION_FAILED;
      else
        *adj_entity_indices_allocated = size;
    }
    else if (*adj_entity_indices_allocated < size) {
      *err = iBase_BAD_ARRAY_DIMENSION;
    }
    if (iBase_SUCCESS != *err) {
      free( all_adj_handles );
      if (allocated_entity_handles) {
        free( *entity_handles );
        *entity_handles = 0;
        *entity_handles_allocated = 0;
      }
      if (allocated_offset) {
        free( *offset );
        *offset = 0;
        *offset_allocated = 0;
      }
      return;
    }

    // Now create an array of unique sorted handles from all_adj_handles.
    // We need to create a copy because we still need all_adj_handles.  We
    // will eventually need to copy the resulting unique list into
    // adj_entity_handles, so if adj_entity_handles is already allocated and
    // of sufficient size, use it rather than allocating another temporary.
    iBase_EntityHandle* unique_adj = 0;
    if (*adj_entity_handles_allocated >= size) {
      unique_adj = *adj_entity_handles;
    }
    else {
      unique_adj = (iBase_EntityHandle*)malloc(sizeof(iBase_EntityHandle) * size);
    }
    std::copy( all_adj_handles, all_adj_handles+size, unique_adj );
    std::sort( unique_adj, unique_adj + size );
    *adj_entity_handles_size = std::unique( unique_adj, unique_adj + size ) - unique_adj;

    // If we created a temporary array for unique_adj rather than using
    // already allocated space in adj_entity_handles, allocate adj_entity_handles
    // and copy the unique handle list into it
    if (*adj_entity_handles != unique_adj) {
      if (!*adj_entity_handles_allocated) {
        *adj_entity_handles = (iBase_EntityHandle*)malloc(
                              sizeof(iBase_EntityHandle) * *adj_entity_handles_size);
        if (!*adj_entity_handles)
          *err = iBase_MEMORY_ALLOCATION_FAILED;
        else
          *adj_entity_handles_allocated = *adj_entity_handles_size;
      }
      else if (*adj_entity_handles_allocated < *adj_entity_handles_size)
        *err = iBase_BAD_ARRAY_DIMENSION;
      if (iBase_SUCCESS != *err) {
        free( unique_adj );
        free( all_adj_handles );
        if (allocated_entity_handles) {
          free( *entity_handles );
          *entity_handles = 0;
          *entity_handles_allocated = 0;
        }
        if (allocated_offset) {
          free( *offset );
          *offset = 0;
          *offset_allocated = 0;
        }
        if (allocated_indices) {
          free( *adj_entity_indices );
          *adj_entity_indices = 0;
          *adj_entity_indices_allocated = 0;
        }
        return;
      }

      std::copy( unique_adj, unique_adj + *adj_entity_handles_size, *adj_entity_handles );
      free( unique_adj );
      unique_adj = *adj_entity_handles;
    }

    // convert from adjacency list to indices into unique_adj
    for (int i = 0; i < *adj_entity_indices_size; ++i)
      (*adj_entity_indices)[i] = std::lower_bound( unique_adj,
        unique_adj + *adj_entity_handles_size, all_adj_handles[i] ) - unique_adj;
    free( all_adj_handles );
  }
void iMesh_getAdjTable ( iMesh_Instance  instance,
int **  adjacency_table,
int *  adjacency_table_allocated,
int *  adjacency_table_size,
int *  err 
)

Get the adjacency table for this implementation.

Get the adjacency table for this implementation. This table is a 4x4 array whose entries characterize how the implementation behaves when adjacent and intermediate entities are queried. Entry [i,j] (i=row, j=col, 0-based indexing) will have one of the values in the iBase_AdjacencyCost enumeration. Off-diagonal entres, i!=j, represents the relative cost of retrieving entities of dimension i adjacent to entities of dimension j. Diagonal entries, i==j, indicate whether or not handles to ALL entities of that dimension are obtainable from calls that return entity handles. This is always true by definition for i==j==0 as well as for i==j==2 in a 2D mesh and i==j==3 in a 3D mesh. However, diagonal entries [1,1] for a 2D mesh and both [1,1] and [2,2] for a 3D mesh indicate whether or not handles to ALL intermediate dimensioned entities (ALL edges in a 2D mesh or ALL edges and ALL faces in a 3D mesh) are obtainable from calls that return entity handles. A value of iBase_AVAILABLE for a diagonal entry indicates that handles are obtainable for ALL entities of that dimension while a value of iBase_UNAVAILABLE indicates that handles are not obtainable for ALL entities of that dimension.

Parameters:
[in]instanceiMesh instance handle
[out]adjacency_tablePointer to array representing adjacency table Array pointer, allocated and occupied sizes argument trio)
[in,out]adjacency_table_allocatedPointer to allocated size of adjacency_table
[out]adjacency_table_sizePointer to occupied size of adjacency_table
[out]errReturned Error status (see iBase_ErrorType)

Definition at line 411 of file iMesh_MOAB.cpp.

References ALLOC_CHECK_ARRAY_NOFAIL, iBase_SUCCESS, iMesh_getGeometricDimension, MBIMESHI, munge_adj_table(), and RETURN.

  {
    int geom_dim;
    iMesh_getGeometricDimension(instance, &geom_dim, err);

    ALLOC_CHECK_ARRAY_NOFAIL(adjacency_table, 16);
    memcpy(*adjacency_table, MBIMESHI->AdjTable, 16*sizeof(int));
    munge_adj_table(*adjacency_table, geom_dim);
    RETURN(iBase_SUCCESS);
  }
void iMesh_getEnt2ndAdj ( iMesh_Instance  instance,
iBase_EntityHandle  entity_handle,
int  bridge_entity_type,
int  requested_entity_type,
iBase_EntityHandle **  adjacent_entities,
int *  adjacent_entities_allocated,
int *  adjacent_entities_size,
int *  err 
)

Get "2nd order" adjacencies to an entity.

Get "2nd order" adjacencies to an entity, that is, from an entity, through other entities of a specified "bridge" dimension, to other entities of another specified "to" dimension. Note 1: If the "bridge" dimension is the same as the "to" dimension or the dimension of the input entity, the output will be empty (and an error code of iBase_INVALID_ARGUMENT returned). This is consistent with the definition of adjacencies and the behavior of iMesh first adjacency calls. Note 2: An entity will never be returned as a second adjacency of itself, on the grounds that this is the most likely expectation of applications, and that it is easier for an application to add the original entity to the returned data than to find and remove it.

Parameters:
[in]instanceiMesh instance handle
[in]entity_handleEntity from which adjacencies are requested
[in]bridge_entity_typeType of bridge entity for 2nd order adjacencies
[in]requested_entity_typeType of adjacent entities returned
[in,out]adjacent_entitiesAdjacent entities Array pointer, allocated and occupied sizes argument trio)
[in,out]adjacent_entities_allocatedAllocated size of returned array
[out]adjacent_entities_sizeOccupied size of returned array
[out]errReturned Error status (see iBase_ErrorType)

Definition at line 2379 of file iMesh_MOAB.cpp.

References iMesh_getEntArr2ndAdj.

  {
    int offsets[2];
    int *offsets_ptr = offsets;
    int offset_size, offset_allocated = 2;

    iMesh_getEntArr2ndAdj(instance,
                          &entity_handle, 1, order_adjacent_key,
                          requested_entity_type,
                          adj_entities, adj_entities_allocated,
                          adj_entities_size, &offsets_ptr, &offset_allocated,
                          &offset_size, err);
  }
void iMesh_getEntAdj ( iMesh_Instance  instance,
const iBase_EntityHandle  entity_handle,
const int  entity_type_requested,
iBase_EntityHandle **  adj_entity_handles,
int *  adj_entity_handles_allocated,
int *  adj_entity_handles_size,
int *  err 
)

Get entities of specified type adjacent to an entity.

Get entities of specified type adjacent to an entity. Specified type must be value in the iBase_EntityType enumeration.

Note 1: Because 'adjacent' as defined by the iMesh data model refers to those entities that bound another, the entity being queried here (in entity_handle arg) is NEVER ALSO returned in adj_entity_handles even if the entity_type_requested matches the entity type in entity_handle.

Parameters:
[in]instanceiMesh instance handle
[in]entity_handleEntity handle being queried
[in]entity_type_requestedType of adjacent entities requested
[in,out]adj_entity_handlesPointer to array of adjacent entities Array pointer, allocated and occupied sizes argument trio)
[in,out]adj_entity_handles_allocatedPointer to allocated size of adj_entity_handles
[out]adj_entity_handles_sizePointer to occupied size of adj_entity_handles
[out]errReturned Error status (see iBase_ErrorType)

Definition at line 2361 of file iMesh_MOAB.cpp.

References iMesh_getEntArrAdj.

  {
    int offsets[2];
    int *offsets_ptr = offsets;
    int offset_size, offset_allocated = 2;

    iMesh_getEntArrAdj(instance,
                       &entity_handle, 1, entity_type_requested,
                       adj_entity_handles, adj_entity_handles_allocated,
                       adj_entity_handles_size, &offsets_ptr, &offset_allocated,
                       &offset_size, err);
  }
void iMesh_getEntArr2ndAdj ( iMesh_Instance  instance,
iBase_EntityHandle const *  entity_handles,
int  entity_handles_size,
int  bridge_entity_type,
int  requested_entity_type,
iBase_EntityHandle **  adj_entity_handles,
int *  adj_entity_handles_allocated,
int *  adj_entity_handles_size,
int **  offset,
int *  offset_allocated,
int *  offset_size,
int *  err 
)

Get "2nd order" adjacencies to an array of entities.

Get "2nd order" adjacencies to an array of entities, that is, from each entity, through other entities of a specified "bridge" dimension, to other entities of another specified "to" dimension. Note 1: If the "bridge" dimension is the same as the "to" dimension, the output will be empty (and an error code of iBase_INVALID_ARGUMENT returned). If the type of a particular entity matches the "bridge" dimension, there will be no entities returned for that input entity. This is consistent with the definition of adjacencies and the behavior of iMesh first adjacency calls. Note 2: An entity will never be returned as a second adjacency of itself, on the grounds that this is the most likely expectation of applications, and that it is easier for an application to add the original entity to the returned data than to find and remove it. Note 3: The entities adjacent to the ith entity in entity_handles are found in adj_entity_handles running from offset[i] to offset[i+1] - 1. This implies that the offset_size will be entity_handles_size + 1.

Parameters:
[in]instanceiMesh instance handle
[in]entity_handlesEntities from which adjacencies are requested
[in]entity_handles_sizeNumber of entities whose adjacencies are requested
[in]bridge_entity_typeType of bridge entity for 2nd order adjacencies
[in]requested_entity_typeType of adjacent entities returned
[in,out]adj_entity_handlesAdjacent entities. Note that the implicit INTERLEAVED storage order rule applies (see iBase_StorageOrder) Array pointer, allocated and occupied sizes argument trio)
[in,out]adj_entity_handles_allocatedAllocated size of returned array
[out]adj_entity_handles_sizeOccupied size of returned array
[in,out]offsetOffset[i] is offset into adj_entity_handles of 2nd order adjacencies of ith entity in entity_handles Array pointer, allocated and occupied sizes argument trio)
[in,out]offset_allocatedAllocated size of offset array
[out]offset_sizeOccupied size of offset array
[out]errReturned Error status (see iBase_ErrorType)

Definition at line 755 of file iMesh_MOAB.cpp.

References ALLOC_CHECK_ARRAY, ALLOC_CHECK_ARRAY_NOFAIL, moab::Range::begin(), CHKENUM, CHKERR, moab::Range::end(), ERROR, moab::MeshTopoUtil::get_bridge_adjacencies(), iBase_ALL_TYPES, iBase_INVALID_ARGUMENT, iBase_INVALID_ENTITY_TYPE, iBase_REGION, iBase_SUCCESS, iBase_VERTEX, KEEP_ARRAY, MB_SUCCESS, MOABI, RETURN, moab::Range::size(), and moab::TYPE_FROM_HANDLE().

  {
    CHKENUM(bridge_entity_type, iBase_EntityType, iBase_INVALID_ENTITY_TYPE);
    CHKENUM(requested_entity_type, iBase_EntityType, iBase_INVALID_ENTITY_TYPE);

    ErrorCode result = MB_SUCCESS;

    ALLOC_CHECK_ARRAY(offset, entity_handles_size+1);

    const EntityHandle* entity_iter = (const EntityHandle*)entity_handles;
    const EntityHandle* const entity_end = entity_iter + entity_handles_size;
    int* off_iter = *offset;
    int prev_off = 0;

    std::vector<EntityHandle> all_adj_ents;
    MeshTopoUtil mtu(MOABI);

    int min_bridge = iBase_VERTEX, max_bridge = iBase_REGION;
    int min_req    = iBase_VERTEX, max_req    = iBase_REGION;
    if (iBase_ALL_TYPES != bridge_entity_type)
      min_bridge = max_bridge = bridge_entity_type;
    if (iBase_ALL_TYPES != requested_entity_type)
      min_req = max_req = requested_entity_type;

    for ( ; entity_iter != entity_end; ++entity_iter)
    {
      *off_iter = prev_off;
      off_iter++;
      Range adj_ents;

      int source = CN::Dimension(TYPE_FROM_HANDLE(*entity_iter));
      for (int bridge = min_bridge; bridge <= max_bridge; ++bridge) {
        if (source == bridge)
          continue;
        for (int requested = min_req; requested <= max_req; ++requested) {
          if (bridge == requested)
            continue;
          result = mtu.get_bridge_adjacencies( *entity_iter,
                                               bridge,
                                               requested, adj_ents );
          CHKERR(result, "iMesh_getEntArr2ndAdj: trouble getting adjacency list.");
        }
      }

      std::copy(adj_ents.begin(), adj_ents.end(), std::back_inserter(all_adj_ents));
      prev_off += adj_ents.size();
    }
    *off_iter = prev_off;

    ALLOC_CHECK_ARRAY_NOFAIL(adj_entity_handles, all_adj_ents.size() );
    memcpy(*adj_entity_handles, &all_adj_ents[0],
           sizeof(EntityHandle)*all_adj_ents.size() );

    KEEP_ARRAY(offset);

      // Return an error if the bridge and requested entity types are different
    if (iBase_ALL_TYPES != bridge_entity_type &&
        bridge_entity_type == requested_entity_type)
      ERROR(iBase_INVALID_ARGUMENT, "iMesh_getEntArr2ndAdj: bridge and "
            "requested entity types must be different.");
    else
      RETURN(iBase_SUCCESS);
  }
void iMesh_getEntArrAdj ( iMesh_Instance  instance,
const iBase_EntityHandle entity_handles,
const int  entity_handles_size,
const int  entity_type_requested,
iBase_EntityHandle **  adjacentEntityHandles,
int *  adjacentEntityHandles_allocated,
int *  adj_entity_handles_size,
int **  offset,
int *  offset_allocated,
int *  offset_size,
int *  err 
)

Get entities of specified type adjacent to entities.

Get entities of specified type adjacent to entities. Specified type must be value in the iBase_EntityType enumeration. offset(i) is index of first entity in adjacentEntityHandles array adjacent to entity_handles[i]. More precisely, the entities adjacent to the ith entity in entity_handles are found in adjacentEntityHandles running from offset[i] to offset[i+1] - 1. This implies that the offset_size will be entity_handles_size + 1.

Note 1: Because 'adjacent' as defined by the iMesh data model refers to those entities that bound another, the entities being queried here (in entity_handles arg) are NEVER ALSO returned in adjacentEntityHandles even if the entity_type_requested matches the entity type(s) in entity_handles.

Parameters:
[in]instanceiMesh instance handle
[in]entity_handlesArray of entity handles being queried
[in]entity_handles_sizeNumber of entities in entity_handles array
[in]entity_type_requestedType of adjacent entities requested
[in,out]adjacentEntityHandlesPointer to array of adjacentEntityHandles Array pointer, allocated and occupied sizes argument trio) returned from function. Note that the implicit INTERLEAVED storage order rule applies (see iBase_StorageOrder)
[in,out]adjacentEntityHandles_allocatedPointer to allocated size of adjacentEntityHandles array
[out]adj_entity_handles_sizePointer to occupied size of adjacentEntityHandles array
[in,out]offsetPointer to array of offsets returned from function Array pointer, allocated and occupied sizes argument trio)
[in,out]offset_allocatedPointer to allocated size of offset array
[out]offset_sizePointer to occupied size of offset array
[out]errReturned Error status (see iBase_ErrorType)

Definition at line 626 of file iMesh_MOAB.cpp.

References ALLOC_CHECK_ARRAY, dim, ERROR, iBase_ALL_TYPES, iBase_BAD_ARRAY_SIZE, iBase_MEMORY_ALLOCATION_FAILED, iBase_SUCCESS, iBase_VERTEX, KEEP_ARRAY, MB_SUCCESS, MBPOLYHEDRON, MOABI, RETURN, moab::TYPE_FROM_HANDLE(), and smoab::UNION.

  {
    ErrorCode result = MB_SUCCESS;

    ALLOC_CHECK_ARRAY(offset, entity_handles_size+1);

    const EntityHandle* entity_iter = (const EntityHandle*)entity_handles;
    const EntityHandle* const entity_end = entity_iter + entity_handles_size;
    int* off_iter = *offset;
    int prev_off = 0;

    std::vector<EntityHandle> conn_storage;
    std::vector<EntityHandle> adj_ents;
    const EntityHandle *connect;
    int num_connect;

    EntityHandle* array; // ptr to working array of result handles
    int array_alloc;       // allocated size of 'array'
    const bool allocated_array = !*adjacentEntityHandles_allocated || !*adjacentEntityHandles;
    if (allocated_array) {
      array = 0;
      array_alloc = 0;
    }
    else {
      array = reinterpret_cast<EntityHandle*>(*adjacentEntityHandles);
      array_alloc = *adjacentEntityHandles_allocated;
    }

    for ( ; entity_iter != entity_end; ++entity_iter)
    {
      *off_iter = prev_off;
      off_iter++;

      if (iBase_VERTEX == entity_type_requested &&
          TYPE_FROM_HANDLE(*entity_iter) != MBPOLYHEDRON) {
        if (CN::Dimension(TYPE_FROM_HANDLE(*entity_iter)) == 0)
          continue;
        result = MOABI->get_connectivity(*entity_iter, connect, num_connect, false, &conn_storage);
        if (MB_SUCCESS != result) {
          if (allocated_array)
            free(array);
          ERROR(result, "iMesh_getEntArrAdj: trouble getting adjacency list.");
        }
      }
      else if (iBase_ALL_TYPES == entity_type_requested) {
        adj_ents.clear();
        for (int dim = 0; dim < 4; ++dim) {
          if (CN::Dimension(TYPE_FROM_HANDLE(*entity_iter)) == dim)
            continue;
          result = MOABI->get_adjacencies( entity_iter, 1, dim, false, adj_ents, Interface::UNION );
          if (MB_SUCCESS != result) {
            if (allocated_array)
              free(array);
            ERROR(result, "iMesh_getEntArrAdj: trouble getting adjacency list.");
          }
        }
        connect = &adj_ents[0];
        num_connect = adj_ents.size();
      }
      else {
        if (CN::Dimension(TYPE_FROM_HANDLE(*entity_iter)) == entity_type_requested)
          continue;
        adj_ents.clear();
        result = MOABI->get_adjacencies( entity_iter, 1,
                                       entity_type_requested, false, adj_ents );
        if (MB_SUCCESS != result) {
          if (allocated_array)
            free(array);
          ERROR(result, "iMesh_getEntArrAdj: trouble getting adjacency list.");
        }
        connect = &adj_ents[0];
        num_connect = adj_ents.size();
      }

      if (prev_off + num_connect <= array_alloc) {
        std::copy(connect, connect+num_connect, array+prev_off);
      }
      else if (allocated_array) {
          // if array is not allocated yet, guess at initial size
          // as the number of adjacencies for the first entity times
          // the number of input entities.  This will result in a single
          // exact allocation if one input entity or typical queries
          // such as connectivity of a non-mixed mesh or regions adjacent
          // to faces.
        if (!array_alloc)
          array_alloc = entity_handles_size * num_connect;
        else
          array_alloc = std::max(array_alloc*2, prev_off+num_connect);
        EntityHandle* new_array = (EntityHandle*)realloc( array, array_alloc*sizeof(EntityHandle) );
        if (!new_array) {
          free(array);
          RETURN(iBase_MEMORY_ALLOCATION_FAILED);
        }
        else
          array = new_array;
        std::copy(connect, connect+num_connect, array+prev_off);
      }
      // else do nothing.  Will catch error later when comparing
      //  occupied to allocated sizes.  Continue here because
      //  must pass back required size.

      prev_off += num_connect;
    }
    *off_iter = prev_off;
    *adjacentEntityHandles_size = prev_off;

    if (*adjacentEntityHandles_size > array_alloc) {
      if (allocated_array)
        free(array);
      RETURN(iBase_BAD_ARRAY_SIZE);
    }
    else if (allocated_array) {
      *adjacentEntityHandles = reinterpret_cast<iBase_EntityHandle*>(array);
      *adjacentEntityHandles_allocated = array_alloc;
    }

    KEEP_ARRAY(offset);
    RETURN(iBase_SUCCESS);
  }
void iMesh_setAdjTable ( iMesh_Instance  instance,
int *  adjacency_table,
int  adjacency_table_size,
int *  err 
)

Set the adjacency table as requested by the application.

Set the adjacency table as requested by the application. See iMesh_getAdjTable for a description of the meaning of the entries in this table. This call to set the adjacency table allows an application to request how it would like an implementation to behave when adjacent and/or intermediate entities are queried. If the implementation is not able to accommodate the requested query behavior associated with ANY entry in the table, the call will fail and return an error of iBase_NOT_SUPPORTED. Otherwise, the implementation is able to accommodate the requested query behavior associated with ALL table entries and the call will succeed. In either case, however, the implementation will over-write the entries in the adjacency_table argument with the same values that would be obtained in a succeeding call to iMesh_getAdjTable.

Parameters:
[in]instanceiMesh instance handle
[in,out]adjacency_tableArray representing adjacency table requested by application
[in]adjacency_table_sizeSize of adj table array
[out]errReturned Error status (see iBase_ErrorType)

Definition at line 425 of file iMesh_MOAB.cpp.

References iBase_INVALID_ARGUMENT, iBase_SUCCESS, iMesh_getGeometricDimension, MBIMESHI, munge_adj_table(), and RETURN.

  {
    if (16 != adj_table_size) {
      RETURN(iBase_INVALID_ARGUMENT);
    }

    int geom_dim;
    iMesh_getGeometricDimension(instance, &geom_dim, err);

    memcpy(MBIMESHI->AdjTable, adj_table, 16*sizeof(int));
    munge_adj_table(adj_table, geom_dim);
    RETURN(iBase_SUCCESS);
  }
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines