MOAB: Mesh Oriented datABase  (version 5.4.1)
MBMesquite::ArrayMesh Class Reference

#include <ArrayMesh.hpp>

+ Inheritance diagram for MBMesquite::ArrayMesh:
+ Collaboration diagram for MBMesquite::ArrayMesh:

Classes

struct  Tag

Public Member Functions

 ArrayMesh (int coords_per_vertex, unsigned long num_vertices, double *interleaved_vertex_coords, const int *vertex_fixed_flags, unsigned long num_elements, EntityTopology element_type, const unsigned long *element_connectivity_array, bool one_based_conn_indices=false, unsigned nodes_per_element=0, const int *vertex_slaved_flags=0)
 ArrayMesh (int coords_per_vertex, unsigned long num_vertices, double *interleaved_vertex_coords, const int *vertex_fixed_flags, unsigned long num_elements, const EntityTopology *element_types, const unsigned long *element_connectivity_array, const unsigned long *element_connectivity_offsets=NULL, bool one_based_conn_indices=false, const int *vertex_slaved_flags=0)
 ArrayMesh ()
 ~ArrayMesh ()
void clear_mesh ()
void set_mesh (int coords_per_vertex, unsigned long num_vertices, double *interleaved_vertex_coords, const int *vertex_fixed_flags, unsigned long num_elements, EntityTopology element_type, const unsigned long *element_connectivity_array, bool one_based_conn_indices=false, unsigned nodes_per_element=0, const int *vertex_slaved_flags=0)
TagHandle add_read_only_tag_data (const char *tag_name, TagType data_type, int vals_per_entity, const void *vertex_data, const void *element_data, const void *default_value, MsqError &err)
 Give mesquite access to per-entity application data via a tag.
TagHandle add_writable_tag_data (const char *tag_name, TagType tag_data_type, int vals_per_entity, void *vertex_data, void *element_data, const void *default_value, MsqError &err)
 Give mesquite access to per-entity application data via a tag.
virtual int get_geometric_dimension (MsqError &err)
 Returns whether this mesh lies in a 2D or 3D coordinate system.
virtual void get_all_elements (std::vector< ElementHandle > &elements, MsqError &err)
 Get all elements in mesh.
virtual void get_all_vertices (std::vector< VertexHandle > &vertices, MsqError &err)
 Get all vertices in mesh.
virtual VertexIteratorvertex_iterator (MsqError &err)
virtual ElementIteratorelement_iterator (MsqError &err)
virtual void vertices_get_fixed_flag (const VertexHandle vert_array[], std::vector< bool > &fixed_flag_array, size_t num_vtx, MsqError &err)
virtual void vertices_get_slaved_flag (const VertexHandle vert_array[], std::vector< bool > &slaved_flag_array, size_t num_vtx, MsqError &err)
virtual void vertices_get_coordinates (const VertexHandle vert_array[], MsqVertex *coordinates, size_t num_vtx, MsqError &err)
 Get/set location of a vertex.
virtual void vertex_set_coordinates (VertexHandle vertex, const Vector3D &coordinates, MsqError &err)
virtual void vertex_set_byte (VertexHandle vertex, unsigned char byte, MsqError &err)
virtual void vertices_set_byte (const VertexHandle *vert_array, const unsigned char *byte_array, size_t array_size, MsqError &err)
virtual void vertex_get_byte (const VertexHandle vertex, unsigned char *byte, MsqError &err)
virtual void vertices_get_byte (const VertexHandle *vertex, unsigned char *byte_array, size_t array_size, MsqError &err)
virtual void vertices_get_attached_elements (const VertexHandle *vertex_array, size_t num_vertex, std::vector< ElementHandle > &elements, std::vector< size_t > &offsets, MsqError &err)
 get elements adjacent to vertices
virtual void elements_get_attached_vertices (const ElementHandle *elem_handles, size_t num_elems, std::vector< VertexHandle > &vert_handles, std::vector< size_t > &offsets, MsqError &err)
 Get element connectivity.
virtual void elements_get_topologies (const ElementHandle *element_handle_array, EntityTopology *element_topologies, size_t num_elements, MsqError &err)
virtual TagHandle tag_create (const std::string &tag_name, TagType type, unsigned length, const void *default_value, MsqError &err)
 Create a tag.
virtual void tag_delete (TagHandle handle, MsqError &err)
 Remove a tag and all corresponding data.
virtual TagHandle tag_get (const std::string &name, MsqError &err)
 Get handle for existing tag, by name.
virtual void tag_properties (TagHandle handle, std::string &name_out, TagType &type_out, unsigned &length_out, MsqError &err)
 Get properites of tag.
virtual void tag_set_element_data (TagHandle handle, size_t num_elems, const ElementHandle *elem_array, const void *tag_data, MsqError &err)
 Set tag values on elements.
virtual void tag_set_vertex_data (TagHandle handle, size_t num_elems, const VertexHandle *node_array, const void *tag_data, MsqError &err)
 Set tag values on vertices.
virtual void tag_get_element_data (TagHandle handle, size_t num_elems, const ElementHandle *elem_array, void *tag_data, MsqError &err)
 Get tag values on elements.
virtual void tag_get_vertex_data (TagHandle handle, size_t num_elems, const VertexHandle *node_array, void *tag_data, MsqError &err)
 Get tag values on vertices.
virtual void release_entity_handles (const EntityHandle *handle_array, size_t num_handles, MsqError &err)
virtual void release ()

Private Member Functions

const unsigned long * elem_verts (size_t elem_index, int &num_vertex) const
void build_vertex_adjacency_list ()
Tagallocate_tag (const char *name, bool owned, TagType type, unsigned size, const void *vertex_ro_data, void *vertex_rw_data, const void *element_ro_data, void *element_rw_data, const void *default_value, MsqError &err)
bool valid () const

Static Private Member Functions

static unsigned bytes (TagType type)
static void fill (unsigned char *buffer, const unsigned char *value, size_t size, size_t count)

Private Attributes

int mDimension
 Coordinates per vertex.
unsigned long vertexCount
 Number of vertices.
double * coordArray
 Interleaved vertex coordinates.
const int * fixedFlags
 Vertex fixed flags.
const int * slavedFlags
 Vertex slaved flags.
unsigned char * vertexByteArray
 Vertex bytes.
unsigned long elementCount
 Number of elements.
const unsigned long * connArray
 Element connectivity.
const unsigned long * connOffsets
unsigned long * allocConnOffsets
EntityTopology elementType
 Type for all elements if connOffsets is NULL.
const EntityTopologyelementTypes
 Type for each element type if connOffsets is not NULL.
unsigned nodesPerElement
 Nodes per element if connOffsets is NULL.
bool oneBasedArrays
 FORTRAN-style array indexing.
unsigned long * vertexAdjacencyList
unsigned long * vertexAdjacencyOffsets
TagtagList

Detailed Description

Definition at line 41 of file ArrayMesh.hpp.


Constructor & Destructor Documentation

MBMesquite::ArrayMesh::ArrayMesh ( int  coords_per_vertex,
unsigned long  num_vertices,
double *  interleaved_vertex_coords,
const int *  vertex_fixed_flags,
unsigned long  num_elements,
EntityTopology  element_type,
const unsigned long *  element_connectivity_array,
bool  one_based_conn_indices = false,
unsigned  nodes_per_element = 0,
const int *  vertex_slaved_flags = 0 
)

Create a MBMesquite::Mesh instance that wraps application-provided arrays.

Note: An instance of this class will reference the passed in arrays. It will not copy them.

Parameters:
coords_per_vertexDimension of the mesh (2 or 3)
num_verticesNumber of vertices in the mesh
interleaved_vertex_coordsVertex coordinates. Ordered as [x1, y1, z1, x2, y2, z2, ...]
vertex_fixed_flagsOne value per vertex. Zero if vertex is free, one if the position is fixed.
num_elementsNumber of elements in the mesh
element_typeThe type of the elements
element_connectivity_arrayElement connectivity, specified as vertex indices such that the location of the vertex coordinates in vertex_coords is at 3 times the value in this array.
one_based_conn_indicesUse one-based (Fortran) array indexing.
nodes_per_elementNumber of nodes in each element. If not specified, number of nodes in a linear element with the type 'element_type' is assumed.
vertex_slaved_flagsOne value per vertex. Zero if vertex is free, one if the vertex is slaved to the logical position of the element mapping/shape function..

Definition at line 77 of file ArrayMesh.cpp.

References coordArray, corners, fixedFlags, mDimension, nodesPerElement, oneBasedArrays, valid(), and vertexByteArray.

    : mDimension( coords_per_vertex ), vertexCount( num_vertices ), coordArray( interleaved_vertex_coords ),
      fixedFlags( vertex_fixed_flags ), slavedFlags( vertex_slaved_flags ),
      vertexByteArray( new unsigned char[num_vertices + one_based_conn_indices] ), elementCount( num_elements ),
      connArray( element_connectivity_array ), connOffsets( 0 ), allocConnOffsets( 0 ), elementType( element_type ),
      elementTypes( 0 ), nodesPerElement( nodes_per_element ), oneBasedArrays( one_based_conn_indices ),
      vertexAdjacencyList( 0 ), vertexAdjacencyOffsets( 0 ), tagList( 0 )
{
    if( oneBasedArrays )
    {
        coordArray -= mDimension;
        --fixedFlags;
    }

    if( nodesPerElement < 2 ) nodesPerElement = TopologyInfo::corners( element_type );

    assert( valid() );
    memset( vertexByteArray, 0, num_vertices + one_based_conn_indices );
}
MBMesquite::ArrayMesh::ArrayMesh ( int  coords_per_vertex,
unsigned long  num_vertices,
double *  interleaved_vertex_coords,
const int *  vertex_fixed_flags,
unsigned long  num_elements,
const EntityTopology element_types,
const unsigned long *  element_connectivity_array,
const unsigned long *  element_connectivity_offsets = NULL,
bool  one_based_conn_indices = false,
const int *  vertex_slaved_flags = 0 
)

Create a MBMesquite::Mesh instance that wraps application-provided arrays.

Note: An instance of this class will reference the passed in arrays. It will not copy them.

Parameters:
coords_per_vertexDimension of the mesh (2 or 3)
num_verticesNumber of vertices in the mesh
interleaved_vertex_coordsVertex coordinates. Ordered as [x1, y1, z1, x2, y2, z2, ...]
vertex_fixed_flagsOne value per vertex. Zero if vertex is free, one if the poistion is fixed.
num_elementsNumber of elements in the mesh
element_typesThe topological type for each element.
element_connectivity_arrayElement connectivity, specified as vertex indices such that the location of the vertex coordinates in vertex_coords is at 3 times the value in this array.
element_connectivity_offsetsAn optional array of length one greater than num_elements. Each entry other than the last should contain the value of the index into element_connectivity_array at which the connectivity data for the corresponding element begins. The last entry in the array must be the next-to-last entry plus the length of the connectivity data for the last element, such that the length of the connectivity data for any element 'e' can be calculated with: n = element_connectivity_offsets[e+1] - element_connectivity_offsets[e]

If this array is not specified, then it will be assumed that the length of the connectivity data for each element is the number of corners necessary to represent its topological type (that it has no higher-order nodes.)

Parameters:
one_based_conn_indicesUse one-based (Fortran) array indexing.
vertex_slaved_flagsOne value per vertex. Zero if vertex is free, one if the vertex is slaved to the logical position of the element mapping/shape function..

Definition at line 106 of file ArrayMesh.cpp.

References allocConnOffsets, connArray, connOffsets, coordArray, corners, elementTypes, fixedFlags, mDimension, oneBasedArrays, valid(), and vertexByteArray.

    : mDimension( coords_per_vertex ), vertexCount( num_vertices ), coordArray( interleaved_vertex_coords ),
      fixedFlags( vertex_fixed_flags ), slavedFlags( vertex_slaved_flags ),
      vertexByteArray( new unsigned char[num_vertices + one_based_conn_indices] ), elementCount( num_elements ),
      connArray( element_connectivity_array ), connOffsets( element_connectivity_offsets ), allocConnOffsets( 0 ),
      elementType( MIXED ), elementTypes( element_types ), nodesPerElement( 0 ),
      oneBasedArrays( one_based_conn_indices ), vertexAdjacencyList( 0 ), vertexAdjacencyOffsets( 0 ), tagList( 0 )
{
    if( oneBasedArrays )
    {
        coordArray -= mDimension;
        --fixedFlags;
        if( element_connectivity_offsets ) --connArray;
    }

    if( !element_connectivity_offsets )
    {
        connOffsets = allocConnOffsets = new unsigned long[num_elements + 1];
        allocConnOffsets[0]            = 0;
        for( unsigned long i = 1; i <= num_elements; ++i )
            allocConnOffsets[i] = allocConnOffsets[i - 1] + TopologyInfo::corners( elementTypes[i - 1] );
    }

    assert( valid() );
    memset( vertexByteArray, 0, num_vertices + one_based_conn_indices );
}

Definition at line 241 of file ArrayMesh.cpp.

References clear_mesh().

{
    clear_mesh();
}

Member Function Documentation

TagHandle MBMesquite::ArrayMesh::add_read_only_tag_data ( const char *  tag_name,
TagType  data_type,
int  vals_per_entity,
const void *  vertex_data,
const void *  element_data,
const void *  default_value,
MsqError err 
)

Give mesquite access to per-entity application data via a tag.

Allow mesquite to access application-owned per-element and/or per-vertex data by assigning a tag name to said data. Mesquite will be allowed to read application-owned data. Any attempt by Mesquite to modify the data will result in an internal error condition.

Parameters:
tag_nameTag name through which Mesquite can access the data
data_typeThe type of the data
vals_per_entityNumber of values of type data_type that each entity has Data is assumed to be interleaved such that all values associated with a a single entity are adjacent in memory
vertex_dataA pointer to the region of memory in which per-vertex values are stored. May be NULL if per-vertex data is not available. Values are assumed to be in the same order as the vertex coordinates passed to the constructor or set_mesh.
vertex_dataA pointer to the region of memory in which per-element values are stored. May be NULL if per-element data is not available. Values are assumed to be in the same order as the element connetivity passed to the constructor or set_mesh.
default_valueValue to return for all vertices and/or all elements if the cooresponding data array is null. May be NULL if no default value.

Definition at line 597 of file ArrayMesh.cpp.

References allocate_tag(), and MSQ_ERRZERO.

Referenced by ArrayMeshTest::test_delete_tag(), and ArrayMeshTest::test_tag_data().

{
    Tag* tag = allocate_tag( tag_name, false, data_type, vals_per_entity, vertex_data, 0, element_data, 0,
                             default_value, err );
    MSQ_ERRZERO( err );
    return reinterpret_cast< TagHandle >( tag );
}
TagHandle MBMesquite::ArrayMesh::add_writable_tag_data ( const char *  tag_name,
TagType  tag_data_type,
int  vals_per_entity,
void *  vertex_data,
void *  element_data,
const void *  default_value,
MsqError err 
)

Give mesquite access to per-entity application data via a tag.

Allow mesquite to access and/or set application-owned per-element and/or per-vertex data by assigning a tag name to said data.

Parameters:
tag_nameTag name through which Mesquite can access the data
data_typeThe type of the data
vals_per_entityNumber of values of type data_type that each entity has Data is assumed to be interleaved such that all values associated with a a single entity are adjacent in memory
vertex_dataA pointer to the region of memory in which per-vertex values are stored. May be NULL if per-vertex data is not available. Values are assumed to be in the same order as the vertex coordinates passed to the constructor or set_mesh.
vertex_dataA pointer to the region of memory in which per-element values are stored. May be NULL if per-element data is not available. Values are assumed to be in the same order as the element connetivity passed to the constructor or set_mesh.
default_valueValue to return for all vertices and/or all elements if the cooresponding data array is null. May be NULL if no default value.

Definition at line 611 of file ArrayMesh.cpp.

References allocate_tag(), and MSQ_ERRZERO.

Referenced by ArrayMeshTest::test_delete_tag(), and ArrayMeshTest::test_tag_data().

{
    Tag* tag = allocate_tag( tag_name, false, data_type, vals_per_entity, 0, vertex_data, 0, element_data,
                             default_value, err );
    MSQ_ERRZERO( err );
    return reinterpret_cast< TagHandle >( tag );
}
ArrayMesh::Tag * MBMesquite::ArrayMesh::allocate_tag ( const char *  name,
bool  owned,
TagType  type,
unsigned  size,
const void *  vertex_ro_data,
void *  vertex_rw_data,
const void *  element_ro_data,
void *  element_rw_data,
const void *  default_value,
MsqError err 
) [private]

Definition at line 537 of file ArrayMesh.cpp.

References bytes(), MBMesquite::ArrayMesh::Tag::defaultValue, MBMesquite::ArrayMesh::Tag::eleReadPtr, MBMesquite::ArrayMesh::Tag::eleWritePtr, MSQ_SETERR, MBMesquite::ArrayMesh::Tag::name, MBMesquite::ArrayMesh::Tag::next, MBMesquite::ArrayMesh::Tag::owned, MBMesquite::ArrayMesh::Tag::size, MBMesquite::MsqError::TAG_ALREADY_EXISTS, tagList, MBMesquite::ArrayMesh::Tag::type, MBMesquite::ArrayMesh::Tag::vtxReadPtr, and MBMesquite::ArrayMesh::Tag::vtxWritePtr.

Referenced by add_read_only_tag_data(), add_writable_tag_data(), and tag_create().

{
    // check if name is already in use
    for( Tag* iter = tagList; iter; iter = iter->next )
    {
        if( !strcmp( iter->name, name ) )
        {
            MSQ_SETERR( err )( MsqError::TAG_ALREADY_EXISTS );
            return 0;
        }
    }

    // allocate object
    Tag* result = new Tag;

    // initialize members
    result->type  = type;
    result->size  = size * bytes( type );
    result->owned = owned;
    result->name  = new char[strlen( name ) + 1];
    strcpy( result->name, name );

    result->vtxWritePtr = reinterpret_cast< unsigned char* >( vertex_rw_data );
    if( vertex_rw_data )
        result->vtxReadPtr = reinterpret_cast< unsigned char* >( vertex_rw_data );
    else
        result->vtxReadPtr = reinterpret_cast< const unsigned char* >( vertex_ro_data );

    result->eleWritePtr = reinterpret_cast< unsigned char* >( element_rw_data );
    if( element_rw_data )
        result->eleReadPtr = reinterpret_cast< unsigned char* >( element_rw_data );
    else
        result->eleReadPtr = reinterpret_cast< const unsigned char* >( element_ro_data );

    if( default_value )
    {
        result->defaultValue = new unsigned char[result->size];
        memcpy( result->defaultValue, default_value, result->size );
    }
    else
    {
        result->defaultValue = 0;
    }

    // prepend to tag list
    result->next = tagList;
    tagList      = result;

    return result;
}

Definition at line 470 of file ArrayMesh.cpp.

References conn, elem_verts(), elementCount, n, oneBasedArrays, vertexAdjacencyList, vertexAdjacencyOffsets, and vertexCount.

Referenced by vertices_get_attached_elements().

{
    delete[] vertexAdjacencyList;
    delete[] vertexAdjacencyOffsets;
    vertexAdjacencyOffsets = new unsigned long[vertexCount + oneBasedArrays + 1];

    // for each vertex, store the number of elements the previous
    // vertex occurs in.
    memset( vertexAdjacencyOffsets, 0, sizeof( unsigned long ) * ( vertexCount + oneBasedArrays + 1 ) );
    for( size_t i = 0; i < elementCount; ++i )
    {
        int n;
        const unsigned long* conn = elem_verts( i, n );
        for( int j = 0; j < n; ++j )
            ++vertexAdjacencyOffsets[conn[j] + 1];
    }

    // convert vertexAdjacencyOffsets from a shifted list of counts
    // to a list of offsts
    for( size_t i = 1; i <= vertexCount + oneBasedArrays; ++i )
        vertexAdjacencyOffsets[i] += vertexAdjacencyOffsets[i - 1];

    // allocate space and populate with reverse connectivity
    vertexAdjacencyList = new unsigned long[vertexAdjacencyOffsets[vertexCount + oneBasedArrays]];
    for( size_t i = 0; i < elementCount; ++i )
    {
        int n;
        const unsigned long* conn = elem_verts( i, n );
        for( int j = 0; j < n; ++j )
            vertexAdjacencyList[vertexAdjacencyOffsets[conn[j]]++] = i;
    }

    for( size_t i = vertexCount + oneBasedArrays; i > 0; --i )
        vertexAdjacencyOffsets[i] = vertexAdjacencyOffsets[i - 1];
    vertexAdjacencyOffsets[0] = 0;
}
unsigned MBMesquite::ArrayMesh::bytes ( TagType  type) [static, private]

Definition at line 507 of file ArrayMesh.cpp.

References MBMesquite::Mesh::BOOL, MBMesquite::Mesh::BYTE, MBMesquite::Mesh::DOUBLE, MBMesquite::Mesh::HANDLE, and MBMesquite::Mesh::INT.

Referenced by allocate_tag(), and tag_properties().

{
    switch( type )
    {
        case BYTE:
        case BOOL:
            return 1;
        case INT:
            return sizeof( int );
        case DOUBLE:
            return sizeof( double );
        case HANDLE:
            return sizeof( EntityHandle );
    }
    return 0;
}
const unsigned long * MBMesquite::ArrayMesh::elem_verts ( size_t  elem_index,
int &  num_vertex 
) const [inline, private]

Definition at line 246 of file ArrayMesh.cpp.

References connArray, connOffsets, elementCount, and nodesPerElement.

Referenced by build_vertex_adjacency_list(), and elements_get_attached_vertices().

{
    assert( e < elementCount );
    if( connOffsets )
    {
        n = connOffsets[e + 1] - connOffsets[e];
        return connArray + connOffsets[e];
    }
    else
    {
        n = nodesPerElement;
        return connArray + nodesPerElement * e;
    }
}

Definition at line 285 of file ArrayMesh.cpp.

References elementCount.

{
    return new IndexIterator( 0, elementCount );
}
void MBMesquite::ArrayMesh::elements_get_attached_vertices ( const ElementHandle elem_handles,
size_t  num_elems,
std::vector< VertexHandle > &  vert_handles,
std::vector< size_t > &  offsets,
MsqError err 
) [virtual]

Get element connectivity.

Get the connectivity (ordered list of vertex handles) for each element in the input array.

Parameters:
elem_handlesThe array of element handles for which to retrieve the connectivity list.
num_elemsThe length of #elem_handles
vert_handlesArray in which to place the vertex handles in each elements connectivity.
offsetsFor each element in #elem_handles, the value in the same position in this array is the index into #vert_handles at which the connectivity list for that element begins.

Implements MBMesquite::Mesh.

Definition at line 421 of file ArrayMesh.cpp.

References conn, elem_verts(), and elementCount.

{
    const size_t* indices = (const size_t*)elem_handles;
    offsets.resize( num_elems + 1 );
    vert_handles.clear();
    for( size_t i = 0; i < num_elems; ++i )
    {
        assert( indices[i] < elementCount );
        int count;
        const unsigned long* conn = elem_verts( indices[i], count );
        size_t prev_size          = vert_handles.size();
        offsets[i]                = prev_size;
        vert_handles.resize( prev_size + count );
        std::copy( conn, conn + count, (size_t*)( &vert_handles[prev_size] ) );
    }
    offsets[num_elems] = vert_handles.size();
}
void MBMesquite::ArrayMesh::elements_get_topologies ( const ElementHandle element_handle_array,
EntityTopology element_topologies,
size_t  num_elements,
MsqError err 
) [virtual]

Returns the topologies of the given entities. The "entity_topologies" array must be at least "num_elements" in size.

Implements MBMesquite::Mesh.

Definition at line 443 of file ArrayMesh.cpp.

References elementCount, elementType, elementTypes, and MBMesquite::MIXED.

{
    const size_t* indices = (const size_t*)handles;
    if( elementType == MIXED )
        for( size_t i = 0; i < num_elements; ++i )
        {
            assert( indices[i] < elementCount );
            element_topologies[i] = elementTypes[indices[i]];
        }
    else
        for( size_t i = 0; i < num_elements; ++i )
        {
            assert( indices[i] < elementCount );
            element_topologies[i] = elementType;
        }
}
void MBMesquite::ArrayMesh::fill ( unsigned char *  buffer,
const unsigned char *  value,
size_t  size,
size_t  count 
) [static, private]

Definition at line 524 of file ArrayMesh.cpp.

References size.

Referenced by tag_get_element_data(), tag_get_vertex_data(), tag_set_element_data(), and tag_set_vertex_data().

{
    if( !value )
    {
        memset( buffer, 0, size * count );
        return;
    }

    unsigned char* const end = buffer + size * count;
    for( unsigned char* iter = buffer; iter != end; iter += size )
        memcpy( iter, value, size );
}
void MBMesquite::ArrayMesh::get_all_elements ( std::vector< ElementHandle > &  elements,
MsqError err 
) [virtual]

Get all elements in mesh.

Get the handles of every element in the active mesh.

Implements MBMesquite::Mesh.

Definition at line 266 of file ArrayMesh.cpp.

References elementCount.

Referenced by ArrayMeshTest::test_tag_data(), and HigherOrderTest::test_tri_open_domain().

{
    elements.resize( elementCount );
    for( unsigned long i = 0; i < elementCount; ++i )
        elements[i] = (Mesh::ElementHandle)i;
}
void MBMesquite::ArrayMesh::get_all_vertices ( std::vector< VertexHandle > &  vertices,
MsqError err 
) [virtual]

Get all vertices in mesh.

Get the handles of every vertex in the active mesh

Implements MBMesquite::Mesh.

Definition at line 273 of file ArrayMesh.cpp.

References oneBasedArrays, and vertexCount.

Referenced by TerminationCriterionTest::test_abs_vtx_movement_culling(), ArrayMeshTest::test_tag_data(), and HigherOrderTest::test_tri_open_domain().

{
    vertices.resize( vertexCount );
    for( unsigned long i = 0; i < vertexCount; ++i )
        vertices[i] = ( Mesh::VertexHandle )( i + oneBasedArrays );
}

Returns whether this mesh lies in a 2D or 3D coordinate system.

Implements MBMesquite::Mesh.

Definition at line 261 of file ArrayMesh.cpp.

References mDimension.

{
    return mDimension;
}
void MBMesquite::ArrayMesh::release ( ) [virtual]

Instead of deleting a Mesh when you think you are done, call release(). In simple cases, the implementation could just call the destructor. More sophisticated implementations may want to keep the Mesh object to live longer than Mesquite is using it.

Implements MBMesquite::Mesh.

Definition at line 468 of file ArrayMesh.cpp.

{}
void MBMesquite::ArrayMesh::release_entity_handles ( const EntityHandle handle_array,
size_t  num_handles,
MsqError err 
) [virtual]

Tells the mesh that the client is finished with a given entity handle.

Implements MBMesquite::Mesh.

Definition at line 463 of file ArrayMesh.cpp.

References MSQ_SETERR, and NOT_IMPLEMENTED.

void MBMesquite::ArrayMesh::set_mesh ( int  coords_per_vertex,
unsigned long  num_vertices,
double *  interleaved_vertex_coords,
const int *  vertex_fixed_flags,
unsigned long  num_elements,
EntityTopology  element_type,
const unsigned long *  element_connectivity_array,
bool  one_based_conn_indices = false,
unsigned  nodes_per_element = 0,
const int *  vertex_slaved_flags = 0 
)

Definition at line 203 of file ArrayMesh.cpp.

References clear_mesh(), connArray, coordArray, corners, elementCount, elementType, fixedFlags, mDimension, nodesPerElement, oneBasedArrays, slavedFlags, valid(), vertexByteArray, and vertexCount.

{
    clear_mesh();
    mDimension     = coords_per_vertex;
    vertexCount    = num_vertices;
    coordArray     = interleaved_vertex_coords;
    fixedFlags     = vertex_fixed_flags;
    slavedFlags    = vertex_slaved_flags;
    elementCount   = num_elements;
    connArray      = element_connectivity_array;
    elementType    = element_type;
    oneBasedArrays = one_based_conn_indices;

    if( oneBasedArrays )
    {
        coordArray -= mDimension;
        --fixedFlags;
    }

    if( nodes_per_element < 2 )
        nodesPerElement = TopologyInfo::corners( element_type );
    else
        nodesPerElement = nodes_per_element;

    vertexByteArray = new unsigned char[num_vertices + one_based_conn_indices];
    assert( valid() );
    memset( vertexByteArray, 0, num_vertices + one_based_conn_indices );
}
TagHandle MBMesquite::ArrayMesh::tag_create ( const std::string &  tag_name,
TagType  type,
unsigned  length,
const void *  default_value,
MsqError err 
) [virtual]

Create a tag.

Create a user-defined data type that can be attached to any element or vertex in the mesh. For an opaque or undefined type, use type=BYTE and length=sizeof(..).

Parameters:
tag_nameA unique name for the data object
typeThe type of the data
lengthNumber of values per entity (1->scalar, >1 ->vector)
default_valueDefault value to assign to all entities - may be NULL
Returns:
- Handle for tag definition

Implements MBMesquite::Mesh.

Definition at line 625 of file ArrayMesh.cpp.

References allocate_tag(), MSQ_ERRZERO, and size.

Referenced by ArrayMeshTest::test_delete_tag(), and ArrayMeshTest::test_tag_data().

{
    Tag* tag = allocate_tag( tag_name.c_str(), true, data_type, size, 0, 0, 0, 0, default_value, err );
    MSQ_ERRZERO( err );
    return reinterpret_cast< TagHandle >( tag );
}
void MBMesquite::ArrayMesh::tag_delete ( TagHandle  handle,
MsqError err 
) [virtual]

Remove a tag and all corresponding data.

Delete a tag.

Implements MBMesquite::Mesh.

Definition at line 636 of file ArrayMesh.cpp.

References MBMesquite::ArrayMesh::Tag::defaultValue, MBMesquite::ArrayMesh::Tag::eleWritePtr, MSQ_SETERR, MBMesquite::ArrayMesh::Tag::name, MBMesquite::ArrayMesh::Tag::next, MBMesquite::ArrayMesh::Tag::owned, MBMesquite::MsqError::TAG_NOT_FOUND, tagList, and MBMesquite::ArrayMesh::Tag::vtxWritePtr.

Referenced by ArrayMeshTest::test_delete_tag().

{
    Tag* ptr = reinterpret_cast< Tag* >( handle );
    // find previous tag pointer in list
    if( !tagList )
    {
        MSQ_SETERR( err )( "Invalid tag handle", MsqError::TAG_NOT_FOUND );
        return;
    }
    Tag* prev = 0;
    for( prev = tagList; prev && prev->next != ptr; prev = prev->next )
        ;
    if( !prev && tagList != ptr )
    {
        MSQ_SETERR( err )( "Invalid tag handle", MsqError::TAG_NOT_FOUND );
        return;
    }
    delete[] ptr->name;
    delete[] ptr->defaultValue;
    if( ptr->owned )
    {
        delete[] ptr->vtxWritePtr;
        delete[] ptr->eleWritePtr;
    }
    if( prev )
        prev->next = ptr->next;
    else
        tagList = ptr->next;
    delete ptr;
}
TagHandle MBMesquite::ArrayMesh::tag_get ( const std::string &  name,
MsqError err 
) [virtual]

Get handle for existing tag, by name.

Check for the existance of a tag given it's name and if it exists return a handle for it. If the specified tag does not exist, zero should be returned WITHOUT flagging an error.

Implements MBMesquite::Mesh.

Definition at line 667 of file ArrayMesh.cpp.

References MSQ_SETERR, MBMesquite::ArrayMesh::Tag::next, MBMesquite::MsqError::TAG_NOT_FOUND, and tagList.

Referenced by ArrayMeshTest::test_delete_tag(), and ArrayMeshTest::test_tag_data().

{
    for( Tag* iter = tagList; iter; iter = iter->next )
        if( name == iter->name ) return reinterpret_cast< TagHandle >( iter );
    MSQ_SETERR( err )( MsqError::TAG_NOT_FOUND );
    return 0;
}
void MBMesquite::ArrayMesh::tag_get_element_data ( TagHandle  handle,
size_t  num_elems,
const ElementHandle elem_array,
void *  tag_data,
MsqError err 
) [virtual]

Get tag values on elements.

Get the value of a tag for a list of mesh elements.

Parameters:
handleThe tag
num_elemsLength of elem_array
elem_arrayArray of elements for which to get the tag value.
tag_dataReturn buffer in which to copy tag data, contiguous in memory. This data is expected to be num_elems*tag_length*sizeof(tag_type) bytes.

Implements MBMesquite::Mesh.

Definition at line 749 of file ArrayMesh.cpp.

References MBMesquite::ArrayMesh::Tag::defaultValue, MBMesquite::ArrayMesh::Tag::eleReadPtr, fill(), MSQ_SETERR, MBMesquite::ArrayMesh::Tag::size, and MBMesquite::MsqError::TAG_NOT_FOUND.

Referenced by ArrayMeshTest::test_tag_data().

{
    unsigned char* ptr = reinterpret_cast< unsigned char* >( data );
    const Tag* tag     = reinterpret_cast< const Tag* >( handle );
    if( tag->eleReadPtr )
    {
        for( size_t i = 0; i < count; ++i )
        {
            size_t idx = reinterpret_cast< size_t >( entities[i] );
            memcpy( ptr + i * tag->size, tag->eleReadPtr + idx * tag->size, tag->size );
        }
    }
    else if( tag->defaultValue )
    {
        fill( ptr, tag->defaultValue, tag->size, count );
    }
    else
    {
        MSQ_SETERR( err )( MsqError::TAG_NOT_FOUND );
    }
}
void MBMesquite::ArrayMesh::tag_get_vertex_data ( TagHandle  handle,
size_t  num_elems,
const VertexHandle node_array,
void *  tag_data,
MsqError err 
) [virtual]

Get tag values on vertices.

Get the value of a tag for a list of mesh vertices.

Parameters:
handleThe tag
num_elemsLength of elem_array
elem_arrayArray of vertices for which to get the tag value.
tag_dataReturn buffer in which to copy tag data, contiguous in memory. This data is expected to be num_elems*tag_length*sizeof(tag_type) bytes.

Implements MBMesquite::Mesh.

Definition at line 775 of file ArrayMesh.cpp.

References MBMesquite::ArrayMesh::Tag::defaultValue, fill(), MSQ_SETERR, oneBasedArrays, MBMesquite::ArrayMesh::Tag::size, MBMesquite::MsqError::TAG_NOT_FOUND, and MBMesquite::ArrayMesh::Tag::vtxReadPtr.

Referenced by ArrayMeshTest::test_tag_data().

{
    unsigned char* ptr = reinterpret_cast< unsigned char* >( data );
    const Tag* tag     = reinterpret_cast< const Tag* >( handle );
    if( tag->vtxReadPtr )
    {
        for( size_t i = 0; i < count; ++i )
        {
            size_t idx = reinterpret_cast< size_t >( entities[i] ) - oneBasedArrays;
            memcpy( ptr + i * tag->size, tag->vtxReadPtr + idx * tag->size, tag->size );
        }
    }
    else if( tag->defaultValue )
    {
        fill( ptr, tag->defaultValue, tag->size, count );
    }
    else
    {
        MSQ_SETERR( err )( MsqError::TAG_NOT_FOUND );
    }
}
void MBMesquite::ArrayMesh::tag_properties ( TagHandle  handle,
std::string &  name_out,
TagType type_out,
unsigned &  length_out,
MsqError err 
) [virtual]

Get properites of tag.

Get data type and number of values per entity for tag.

Parameters:
handleTag to get properties of.
name_outPassed back tag name.
type_outPassed back tag type.
length_outPassed back number of values per entity.

Implements MBMesquite::Mesh.

Definition at line 675 of file ArrayMesh.cpp.

References bytes(), MBMesquite::ArrayMesh::Tag::name, MBMesquite::ArrayMesh::Tag::size, and MBMesquite::ArrayMesh::Tag::type.

Referenced by ArrayMeshTest::test_tag_data().

{
    const Tag* ptr = reinterpret_cast< const Tag* >( handle );
    name           = ptr->name;
    type           = ptr->type;
    size           = ptr->size / bytes( ptr->type );
}
void MBMesquite::ArrayMesh::tag_set_element_data ( TagHandle  handle,
size_t  num_elems,
const ElementHandle elem_array,
const void *  tag_data,
MsqError err 
) [virtual]

Set tag values on elements.

Set the value of a tag for a list of mesh elements.

Parameters:
handleThe tag
num_elemsLength of elem_array
elem_arrayArray of elements for which to set the tag value.
tag_dataTag data for each element, contiguous in memory. This data is expected to be num_elems*tag_length*sizeof(tag_type) bytes.

Implements MBMesquite::Mesh.

Definition at line 683 of file ArrayMesh.cpp.

References MBMesquite::ArrayMesh::Tag::defaultValue, elementCount, MBMesquite::ArrayMesh::Tag::eleReadPtr, MBMesquite::ArrayMesh::Tag::eleWritePtr, fill(), MSQ_SETERR, MBMesquite::ArrayMesh::Tag::owned, MBMesquite::ArrayMesh::Tag::size, and MBMesquite::MsqError::TAG_ALREADY_EXISTS.

Referenced by ArrayMeshTest::test_tag_data().

{
    Tag* tag = reinterpret_cast< Tag* >( handle );
    if( !tag->eleWritePtr )
    {
        if( !tag->owned )
        {
            MSQ_SETERR( err )
            ( "Attempt to set non-writeable (application owned) tag data", MsqError::TAG_ALREADY_EXISTS );
            return;
        }
        else
        {
            assert( !tag->eleReadPtr );
            tag->eleReadPtr = tag->eleWritePtr = new unsigned char[elementCount * tag->size];
            fill( tag->eleWritePtr, tag->defaultValue, tag->size, elementCount );
        }
    }
    const unsigned char* ptr = reinterpret_cast< const unsigned char* >( data );
    // as we're working with user-supplied arrays, make sure we're not
    // memcpying overlapping regions
    assert( ptr + tag->size * elementCount <= tag->eleWritePtr || tag->eleWritePtr + tag->size * elementCount <= ptr );
    for( size_t i = 0; i < count; ++i )
    {
        size_t idx = reinterpret_cast< size_t >( entities[i] );
        memcpy( tag->eleWritePtr + idx * tag->size, ptr + i * tag->size, tag->size );
    }
}
void MBMesquite::ArrayMesh::tag_set_vertex_data ( TagHandle  handle,
size_t  num_elems,
const VertexHandle node_array,
const void *  tag_data,
MsqError err 
) [virtual]

Set tag values on vertices.

Set the value of a tag for a list of mesh vertices.

Parameters:
handleThe tag
num_elemsLength of node_array
node_arrayArray of vertices for which to set the tag value.
tag_dataTag data for each element, contiguous in memory. This data is expected to be num_elems*tag_length*sizeof(tag_type) bytes.

Implements MBMesquite::Mesh.

Definition at line 716 of file ArrayMesh.cpp.

References MBMesquite::ArrayMesh::Tag::defaultValue, fill(), MSQ_SETERR, oneBasedArrays, MBMesquite::ArrayMesh::Tag::owned, MBMesquite::ArrayMesh::Tag::size, MBMesquite::MsqError::TAG_ALREADY_EXISTS, vertexCount, MBMesquite::ArrayMesh::Tag::vtxReadPtr, and MBMesquite::ArrayMesh::Tag::vtxWritePtr.

Referenced by ArrayMeshTest::test_tag_data().

{
    Tag* tag = reinterpret_cast< Tag* >( handle );
    if( !tag->vtxWritePtr )
    {
        if( !tag->owned )
        {
            MSQ_SETERR( err )
            ( "Attempt to set non-writeable (application owned) tag data", MsqError::TAG_ALREADY_EXISTS );
            return;
        }
        else
        {
            assert( !tag->vtxReadPtr );
            tag->vtxReadPtr = tag->vtxWritePtr = new unsigned char[vertexCount * tag->size];
            fill( tag->vtxWritePtr, tag->defaultValue, tag->size, vertexCount );
        }
    }
    const unsigned char* ptr = reinterpret_cast< const unsigned char* >( data );
    // as we're working with user-supplied arrays, make sure we're not
    // memcpying overlapping regions
    assert( ptr + tag->size * vertexCount <= tag->vtxWritePtr || tag->vtxWritePtr + tag->size * vertexCount <= ptr );
    for( size_t i = 0; i < count; ++i )
    {
        size_t idx = reinterpret_cast< size_t >( entities[i] ) - oneBasedArrays;
        memcpy( tag->vtxWritePtr + idx * tag->size, ptr + i * tag->size, tag->size );
    }
}
bool MBMesquite::ArrayMesh::valid ( ) const [private]

Definition at line 142 of file ArrayMesh.cpp.

References connArray, elementCount, fixedFlags, nodesPerElement, oneBasedArrays, and vertexCount.

Referenced by ArrayMesh(), and set_mesh().

{
    unsigned long off = oneBasedArrays ? 1 : 0;
    for( unsigned long i = off; i < vertexCount + off; ++i )
    {
        if( fixedFlags[i] != 0 && fixedFlags[i] != 1 )
        {
            std::cerr << "Invalid vertex fixed flag at index " << i << std::endl;
            return false;
        }
    }

    for( unsigned long i = 0; i < elementCount * nodesPerElement; ++i )
    {
        unsigned long j = connArray[i] - oneBasedArrays;
        if( j >= vertexCount )
        {
            std::cerr << "Invalid connectivity index at index " << j << "(element " << j / elementCount << " node "
                      << j % elementCount << ')' << std::endl;
            return false;
        }
    }

    return true;
}
void MBMesquite::ArrayMesh::vertex_get_byte ( const VertexHandle  vertex,
unsigned char *  byte,
MsqError err 
) [virtual]

Retrieve the byte value for the specified vertex or vertices. The byte value is 0 if it has not yet been set via one of the _set_byte() functions.

Implements MBMesquite::Mesh.

Definition at line 381 of file ArrayMesh.cpp.

References oneBasedArrays, vertexByteArray, and vertexCount.

{
    assert( (size_t)vertex < vertexCount + oneBasedArrays );
    *byte = vertexByteArray[(size_t)vertex];
}

Definition at line 280 of file ArrayMesh.cpp.

References oneBasedArrays, and vertexCount.

{
    return new IndexIterator( oneBasedArrays, vertexCount + oneBasedArrays );
}
void MBMesquite::ArrayMesh::vertex_set_byte ( VertexHandle  vertex,
unsigned char  byte,
MsqError err 
) [virtual]

Each vertex has a byte-sized flag that can be used to store flags. This byte's value is neither set nor used by the mesh implementation. It is intended to be used by Mesquite algorithms. Until a vertex's byte has been explicitly set, its value is 0.

Implements MBMesquite::Mesh.

Definition at line 362 of file ArrayMesh.cpp.

References oneBasedArrays, vertexByteArray, and vertexCount.

{
    assert( (size_t)vertex < vertexCount + oneBasedArrays );
    vertexByteArray[(size_t)vertex] = byte;
}
void MBMesquite::ArrayMesh::vertex_set_coordinates ( VertexHandle  vertex,
const Vector3D coordinates,
MsqError err 
) [virtual]

Implements MBMesquite::Mesh.

Definition at line 347 of file ArrayMesh.cpp.

References coordArray, MBMesquite::Vector3D::get_coordinates(), MBMesquite::MsqError::INVALID_STATE, mDimension, MSQ_SETERR, oneBasedArrays, and vertexCount.

{
    size_t i = (size_t)vert;
    assert( i < vertexCount + oneBasedArrays );
    if( mDimension == 3 )
        coordinates.get_coordinates( coordArray + 3 * i );
    else if( mDimension == 2 )
    {
        coordArray[2 * i]     = coordinates[0];
        coordArray[2 * i + 1] = coordinates[1];
    }
    else
        MSQ_SETERR( err )( MsqError::INVALID_STATE );
}
void MBMesquite::ArrayMesh::vertices_get_attached_elements ( const VertexHandle vertex_array,
size_t  num_vertex,
std::vector< ElementHandle > &  elements,
std::vector< size_t > &  offsets,
MsqError err 
) [virtual]

get elements adjacent to vertices

Get adjacency data for vertices

Parameters:
vertex_arrayArray of vertex handles specifying the list of vertices to retrieve adjacency data for.
num_vertexNumber of vertex handles in #vertex_array
elementsThe array in which to place the handles of elements adjacent to the input vertices.
offsetsFor each vertex in #vertex_array, the value in the corresponding position in this array is the index into #elem_array at which the adjacency list begins for that vertex.

Implements MBMesquite::Mesh.

Definition at line 400 of file ArrayMesh.cpp.

References build_vertex_adjacency_list(), oneBasedArrays, vertexAdjacencyList, vertexAdjacencyOffsets, and vertexCount.

{
    const size_t* indices = (const size_t*)vertex_array;
    if( !vertexAdjacencyList ) build_vertex_adjacency_list();

    elements.clear();
    offsets.resize( num_vertex + 1 );
    for( size_t i = 0; i < num_vertex; ++i )
    {
        offsets[i] = elements.size();
        assert( indices[i] < vertexCount + oneBasedArrays );
        for( size_t j = vertexAdjacencyOffsets[indices[i]]; j < vertexAdjacencyOffsets[indices[i] + 1]; ++j )
            elements.push_back( (ElementHandle)vertexAdjacencyList[j] );
    }
    offsets[num_vertex] = elements.size();
}
void MBMesquite::ArrayMesh::vertices_get_byte ( const VertexHandle vertex,
unsigned char *  byte_array,
size_t  array_size,
MsqError err 
) [virtual]

Implements MBMesquite::Mesh.

Definition at line 387 of file ArrayMesh.cpp.

References oneBasedArrays, vertexByteArray, and vertexCount.

{
    const size_t* indices = (const size_t*)vert_array;
    for( size_t i = 0; i < array_size; ++i )
    {
        assert( indices[i] < vertexCount + oneBasedArrays );
        byte_array[i] = vertexByteArray[indices[i]];
    }
}
void MBMesquite::ArrayMesh::vertices_get_coordinates ( const VertexHandle  vert_array[],
MsqVertex coordinates,
size_t  num_vtx,
MsqError err 
) [virtual]

Get/set location of a vertex.

Implements MBMesquite::Mesh.

Definition at line 325 of file ArrayMesh.cpp.

References coordArray, MBMesquite::MsqError::INVALID_STATE, mDimension, MSQ_SETERR, oneBasedArrays, MBMesquite::Vector3D::set(), and vertexCount.

{
    const size_t* indices = (const size_t*)vert_array;
    if( mDimension == 3 )
        for( size_t i = 0; i < num_vtx; ++i )
        {
            assert( indices[i] < vertexCount + oneBasedArrays );
            coordinates[i].set( coordArray + 3 * indices[i] );
        }
    else if( mDimension == 2 )
        for( size_t i = 0; i < num_vtx; ++i )
        {
            assert( indices[i] < vertexCount + oneBasedArrays );
            coordinates[i].set( coordArray[2 * indices[i]], coordArray[2 * indices[i] + 1], 0.0 );
        }
    else
        MSQ_SETERR( err )( MsqError::INVALID_STATE );
}
void MBMesquite::ArrayMesh::vertices_get_fixed_flag ( const VertexHandle  vert_array[],
std::vector< bool > &  fixed_flag_array,
size_t  num_vtx,
MsqError err 
) [virtual]

Returns a pointer to an iterator that iterates over the set of all vertices in this mesh. The calling code should delete the returned iterator when it is finished with it. If vertices are added or removed from the Mesh after obtaining an iterator, the behavior of that iterator is undefined. Returns a pointer to an iterator that iterates over the set of all top-level elements in this mesh. The calling code should delete the returned iterator when it is finished with it. If elements are added or removed from the Mesh after obtaining an iterator, the behavior of that iterator is undefined. Returns true or false, indicating whether the vertex is allowed to be repositioned. True indicates that the vertex is fixed and cannot be moved. Note that this is a read-only property; this flag can't be modified by users of the MBMesquite::Mesh interface.

Implements MBMesquite::Mesh.

Definition at line 290 of file ArrayMesh.cpp.

References fixedFlags, and vertexCount.

{
    fixed_flag_array.resize( num_vtx );
    const size_t* indices = (const size_t*)vert_array;
    for( size_t i = 0; i < num_vtx; ++i )
    {
        assert( indices[i] < vertexCount );
        fixed_flag_array[i] = !!fixedFlags[indices[i]];
    }
}
void MBMesquite::ArrayMesh::vertices_get_slaved_flag ( const VertexHandle  vert_array[],
std::vector< bool > &  slaved_flag_array,
size_t  num_vtx,
MsqError err 
) [virtual]

Returns true or false, indicating whether the vertex is a higher-order node that should be slaved to the logical mid-point of the element side it lies on or not, respectively.

Note: This function will never be called unless this behavior is requested by calling: InstructionQueue::set_slaved_ho_node_mode( Settings::SLAVE_FLAG )

Implements MBMesquite::Mesh.

Definition at line 304 of file ArrayMesh.cpp.

References MBMesquite::MsqError::INVALID_STATE, MSQ_SETERR, slavedFlags, and vertexCount.

{
    if( !slavedFlags )
    {
        MSQ_SETERR( err )
        ( "No data provided to ArrayMesh for Settings::SLAVE_FLAG", MsqError::INVALID_STATE );
        return;
    }

    slaved_flags.resize( num_vtx );
    const size_t* indices = (const size_t*)vert_array;
    for( size_t i = 0; i < num_vtx; ++i )
    {
        assert( indices[i] < vertexCount );
        slaved_flags[i] = !!slavedFlags[indices[i]];
    }
}
void MBMesquite::ArrayMesh::vertices_set_byte ( const VertexHandle vert_array,
const unsigned char *  byte_array,
size_t  array_size,
MsqError err 
) [virtual]

Implements MBMesquite::Mesh.

Definition at line 368 of file ArrayMesh.cpp.

References oneBasedArrays, vertexByteArray, and vertexCount.

Referenced by TerminationCriterionTest::test_abs_vtx_movement_culling().

{
    const size_t* indices = (const size_t*)vert_array;
    for( size_t i = 0; i < array_size; ++i )
    {
        assert( indices[i] < vertexCount + oneBasedArrays );
        vertexByteArray[indices[i]] = byte_array[i];
    }
}

Member Data Documentation

unsigned long* MBMesquite::ArrayMesh::allocConnOffsets [private]

Same as connOffsets if allocated by constructor. NULL if connOffsets is either NULL or application-provided data.

Definition at line 333 of file ArrayMesh.hpp.

Referenced by ArrayMesh(), and clear_mesh().

const unsigned long* MBMesquite::ArrayMesh::connArray [private]

Element connectivity.

Definition at line 328 of file ArrayMesh.hpp.

Referenced by ArrayMesh(), elem_verts(), set_mesh(), and valid().

const unsigned long* MBMesquite::ArrayMesh::connOffsets [private]

Offsets into connectivity array for each element. If NULL, then all elements are of the same type and have the same number of vertices.

Definition at line 329 of file ArrayMesh.hpp.

Referenced by ArrayMesh(), clear_mesh(), and elem_verts().

Interleaved vertex coordinates.

Definition at line 322 of file ArrayMesh.hpp.

Referenced by ArrayMesh(), clear_mesh(), set_mesh(), vertex_set_coordinates(), and vertices_get_coordinates().

Type for all elements if connOffsets is NULL.

Definition at line 337 of file ArrayMesh.hpp.

Referenced by clear_mesh(), elements_get_topologies(), and set_mesh().

Type for each element type if connOffsets is not NULL.

Definition at line 338 of file ArrayMesh.hpp.

Referenced by ArrayMesh(), clear_mesh(), and elements_get_topologies().

const int* MBMesquite::ArrayMesh::fixedFlags [private]

Vertex fixed flags.

Definition at line 323 of file ArrayMesh.hpp.

Referenced by ArrayMesh(), clear_mesh(), set_mesh(), valid(), and vertices_get_fixed_flag().

Nodes per element if connOffsets is NULL.

Definition at line 339 of file ArrayMesh.hpp.

Referenced by ArrayMesh(), clear_mesh(), elem_verts(), set_mesh(), and valid().

const int* MBMesquite::ArrayMesh::slavedFlags [private]

Vertex slaved flags.

Definition at line 324 of file ArrayMesh.hpp.

Referenced by set_mesh(), and vertices_get_slaved_flag().

Definition at line 374 of file ArrayMesh.hpp.

Referenced by allocate_tag(), clear_mesh(), tag_delete(), and tag_get().

List of all members.


The documentation for this class was generated from the following files:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines