MOAB: Mesh Oriented datABase  (version 5.4.1)
moab::ScdBox Class Reference

#include <ScdInterface.hpp>

+ Collaboration diagram for moab::ScdBox:

Public Member Functions

 ~ScdBox ()
 Destructor.
ScdInterfacesc_impl () const
 Return the ScdInterface responsible for this box.
ErrorCode add_vbox (ScdBox *vbox, HomCoord from1, HomCoord to1, HomCoord from2, HomCoord to2, HomCoord from3, HomCoord to3, bool bb_input=false, const HomCoord &bb_min=HomCoord::getUnitv(0), const HomCoord &bb_max=HomCoord::getUnitv(0))
 Add a vertex box to this box.
bool boundary_complete () const
 Return whether this box has all its vertices defined.
int box_dimension () const
 Return highest topological dimension of box.
EntityHandle start_vertex () const
 Starting vertex handle for this box.
EntityHandle start_element () const
 Starting entity handle for this box.
int num_elements () const
 Return the number of elements in the box.
int num_vertices () const
 Return the number of vertices in the box.
const int * box_dims () const
 Return the parametric coordinates for this box.
HomCoord box_min () const
 Return the lower corner parametric coordinates for this box.
HomCoord box_max () const
 Return the upper corner parametric coordinates for this box.
HomCoord box_size () const
 Return the parameter extents for this box.
void box_size (int *ijk) const
 Return the parametric extents for this box.
void box_size (int &i, int &j, int &k) const
 Return the parametric extents for this box.
EntityHandle get_element (const HomCoord &ijk) const
 Get the element at the specified coordinates.
EntityHandle get_element (int i, int j=0, int k=0) const
 Get the element at the specified coordinates.
EntityHandle get_vertex (const HomCoord &ijk) const
 Get the vertex at the specified coordinates.
EntityHandle get_vertex (int i, int j=0, int k=0) const
 Get the vertex at the specified coordinates.
ErrorCode get_params (EntityHandle ent, int &i, int &j, int &k, int &dir) const
 Get parametric coordinates of the specified entity.
ErrorCode get_params (EntityHandle ent, int &i, int &j, int &k) const
ErrorCode get_params (EntityHandle ent, HomCoord &ijkd) const
 Get parametric coordinates of the specified entity.
ErrorCode get_adj_edge_or_face (int dim, int i, int j, int k, int dir, EntityHandle &ent, bool create_if_missing=true) const
 Get the adjacent edge or face at a parametric location This function gets the left (i=0), front (j=0), or bottom (k=0) edge or face for a parametric element. Left, front, or bottom is indicated by dir = 0, 1, or 2, resp. All edges and faces in a structured mesh block can be accessed using these parameters.
bool contains (int i, int j, int k) const
 Return whether the box contains the parameters passed in.
bool contains (const HomCoord &ijk) const
 Return whether the box contains the parameters passed in.
void box_set (EntityHandle this_set)
 Set/Get the entity set representing the box.
EntityHandle box_set ()
ErrorCode get_coordinate_arrays (double *&xc, double *&yc, double *&zc)
 Get coordinate arrays for vertex coordinates for a structured block.
ErrorCode get_coordinate_arrays (const double *&xc, const double *&yc, const double *&zc) const
 Get read-only coordinate arrays for vertex coordinates for a structured block.
bool locally_periodic_i () const
 Return whether box is locally periodic in i.
bool locally_periodic_j () const
 Return whether box is locally periodic in j.
bool locally_periodic_k () const
 Return whether box is locally periodic in k.
void locally_periodic (bool lperiodic[3])
 Set local periodicity.
const int * locally_periodic () const
 Get local periodicity.
ScdParDatapar_data ()
 Return parallel data.
const ScdParDatapar_data () const
 Return parallel data.
void par_data (const ScdParData &par_datap)
 set parallel data

Private Member Functions

 ScdBox (ScdInterface *sc_impl, EntityHandle box_set, EntitySequence *seq1, EntitySequence *seq2=NULL)
 Constructor.
EntityHandle get_vertex_from_seq (int i, int j, int k) const
 function to get vertex handle directly from sequence
ErrorCode vert_dat (ScdVertexData *vert_dat)
 set the vertex sequence
ScdVertexDatavert_dat () const
 get the vertex sequence
ErrorCode elem_seq (EntitySequence *elem_seq)
 set the element sequence
StructuredElementSeqelem_seq () const
 get the element sequence
void start_vertex (EntityHandle startv)
 Set the starting vertex handle for this box.
void start_element (EntityHandle starte)
 Set the starting entity handle for this box.

Private Attributes

ScdInterfacescImpl
 interface instance
EntityHandle boxSet
 entity set representing this box
ScdVertexDatavertDat
StructuredElementSeqelemSeq
 element sequence this box represents
EntityHandle startVertex
 starting vertex handle for this box
EntityHandle startElem
 starting element handle for this box
int boxDims [6]
 lower and upper corners
int locallyPeriodic [3]
 is locally periodic in i or j or k
ScdParData parData
 parallel data associated with this box, if any
HomCoord boxSize
 parameter extents
int boxSizeIJ
 convenience parameters, (boxSize[1]-1)*(boxSize[0]-1) and boxSize[0]-1
int boxSizeIJM1
int boxSizeIM1

Friends

class ScdInterface

Detailed Description

Examples:
StructuredMeshSimple.cpp.

Definition at line 506 of file ScdInterface.hpp.


Constructor & Destructor Documentation

Destructor.

Definition at line 548 of file ScdInterface.cpp.

References moab::ScdInterface::box_set_tag(), boxSet, moab::Core::is_valid(), mbcore, moab::ScdInterface::mbImpl, moab::ScdInterface::remove_box(), scImpl, and moab::Interface::tag_set_data().

{
    // Reset the tag on the set
    if( boxSet )
    {
        // It is possible that the box set entity has been deleted (e.g. by
        // Core::clean_up_failed_read)
        Core* mbcore = dynamic_cast< Core* >( scImpl->mbImpl );
        assert( mbcore != NULL );
        if( mbcore->is_valid( boxSet ) )
        {
            ScdBox* tmp_ptr = NULL;
            scImpl->mbImpl->tag_set_data( scImpl->box_set_tag(), &boxSet, 1, &tmp_ptr );
        }
        else
            boxSet = 0;
    }

    scImpl->remove_box( this );
}
moab::ScdBox::ScdBox ( ScdInterface sc_impl,
EntityHandle  box_set,
EntitySequence seq1,
EntitySequence seq2 = NULL 
) [private]

Constructor.

Create a structured box instance; this constructor is private because it should only be called from ScdInterface, a friend class. This constructor takes two sequences, one of which can be NULL. If both sequences come in non-NULL, the first should be a VertexSequence* corresponding to a structured vertex sequence and the second should be a StructuredElementSeq*. If the 2nd is NULL, the first can be either of those types. The other members of this class are taken from the sequences (e.g. parametric space) or the box set argument. Tags on the box set should be set from the caller.

Parameters:
sc_implA ScdInterface instance
box_setEntity set representing this rectangle
seq1An EntitySequence (see ScdBox description)
seq2An EntitySequence (see ScdBox description), or NULL

Definition at line 474 of file ScdInterface.cpp.

References moab::ScdInterface::add_box(), moab::Range::begin(), boxDims, moab::ScdInterface::boxDimsTag, moab::ScdInterface::boxPeriodicTag, boxSize, boxSizeIJ, boxSizeIJM1, boxSizeIM1, moab::EntitySequence::data(), elemSeq, moab::Range::empty(), ErrorCode, moab::Interface::get_entities_by_dimension(), locallyPeriodic, moab::ScdVertexData::max_params(), moab::ScdElementData::max_params(), MB_SUCCESS, moab::ScdInterface::mbImpl, moab::ScdVertexData::min_params(), moab::ScdElementData::min_params(), scImpl, moab::StructuredElementSeq::sdata(), moab::SequenceData::start_handle(), moab::EntitySequence::start_handle(), startElem, startVertex, moab::Interface::tag_get_data(), and vertDat.

    : scImpl( impl ), boxSet( bset ), vertDat( NULL ), elemSeq( NULL ), startVertex( 0 ), startElem( 0 )
{
    for( int i = 0; i < 6; i++ )
        boxDims[i] = 0;
    for( int i = 0; i < 3; i++ )
        locallyPeriodic[i] = false;
    VertexSequence* vseq = dynamic_cast< VertexSequence* >( seq1 );
    if( vseq ) vertDat = dynamic_cast< ScdVertexData* >( vseq->data() );
    if( vertDat )
    {
        // retrieve the parametric space
        for( int i = 0; i < 3; i++ )
        {
            boxDims[i]     = vertDat->min_params()[i];
            boxDims[3 + i] = vertDat->max_params()[i];
        }
        startVertex = vertDat->start_handle();
    }
    else if( impl->boxDimsTag )
    {
        // look for parametric space info on set
        ErrorCode rval = impl->mbImpl->tag_get_data( impl->boxDimsTag, &bset, 1, boxDims );
        if( MB_SUCCESS == rval )
        {
            Range verts;
            impl->mbImpl->get_entities_by_dimension( bset, 0, verts );
            if( !verts.empty() ) startVertex = *verts.begin();
        }
    }

    elemSeq = dynamic_cast< StructuredElementSeq* >( seq2 );
    if( !elemSeq ) elemSeq = dynamic_cast< StructuredElementSeq* >( seq1 );

    if( elemSeq )
    {
        if( vertDat )
        {
            // check the parametric space to make sure it's consistent
            assert( elemSeq->sdata()->min_params() == HomCoord( boxDims, 3 ) &&
                    ( elemSeq->sdata()->max_params() + HomCoord( 1, 1, 1 ) ) == HomCoord( boxDims, 3 ) );
        }
        else
        {
            // get the parametric space from the element sequence
            for( int i = 0; i < 3; i++ )
            {
                boxDims[i]     = elemSeq->sdata()->min_params()[i];
                boxDims[3 + i] = elemSeq->sdata()->max_params()[i];
            }
        }

        startElem = elemSeq->start_handle();
    }
    else
    {
        Range elems;
        impl->mbImpl->get_entities_by_dimension(
            bset, ( boxDims[2] == boxDims[5] ? ( boxDims[1] == boxDims[4] ? 1 : 2 ) : 3 ), elems );
        if( !elems.empty() ) startElem = *elems.begin();
        // call the following w/o looking at return value, since it doesn't really need to be there
        if( impl->boxPeriodicTag ) impl->mbImpl->tag_get_data( impl->boxPeriodicTag, &bset, 1, locallyPeriodic );
    }

    assert( vertDat || elemSeq || boxDims[0] != boxDims[3] || boxDims[1] != boxDims[4] || boxDims[2] != boxDims[5] );

    boxSize     = HomCoord( boxDims + 3, 3 ) - HomCoord( boxDims, 3 ) + HomCoord( 1, 1, 1 );
    boxSizeIJ   = ( boxSize[1] ? boxSize[1] : 1 ) * boxSize[0];
    boxSizeIM1  = boxSize[0] - ( locallyPeriodic[0] ? 0 : 1 );
    boxSizeIJM1 = ( boxSize[1] ? ( boxSize[1] - ( locallyPeriodic[1] ? 0 : 1 ) ) : 1 ) * boxSizeIM1;

    scImpl->add_box( this );
}

Member Function Documentation

ErrorCode moab::ScdBox::add_vbox ( ScdBox vbox,
HomCoord  from1,
HomCoord  to1,
HomCoord  from2,
HomCoord  to2,
HomCoord  from3,
HomCoord  to3,
bool  bb_input = false,
const HomCoord bb_min = HomCoord::getUnitv( 0 ),
const HomCoord bb_max = HomCoord::getUnitv( 0 ) 
)

Add a vertex box to this box.

Definition at line 580 of file ScdInterface.cpp.

References moab::ScdElementData::add_vsequence(), elemSeq, ErrorCode, moab::StructuredElementSeq::sdata(), and vertDat.

Referenced by moab::ScdInterface::construct_box(), create_1d_3_sequences(), create_2d_3_sequences(), create_2dtri_3_sequences(), create_3dtri_3_sequences(), eseq_test1a(), eseq_test1b(), and eseq_test1c().

{
    if( !vbox->vertDat ) return MB_FAILURE;
    ScdVertexData* dum_data = dynamic_cast< ScdVertexData* >( vbox->vertDat );
    ErrorCode rval =
        elemSeq->sdata()->add_vsequence( dum_data, from1, to1, from2, to2, from3, to3, bb_input, bb_min, bb_max );
    return rval;
}

Return whether this box has all its vertices defined.

Tests whether vertex boxs added with add_vbox have completely defined the vertex parametric space for this box.

Definition at line 598 of file ScdInterface.cpp.

References moab::StructuredElementSeq::boundary_complete(), and elemSeq.

Referenced by check_element_sequence(), and moab::ScdNCHelper::create_mesh().

{
    return elemSeq->boundary_complete();
}
int moab::ScdBox::box_dimension ( ) const [inline]

Return highest topological dimension of box.

Definition at line 575 of file ScdInterface.cpp.

References moab::Interface::dimension_from_handle(), moab::ScdInterface::mbImpl, scImpl, and startElem.

Referenced by get_params(), and moab::ScdInterface::tag_shared_vertices().

const int * moab::ScdBox::box_dims ( ) const [inline]

Return the parametric coordinates for this box.

Returns:
IJK parameters of lower and upper corners

Definition at line 1520 of file ScdInterface.hpp.

References boxDims.

Referenced by moab::ScdInterface::assign_global_ids(), and moab::ScdInterface::get_shared_vertices().

{
    return boxDims;
}
HomCoord moab::ScdBox::box_max ( ) const [inline]

Return the upper corner parametric coordinates for this box.

Definition at line 1530 of file ScdInterface.hpp.

References boxDims.

Referenced by access_adjacencies(), check_element_sequence(), check_vertex_sequence(), evaluate_element_sequence(), evaluate_vertex_sequence(), and moab::Skinner::skin_box().

{
    return HomCoord( boxDims + 3, 3 );
}
HomCoord moab::ScdBox::box_min ( ) const [inline]

Return the lower corner parametric coordinates for this box.

Definition at line 1525 of file ScdInterface.hpp.

References boxDims.

Referenced by access_adjacencies(), check_element_sequence(), check_vertex_sequence(), evaluate_element_sequence(), evaluate_vertex_sequence(), and moab::Skinner::skin_box().

{
    return HomCoord( boxDims, 3 );
}
void moab::ScdBox::box_set ( EntityHandle  this_set) [inline]

Set/Get the entity set representing the box.

Definition at line 1589 of file ScdInterface.hpp.

References boxSet.

Referenced by moab::ScdInterface::construct_box(), moab::ScdNCHelper::create_mesh(), create_parallel_mesh(), iMesh_createStructuredMesh(), and moab::ScdInterface::tag_shared_vertices().

{
    boxSet = this_set;
}

Definition at line 1594 of file ScdInterface.hpp.

References boxSet.

{
    return boxSet;
}
HomCoord moab::ScdBox::box_size ( ) const [inline]

Return the parameter extents for this box.

Definition at line 1535 of file ScdInterface.hpp.

References boxSize.

Referenced by access_adjacencies(), check_element_sequence(), check_vertex_sequence(), and moab::ScdInterface::construct_box().

{
    return boxSize;
}
void moab::ScdBox::box_size ( int *  ijk) const [inline]

Return the parametric extents for this box.

Parameters:
ijkIJK extents of this box

Definition at line 1540 of file ScdInterface.hpp.

References boxSize.

{
    ijk[0] = boxSize[0];
    ijk[1] = boxSize[1];
    ijk[2] = boxSize[2];
}
void moab::ScdBox::box_size ( int &  i,
int &  j,
int &  k 
) const [inline]

Return the parametric extents for this box.

Parameters:
iI extent of this box
jJ extent of this box
kK extent of this box

Definition at line 1547 of file ScdInterface.hpp.

References boxSize.

{
    i = boxSize[0];
    j = boxSize[1];
    k = boxSize[2];
}
bool moab::ScdBox::contains ( int  i,
int  j,
int  k 
) const [inline]

Return whether the box contains the parameters passed in.

Parameters:
iParametric coordinates being evaluated
jParametric coordinates being evaluated
kParametric coordinates being evaluated

Definition at line 1584 of file ScdInterface.hpp.

Referenced by evaluate_element_sequence(), and evaluate_vertex_sequence().

{
    return contains( HomCoord( i, j, k ) );
}
bool moab::ScdBox::contains ( const HomCoord ijk) const [inline]

Return whether the box contains the parameters passed in.

Parameters:
iParametric coordinates being evaluated
jParametric coordinates being evaluated
kParametric coordinates being evaluated

Definition at line 1579 of file ScdInterface.hpp.

References boxDims.

{
    return ( ijk >= HomCoord( boxDims, 3 ) && ijk <= HomCoord( boxDims + 3, 3 ) );
}
ErrorCode moab::ScdBox::elem_seq ( EntitySequence elem_seq) [private]

set the element sequence

Definition at line 628 of file ScdInterface.cpp.

References boxSize, boxSizeIJM1, boxSizeIM1, elemSeq, moab::StructuredElementSeq::is_periodic(), locallyPeriodic, and MB_SUCCESS.

Referenced by moab::ScdInterface::construct_box().

{
    elemSeq = dynamic_cast< StructuredElementSeq* >( elem_sq );
    if( elemSeq ) elemSeq->is_periodic( locallyPeriodic );

    if( locallyPeriodic[0] ) boxSizeIM1 = boxSize[0] - ( locallyPeriodic[0] ? 0 : 1 );
    if( locallyPeriodic[0] || locallyPeriodic[1] )
        boxSizeIJM1 = ( boxSize[1] ? ( boxSize[1] - ( locallyPeriodic[1] ? 0 : 1 ) ) : 1 ) * boxSizeIM1;

    return ( elemSeq ? MB_SUCCESS : MB_FAILURE );
}
StructuredElementSeq * moab::ScdBox::elem_seq ( ) const [inline, private]

get the element sequence

Definition at line 1604 of file ScdInterface.hpp.

References elemSeq.

{
    return elemSeq;
}
ErrorCode moab::ScdBox::get_adj_edge_or_face ( int  dim,
int  i,
int  j,
int  k,
int  dir,
EntityHandle ent,
bool  create_if_missing = true 
) const

Get the adjacent edge or face at a parametric location This function gets the left (i=0), front (j=0), or bottom (k=0) edge or face for a parametric element. Left, front, or bottom is indicated by dir = 0, 1, or 2, resp. All edges and faces in a structured mesh block can be accessed using these parameters.

Get the entity of specified dimension adjacent to parametric element.

Parameters:
dimDimension of adjacent entity being requested
iParametric coordinates of cell being evaluated
jParametric coordinates of cell being evaluated
kParametric coordinates of cell being evaluated
dirDirection (0, 1, or 2), for getting adjacent edges (2d, 3d) or faces (3d)
entEntity returned from this function
create_if_missingIf true, creates the entity if it doesn't already exist
dimDimension of adjacent entity being requested
iParametric coordinates of cell being evaluated
jParametric coordinates of cell being evaluated
kParametric coordinates of cell being evaluated
dirDirection (0, 1, or 2), for getting adjacent edges (2d, 3d) or faces (3d)
entEntityHandle of adjacent entity
create_if_missingIf true, creates the entity if it doesn't already exist

Definition at line 669 of file ScdInterface.cpp.

References moab::Range::begin(), boxDims, moab::Interface::create_element(), dim, ErrorCode, moab::Interface::get_adjacencies(), get_vertex(), moab::ScdInterface::impl(), locallyPeriodic, MB_SUCCESS, MBEDGE, MBQUAD, scImpl, and moab::Range::size().

Referenced by access_adjacencies(), and moab::Skinner::skin_box().

{
    // describe connectivity of sub-element in static array
    // subconnect[dim-1][dir][numv][ijk] where dimensions are:
    // [dim-1]: dim=1 or 2, so this is 0 or 1
    // [dir]: one of 0..2, for ijk directions in a hex
    // [numv]: number of vertices describing sub entity = 2*dim <= 4
    // [ijk]: 3 values for i, j, k
    int subconnect[2][3][4][3] = { { { { 0, 0, 0 }, { 1, 0, 0 }, { -1, -1, -1 }, { -1, -1, -1 } },    // i edge
                                     { { 0, 0, 0 }, { 0, 1, 0 }, { -1, -1, -1 }, { -1, -1, -1 } },    // j edge
                                     { { 0, 0, 0 }, { 0, 0, 1 }, { -1, -1, -1 }, { -1, -1, -1 } } },  // k edge

                                   { { { 0, 0, 0 }, { 0, 1, 0 }, { 0, 1, 1 }, { 0, 0, 1 } },      // i face
                                     { { 0, 0, 0 }, { 1, 0, 0 }, { 1, 0, 1 }, { 0, 0, 1 } },      // j face
                                     { { 0, 0, 0 }, { 1, 0, 0 }, { 1, 1, 0 }, { 0, 1, 0 } } } };  // k face

    // check proper input dimensions and lower bound
    if( dim < 1 || dim > 2 || i < boxDims[0] || j < boxDims[1] || k < boxDims[2] ) return MB_FAILURE;

    // now check upper bound; parameters must be <= upper corner, since edges/faces
    // follow element parameterization, not vertex parameterization
    else if( ( boxDims[3] != boxDims[0] && i > ( locallyPeriodic[0] ? boxDims[3] + 1 : boxDims[3] ) ) ||
             ( boxDims[4] != boxDims[1] && j > ( locallyPeriodic[1] ? boxDims[4] + 1 : boxDims[4] ) ) ||
             ( boxDims[5] != boxDims[2] && k > boxDims[5] ) )
        return MB_FAILURE;

    // get the vertices making up this entity
    EntityHandle verts[4];
    for( int ind = 0; ind < 2 * dim; ind++ )
    {
        int i1 = i + subconnect[dim - 1][dir][ind][0];
        int j1 = j + subconnect[dim - 1][dir][ind][1];
        // if periodic in i and i1 is boxDims[3]+1, wrap around
        if( locallyPeriodic[0] && i1 == boxDims[3] + 1 ) i1 = boxDims[0];
        // if periodic in j and j1 is boxDims[4]+1, wrap around
        if( locallyPeriodic[1] && j1 == boxDims[4] + 1 ) j1 = boxDims[1];
        verts[ind] = get_vertex( i1, j1, k + subconnect[dim - 1][dir][ind][2] );
        if( !verts[ind] ) return MB_FAILURE;
    }

    Range ents;
    ErrorCode rval = scImpl->impl()->get_adjacencies( verts, 2 * dim, dim, false, ents );
    if( MB_SUCCESS != rval ) return rval;

    if( ents.size() > 1 )
        return MB_FAILURE;

    else if( ents.size() == 1 )
    {
        ent = *ents.begin();
    }
    else if( create_if_missing )
        rval = scImpl->impl()->create_element( ( 1 == dim ? MBEDGE : MBQUAD ), verts, 2 * dim, ent );

    return rval;
}
ErrorCode moab::ScdBox::get_coordinate_arrays ( double *&  xc,
double *&  yc,
double *&  zc 
)

Get coordinate arrays for vertex coordinates for a structured block.

Returns error if there isn't a single vertex sequence associated with this structured block

Parameters:
xcX coordinate array pointer returned
ycY coordinate array pointer returned
zcZ coordinate array pointer returned

Definition at line 603 of file ScdInterface.cpp.

References moab::SequenceData::get_sequence_data(), MB_SUCCESS, and vertDat.

Referenced by moab::ScdInterface::construct_box(), moab::ScdNCHelper::create_mesh(), and iMesh_createStructuredMesh().

{
    if( !vertDat ) return MB_FAILURE;

    xc = reinterpret_cast< double* >( vertDat->get_sequence_data( 0 ) );
    yc = reinterpret_cast< double* >( vertDat->get_sequence_data( 1 ) );
    zc = reinterpret_cast< double* >( vertDat->get_sequence_data( 2 ) );
    return MB_SUCCESS;
}
ErrorCode moab::ScdBox::get_coordinate_arrays ( const double *&  xc,
const double *&  yc,
const double *&  zc 
) const

Get read-only coordinate arrays for vertex coordinates for a structured block.

Returns error if there isn't a single vertex sequence associated with this structured block

Parameters:
xcX coordinate array pointer returned
ycY coordinate array pointer returned
zcZ coordinate array pointer returned

Definition at line 613 of file ScdInterface.cpp.

References moab::SequenceData::get_sequence_data(), MB_SUCCESS, and vertDat.

{
    if( !vertDat ) return MB_FAILURE;
    xc = reinterpret_cast< const double* >( vertDat->get_sequence_data( 0 ) );
    yc = reinterpret_cast< const double* >( vertDat->get_sequence_data( 1 ) );
    zc = reinterpret_cast< const double* >( vertDat->get_sequence_data( 2 ) );
    return MB_SUCCESS;
}
EntityHandle moab::ScdBox::get_element ( const HomCoord ijk) const [inline]

Get the element at the specified coordinates.

Parameters:
ijkParametric coordinates being evaluated
Examples:
StructuredMeshSimple.cpp.

Definition at line 1561 of file ScdInterface.hpp.

Referenced by moab::Core::create_scd_sequence(), evaluate_element_sequence(), and main().

{
    return get_element( ijk[0], ijk[1], ijk[2] );
}
EntityHandle moab::ScdBox::get_element ( int  i,
int  j = 0,
int  k = 0 
) const [inline]

Get the element at the specified coordinates.

Parameters:
iParametric coordinates being evaluated
jParametric coordinates being evaluated
kParametric coordinates being evaluated

Definition at line 1554 of file ScdInterface.hpp.

References boxDims, boxSizeIJM1, boxSizeIM1, and startElem.

{
    return ( !startElem
                 ? 0
                 : startElem + ( k - boxDims[2] ) * boxSizeIJM1 + ( j - boxDims[1] ) * boxSizeIM1 + i - boxDims[0] );
}
ErrorCode moab::ScdBox::get_params ( EntityHandle  ent,
int &  i,
int &  j,
int &  k,
int &  dir 
) const [inline]

Get parametric coordinates of the specified entity.

This function returns MB_ENTITY_NOT_FOUND if the entity is not in this ScdBox.

Parameters:
entEntity being queried
iParametric coordinates returned
jParametric coordinates returned
kParametric coordinates returned
dirParametric coordinate direction returned (in case of getting adjacent edges (2d, 3d) or faces (3d); not modified otherwise

Definition at line 1609 of file ScdInterface.hpp.

References ErrorCode, and MB_SUCCESS.

Referenced by evaluate_element_sequence(), evaluate_vertex_sequence(), and get_params().

{
    HomCoord hc;
    ErrorCode rval = get_params( ent, hc );
    if( MB_SUCCESS == rval )
    {
        i   = hc[0];
        j   = hc[1];
        k   = hc[2];
        dir = hc[3];
    }

    return rval;
}
ErrorCode moab::ScdBox::get_params ( EntityHandle  ent,
int &  i,
int &  j,
int &  k 
) const [inline]

Get parametric coordinates of the specified entity, intermediate entities not allowed (no dir parameter) This function returns MB_ENTITY_NOT_FOUND if the entity is not in this ScdBox, or MB_FAILURE if the entity is an intermediate-dimension entity.

Parameters:
entEntity being queried
iParametric coordinates returned
jParametric coordinates returned
kParametric coordinates returned

Definition at line 1624 of file ScdInterface.hpp.

References ErrorCode, get_params(), and MB_SUCCESS.

{
    HomCoord hc;
    ErrorCode rval = get_params( ent, hc );
    if( MB_SUCCESS == rval )
    {
        i = hc[0];
        j = hc[1];
        k = hc[2];
    }

    return rval;
}

Get parametric coordinates of the specified entity.

This function returns MB_ENTITY_NOT_FOUND if the entity is not in this ScdBox.

Parameters:
entEntity being queried
ijkdParametric coordinates returned (including direction, in case of getting adjacent edges (2d, 3d) or faces (3d))

Definition at line 640 of file ScdInterface.cpp.

References box_dimension(), moab::Interface::dimension_from_handle(), elemSeq, moab::ScdVertexData::get_params(), moab::StructuredElementSeq::get_params(), moab::ScdInterface::impl(), MB_NOT_IMPLEMENTED, scImpl, and vertDat.

{
    // check first whether this is an intermediate entity, so we know what to do
    int dimension = box_dimension();
    int this_dim  = scImpl->impl()->dimension_from_handle( ent );

    if( ( 0 == this_dim && !vertDat ) || ( this_dim && this_dim == dimension ) )
    {
        assert( elemSeq );
        return elemSeq->get_params( ent, ijkd[0], ijkd[1], ijkd[2] );
    }

    else if( !this_dim && vertDat )
        return vertDat->get_params( ent, ijkd[0], ijkd[1], ijkd[2] );

    else
        return MB_NOT_IMPLEMENTED;
}
EntityHandle moab::ScdBox::get_vertex ( const HomCoord ijk) const [inline]

Get the vertex at the specified coordinates.

Parameters:
ijkParametric coordinates being evaluated

Definition at line 1574 of file ScdInterface.hpp.

Referenced by moab::Core::create_scd_sequence(), evaluate_element_sequence(), evaluate_vertex_sequence(), get_adj_edge_or_face(), and main().

{
    return get_vertex( ijk[0], ijk[1], ijk[2] );
}
EntityHandle moab::ScdBox::get_vertex ( int  i,
int  j = 0,
int  k = 0 
) const [inline]

Get the vertex at the specified coordinates.

Parameters:
iParametric coordinates being evaluated
jParametric coordinates being evaluated
kParametric coordinates being evaluated

Definition at line 1566 of file ScdInterface.hpp.

References boxDims, boxSize, boxSizeIJ, get_vertex_from_seq(), startVertex, and vertDat.

{
    return ( vertDat
                 ? startVertex + ( boxDims[2] == -1 && boxDims[5] == -1 ? 0 : ( k - boxDims[2] ) ) * boxSizeIJ +
                       ( boxDims[1] == -1 && boxDims[4] == -1 ? 0 : ( j - boxDims[1] ) ) * boxSize[0] + i - boxDims[0]
                 : get_vertex_from_seq( i, j, k ) );
}
EntityHandle moab::ScdBox::get_vertex_from_seq ( int  i,
int  j,
int  k 
) const [private]

function to get vertex handle directly from sequence

Parameters:
iParameter being queried
jParameter being queried
kParameter being queried

Definition at line 569 of file ScdInterface.cpp.

References elemSeq, and moab::StructuredElementSeq::get_vertex().

Referenced by get_vertex().

{
    assert( elemSeq );
    return elemSeq->get_vertex( i, j, k );
}
void moab::ScdBox::locally_periodic ( bool  lperiodic[3]) [inline]

Set local periodicity.

Parameters:
lperiodicVector of ijk periodicities to set this box to

Definition at line 1653 of file ScdInterface.hpp.

References locallyPeriodic.

Referenced by moab::ScdInterface::assign_global_ids().

{
    for( int i = 0; i < 3; i++ )
        locallyPeriodic[i] = lperiodic[i];
}
const int * moab::ScdBox::locally_periodic ( ) const [inline]

Get local periodicity.

Returns:
Vector of ijk periodicities for this box

Definition at line 1659 of file ScdInterface.hpp.

References locallyPeriodic.

{
    return locallyPeriodic;
}
bool moab::ScdBox::locally_periodic_i ( ) const [inline]

Return whether box is locally periodic in i.

Return whether box is locally periodic in i

Returns:
True if box is locally periodic in i direction

Definition at line 1638 of file ScdInterface.hpp.

References locallyPeriodic.

Referenced by evaluate_element_sequence().

{
    return locallyPeriodic[0];
}
bool moab::ScdBox::locally_periodic_j ( ) const [inline]

Return whether box is locally periodic in j.

Return whether box is locally periodic in j

Returns:
True if box is locally periodic in j direction

Definition at line 1643 of file ScdInterface.hpp.

References locallyPeriodic.

Referenced by evaluate_element_sequence().

{
    return locallyPeriodic[1];
}
bool moab::ScdBox::locally_periodic_k ( ) const [inline]

Return whether box is locally periodic in k.

Return whether box is locally periodic in k

Returns:
True if box is locally periodic in k direction

Definition at line 1648 of file ScdInterface.hpp.

References locallyPeriodic.

{
    return locallyPeriodic[2];
}
int moab::ScdBox::num_elements ( ) const [inline]

Return the number of elements in the box.

Definition at line 1493 of file ScdInterface.hpp.

References boxSize, locallyPeriodic, and startElem.

Referenced by moab::ScdInterface::construct_box(), moab::ScdNCHelper::create_mesh(), moab::ScdInterface::create_scd_sequence(), eseq_test1a(), iMesh_createStructuredMesh(), mb_skin_scd_test(), moab::ScdInterface::tag_shared_vertices(), and test_bound_box().

{
    if( !startElem ) return 0;  // not initialized yet

    /* for a structured mesh, total number of elements is obtained by multiplying
        number of elements in each direction
      number of elements in each direction is given by number of vertices in that direction minus 1
      if periodic in that direction, the last vertex is the same as first one, count one more
      element
      */
    int num_e_i = ( -1 == boxSize[0] || 1 == boxSize[0] ) ? 1 : boxSize[0] - 1;
    if( locallyPeriodic[0] ) ++num_e_i;

    int num_e_j = ( -1 == boxSize[1] || 1 == boxSize[1] ) ? 1 : boxSize[1] - 1;
    if( locallyPeriodic[1] ) ++num_e_j;

    int num_e_k = ( -1 == boxSize[2] || 1 == boxSize[2] ) ? 1 : boxSize[2] - 1;
    if( locallyPeriodic[2] ) ++num_e_k;

    return num_e_i * num_e_j * num_e_k;
}
int moab::ScdBox::num_vertices ( ) const [inline]

Return the number of vertices in the box.

Definition at line 1515 of file ScdInterface.hpp.

References boxSize.

Referenced by moab::ScdInterface::assign_global_ids(), moab::ScdNCHelper::create_mesh(), moab::ScdInterface::create_scd_sequence(), and iMesh_createStructuredMesh().

{
    return boxSize[0] * ( !boxSize[1] ? 1 : boxSize[1] ) * ( !boxSize[2] ? 1 : boxSize[2] );
}

Return parallel data.

Return parallel data, if there is any

Returns:
par_data Parallel data set on this box

Definition at line 749 of file ScdInterface.hpp.

Referenced by moab::ScdInterface::assign_global_ids(), moab::ScdInterface::construct_box(), moab::ScdInterface::get_shared_vertices(), iMesh_createStructuredMesh(), and moab::ScdInterface::tag_shared_vertices().

    {
        return parData;
    }
const ScdParData& moab::ScdBox::par_data ( ) const [inline]

Return parallel data.

Return parallel data, if there is any

Returns:
par_data Parallel data set on this box

Definition at line 758 of file ScdInterface.hpp.

    {
        return parData;
    }
void moab::ScdBox::par_data ( const ScdParData par_datap) [inline]

set parallel data

Set parallel data for this box

Parameters:
par_dataParallel data to be set on this box

Definition at line 767 of file ScdInterface.hpp.

    {
        parData = par_datap;
    }
ScdInterface * moab::ScdBox::sc_impl ( ) const [inline]

Return the ScdInterface responsible for this box.

Definition at line 1468 of file ScdInterface.hpp.

References scImpl.

Referenced by access_adjacencies().

{
    return scImpl;
}
void moab::ScdBox::start_element ( EntityHandle  starte) [inline, private]

Set the starting entity handle for this box.

Definition at line 1488 of file ScdInterface.hpp.

References startElem.

{
    startElem = starte;
}
void moab::ScdBox::start_vertex ( EntityHandle  startv) [inline, private]

Set the starting vertex handle for this box.

Definition at line 1478 of file ScdInterface.hpp.

References startVertex.

{
    startVertex = startv;
}
ErrorCode moab::ScdBox::vert_dat ( ScdVertexData vert_dat) [private]

set the vertex sequence

Definition at line 622 of file ScdInterface.cpp.

References MB_SUCCESS, and vertDat.

{
    vertDat = vert_dt;
    return MB_SUCCESS;
}
ScdVertexData * moab::ScdBox::vert_dat ( ) const [inline, private]

get the vertex sequence

Definition at line 1599 of file ScdInterface.hpp.

References vertDat.

{
    return vertDat;
}

Friends And Related Function Documentation

friend class ScdInterface [friend]

Definition at line 508 of file ScdInterface.hpp.


Member Data Documentation

int moab::ScdBox::boxDims[6] [private]

lower and upper corners

Definition at line 832 of file ScdInterface.hpp.

Referenced by box_dims(), box_max(), box_min(), contains(), get_adj_edge_or_face(), get_element(), get_vertex(), and ScdBox().

entity set representing this box

Definition at line 816 of file ScdInterface.hpp.

Referenced by box_set(), and ~ScdBox().

parameter extents

Definition at line 841 of file ScdInterface.hpp.

Referenced by box_size(), elem_seq(), get_vertex(), num_elements(), num_vertices(), and ScdBox().

int moab::ScdBox::boxSizeIJ [private]

convenience parameters, (boxSize[1]-1)*(boxSize[0]-1) and boxSize[0]-1

Definition at line 844 of file ScdInterface.hpp.

Referenced by get_vertex(), and ScdBox().

Definition at line 845 of file ScdInterface.hpp.

Referenced by elem_seq(), get_element(), and ScdBox().

int moab::ScdBox::boxSizeIM1 [private]

Definition at line 846 of file ScdInterface.hpp.

Referenced by elem_seq(), get_element(), and ScdBox().

element sequence this box represents

Definition at line 823 of file ScdInterface.hpp.

Referenced by add_vbox(), boundary_complete(), elem_seq(), get_params(), get_vertex_from_seq(), and ScdBox().

parallel data associated with this box, if any

Definition at line 838 of file ScdInterface.hpp.

interface instance

Definition at line 813 of file ScdInterface.hpp.

Referenced by box_dimension(), get_adj_edge_or_face(), get_params(), sc_impl(), ScdBox(), and ~ScdBox().

starting element handle for this box

Definition at line 829 of file ScdInterface.hpp.

Referenced by box_dimension(), get_element(), num_elements(), ScdBox(), and start_element().

starting vertex handle for this box

Definition at line 826 of file ScdInterface.hpp.

Referenced by get_vertex(), ScdBox(), and start_vertex().

vertex sequence this box represents, if there's only one, otherwise they're retrieved from the element sequence

Definition at line 820 of file ScdInterface.hpp.

Referenced by add_vbox(), get_coordinate_arrays(), get_params(), get_vertex(), ScdBox(), and vert_dat().

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