![]() |
Mesh Oriented datABase
(version 5.4.1)
Array-based unstructured mesh datastructure
|
#include <ScdElementData.hpp>
Classes | |
class | VertexDataRef |
structure to hold references to bounding vertex blocks More... | |
Public Member Functions | |
ScdElementData (EntityHandle start_handle, const int imin, const int jmin, const int kmin, const int imax, const int jmax, const int kmax, int *is_periodic) | |
constructor | |
virtual | ~ScdElementData () |
EntityHandle | get_vertex (const HomCoord &coords) const |
get handle of vertex at homogeneous coords | |
EntityHandle | get_vertex (int i, int j, int k) const |
EntityHandle | get_element (const int i, const int j, const int k) const |
get handle of element at i, j, k | |
const HomCoord & | min_params () const |
get min params for this element | |
const HomCoord & | max_params () const |
get max params for this element | |
void | param_extents (int &di, int &dj, int &dk) const |
get the number of vertices in each direction, inclusive | |
ErrorCode | get_params (const EntityHandle ehandle, int &i, int &j, int &k) const |
given a handle, get the corresponding parameters | |
int | is_periodic_i () const |
return whether rectangle is periodic in i | |
int | is_periodic_j () const |
return whether rectangle is periodic in j | |
void | is_periodic (int is_p[2]) const |
int | i_min () const |
convenience functions for parameter extents | |
int | j_min () const |
int | k_min () const |
int | i_max () const |
int | j_max () const |
int | k_max () const |
bool | boundary_complete () const |
test the bounding vertex sequences and determine whether they fully define the vertices covering this element block's parameter space | |
bool | contains (const HomCoord &coords) const |
test whether this sequence contains these parameters | |
bool | contains_vertex (const HomCoord &coords) const |
test whether *vertex parameterization* in this sequence contains these parameters | |
ErrorCode | get_params_connectivity (const int i, const int j, const int k, std::vector< EntityHandle > &connectivity) const |
get connectivity of an entity given entity's parameters | |
ErrorCode | add_vsequence (ScdVertexData *vseq, const HomCoord &p1, const HomCoord &q1, const HomCoord &p2, const HomCoord &q2, const HomCoord &p3, const HomCoord &q3, bool bb_input=false, const HomCoord &bb_min=HomCoord::unitv[0], const HomCoord &bb_max=HomCoord::unitv[0]) |
add a vertex seq ref to this element sequence; if bb_input is true, bounding box (in eseq-local coords) of vseq being added is input in bb_min and bb_max (allows partial sharing of vseq rather than the whole vseq); if it's false, the whole vseq is referenced and the eseq-local coordinates is computed from the transformed bounding box of the vseq | |
SequenceData * | subset (EntityHandle start, EntityHandle end, const int *sequence_data_sizes, const int *tag_data_sizes) const |
unsigned long | get_memory_use () const |
Static Public Member Functions | |
static EntityID | calc_num_entities (EntityHandle start_handle, int irange, int jrange, int krange, int *is_periodic=NULL) |
Private Member Functions | |
ScdElementData () | |
bare constructor, so compiler doesn't create one for me | |
Private Attributes | |
HomCoord | boxParams [3] |
parameter min/max/stride for vertices, in homogeneous coords ijkh; elements are one less than max | |
int | dIJK [3] |
difference between max and min params plus one (i.e. # VERTICES in each parametric direction) | |
int | dIJKm1 [3] |
difference between max and min params (i.e. # ELEMENTS in each parametric direction) | |
int | isPeriodic [2] |
whether scd element rectangle is periodic in i and possibly j | |
std::vector< VertexDataRef > | vertexSeqRefs |
list of bounding vertex blocks |
Definition at line 43 of file ScdElementData.hpp.
moab::ScdElementData::ScdElementData | ( | ) | [private] |
bare constructor, so compiler doesn't create one for me
moab::ScdElementData::ScdElementData | ( | EntityHandle | start_handle, |
const int | imin, | ||
const int | jmin, | ||
const int | kmin, | ||
const int | imax, | ||
const int | jmax, | ||
const int | kmax, | ||
int * | is_periodic | ||
) |
constructor
Definition at line 54 of file ScdElementData.cpp.
References boxParams, dIJK, dIJKm1, and isPeriodic.
: SequenceData( 0,
shandle,
shandle + calc_num_entities( shandle, imax - imin, jmax - jmin, kmax - kmin, is_p ) - 1 )
{
// need to have meaningful parameters
assert( imax >= imin && jmax >= jmin && kmax >= kmin );
isPeriodic[0] = ( is_p ? is_p[0] : 0 );
isPeriodic[1] = ( is_p ? is_p[1] : 0 );
boxParams[0] = HomCoord( imin, jmin, kmin );
boxParams[1] = HomCoord( imax, jmax, kmax );
boxParams[2] = HomCoord( 1, 1, 1 );
// assign and compute parameter stuff
dIJK[0] = boxParams[1][0] - boxParams[0][0] + 1;
dIJK[1] = boxParams[1][1] - boxParams[0][1] + 1;
dIJK[2] = boxParams[1][2] - boxParams[0][2] + 1;
dIJKm1[0] = dIJK[0] - ( isPeriodic[0] ? 0 : 1 );
dIJKm1[1] = dIJK[1] - ( isPeriodic[1] ? 0 : 1 );
dIJKm1[2] = dIJK[2] - 1;
}
moab::ScdElementData::~ScdElementData | ( | ) | [virtual] |
Definition at line 85 of file ScdElementData.cpp.
{}
ErrorCode moab::ScdElementData::add_vsequence | ( | ScdVertexData * | vseq, |
const HomCoord & | p1, | ||
const HomCoord & | q1, | ||
const HomCoord & | p2, | ||
const HomCoord & | q2, | ||
const HomCoord & | p3, | ||
const HomCoord & | q3, | ||
bool | bb_input = false , |
||
const HomCoord & | bb_min = HomCoord::unitv[0] , |
||
const HomCoord & | bb_max = HomCoord::unitv[0] |
||
) | [inline] |
add a vertex seq ref to this element sequence; if bb_input is true, bounding box (in eseq-local coords) of vseq being added is input in bb_min and bb_max (allows partial sharing of vseq rather than the whole vseq); if it's false, the whole vseq is referenced and the eseq-local coordinates is computed from the transformed bounding box of the vseq
Definition at line 312 of file ScdElementData.hpp.
References moab::ScdVertexData::max_params(), MB_SUCCESS, moab::ScdVertexData::min_params(), moab::HomXform::three_pt_xform(), and vertexSeqRefs.
Referenced by moab::ScdBox::add_vbox(), and moab::SequenceManager::add_vsequence().
{
// compute the transform given the vseq-local parameters and the mapping to
// this element sequence's parameters passed in minmax
HomXform M;
M.three_pt_xform( p1, q1, p2, q2, p3, q3 );
// min and max in element seq's parameter system may not be same as those in
// vseq's system, so need to take min/max
HomCoord minmax[2];
if( bb_input )
{
minmax[0] = bb_min;
minmax[1] = bb_max;
}
else
{
minmax[0] = vseq->min_params() * M;
minmax[1] = vseq->max_params() * M;
}
// check against other vseq's to make sure they don't overlap
for( std::vector< VertexDataRef >::const_iterator vsit = vertexSeqRefs.begin(); vsit != vertexSeqRefs.end();
++vsit )
if( ( *vsit ).contains( minmax[0] ) || ( *vsit ).contains( minmax[1] ) ) return MB_FAILURE;
HomCoord tmp_min( std::min( minmax[0].i(), minmax[1].i() ), std::min( minmax[0].j(), minmax[1].j() ),
std::min( minmax[0].k(), minmax[1].k() ) );
HomCoord tmp_max( std::max( minmax[0].i(), minmax[1].i() ), std::max( minmax[0].j(), minmax[1].j() ),
std::max( minmax[0].k(), minmax[1].k() ) );
// set up a new vertex sequence reference
VertexDataRef tmp_seq_ref( tmp_min, tmp_max, M, vseq );
// add to the list
vertexSeqRefs.push_back( tmp_seq_ref );
return MB_SUCCESS;
}
bool moab::ScdElementData::boundary_complete | ( | ) | const |
test the bounding vertex sequences and determine whether they fully define the vertices covering this element block's parameter space
Definition at line 87 of file ScdElementData.cpp.
References boxParams, moab::HomCoord::unitv, and vertexSeqRefs.
Referenced by moab::SweptElementSeq::boundary_complete(), moab::StructuredElementSeq::boundary_complete(), and get_vertex().
{
// test the bounding vertex sequences to see if they fully define the
// vertex parameter space for this rectangular block of elements
int p;
std::vector< VertexDataRef > minlist, maxlist;
// pseudo code:
// for each vertex sequence v:
for( std::vector< VertexDataRef >::const_iterator vseq = vertexSeqRefs.begin(); vseq != vertexSeqRefs.end();
++vseq )
{
// test min corner mincorner:
bool mincorner = true;
// for each p = (i-1,j,k), (i,j-1,k), (i,j,k-1):
for( p = 0; p < 3; p++ )
{
// for each vsequence v' != v:
for( std::vector< VertexDataRef >::const_iterator othervseq = vertexSeqRefs.begin();
othervseq != vertexSeqRefs.end(); ++othervseq )
{
if( othervseq == vseq ) continue;
// if v.min-p contained in v'
if( ( *othervseq ).contains( ( *vseq ).minmax[0] - HomCoord::unitv[p] ) )
{
// mincorner = false
mincorner = false;
break;
}
}
if( !mincorner ) break;
}
bool maxcorner = true;
// for each p = (i-1,j,k), (i,j-1,k), (i,j,k-1):
for( p = 0; p < 3; p++ )
{
// for each vsequence v' != v:
for( std::vector< VertexDataRef >::const_iterator othervseq = vertexSeqRefs.begin();
othervseq != vertexSeqRefs.end(); ++othervseq )
{
if( othervseq == vseq ) continue;
// if v.max+p contained in v'
if( ( *othervseq ).contains( ( *vseq ).minmax[1] + HomCoord::unitv[p] ) )
{
// maxcorner = false
maxcorner = false;
break;
}
}
if( !maxcorner ) break;
}
// if mincorner add to min corner list minlist
if( mincorner ) minlist.push_back( *vseq );
// if maxcorner add to max corner list maxlist
if( maxcorner ) maxlist.push_back( *vseq );
}
//
// if minlist.size = 1 & maxlist.size = 1 & minlist[0] = esequence.min &
// maxlist[0] = esequence.max+(1,1,1)
if( minlist.size() == 1 && maxlist.size() == 1 && minlist[0].minmax[0] == boxParams[0] &&
maxlist[0].minmax[1] == boxParams[1] )
// complete
return true;
// else
return false;
}
EntityID moab::ScdElementData::calc_num_entities | ( | EntityHandle | start_handle, |
int | irange, | ||
int | jrange, | ||
int | krange, | ||
int * | is_periodic = NULL |
||
) | [static] |
Definition at line 27 of file ScdElementData.cpp.
References dim, moab::CN::Dimension(), and moab::TYPE_FROM_HANDLE().
{
size_t result = 1;
auto dim = CN::Dimension( TYPE_FROM_HANDLE( start_handle ) );
switch( dim )
{
case 3:
result *= krange;
// fall through
case 2:
result *= ( is_periodic && is_periodic[1] ? ( jrange + 1 ) : jrange );
// fall through
case 1:
result *= ( is_periodic && is_periodic[0] ? ( irange + 1 ) : irange );
break;
default:
result = 0;
assert( false );
break;
}
return result;
}
bool moab::ScdElementData::contains | ( | const HomCoord & | coords | ) | const [inline] |
test whether this sequence contains these parameters
Definition at line 254 of file ScdElementData.hpp.
References boxParams, dIJKm1, moab::HomCoord::i(), moab::HomCoord::j(), and moab::HomCoord::k().
Referenced by moab::SweptElementSeq::contains(), moab::StructuredElementSeq::contains(), and get_params_connectivity().
{
// upper bound is < instead of <= because element params max is one less
// than vertex params max, except in case of 2d or 1d sequence
return ( ( dIJKm1[0] && temp.i() >= boxParams[0].i() && temp.i() < boxParams[0].i() + dIJKm1[0] ) &&
( ( !dIJKm1[1] && temp.j() == boxParams[1].j() ) ||
( dIJKm1[1] && temp.j() >= boxParams[0].j() && temp.j() < boxParams[0].j() + dIJKm1[1] ) ) &&
( ( !dIJKm1[2] && temp.k() == boxParams[1].k() ) ||
( dIJKm1[2] && temp.k() >= boxParams[0].k() && temp.k() < boxParams[0].k() + dIJKm1[2] ) ) );
}
bool moab::ScdElementData::contains_vertex | ( | const HomCoord & | coords | ) | const [inline] |
test whether *vertex parameterization* in this sequence contains these parameters
Definition at line 265 of file ScdElementData.hpp.
References boxParams, dIJK, moab::HomCoord::i(), moab::HomCoord::j(), and moab::HomCoord::k().
{
// upper bound is < instead of <= because element params max is one less
// than vertex params max, except in case of 2d or 1d sequence
return ( ( dIJK[0] && temp.i() >= boxParams[0].i() && temp.i() < boxParams[0].i() + dIJK[0] ) &&
( ( !dIJK[1] && temp.j() == boxParams[1].j() ) ||
( dIJK[1] && temp.j() >= boxParams[0].j() && temp.j() < boxParams[0].j() + dIJK[1] ) ) &&
( ( !dIJK[2] && temp.k() == boxParams[1].k() ) ||
( dIJK[2] && temp.k() >= boxParams[0].k() && temp.k() < boxParams[0].k() + dIJK[2] ) ) );
}
EntityHandle moab::ScdElementData::get_element | ( | const int | i, |
const int | j, | ||
const int | k | ||
) | const [inline] |
get handle of element at i, j, k
Definition at line 209 of file ScdElementData.hpp.
References dIJKm1, i_min(), j_min(), k_min(), and moab::SequenceData::start_handle().
Referenced by moab::SweptElementSeq::get_element(), and moab::StructuredElementSeq::get_element().
{
return start_handle() + ( i - i_min() ) + ( j - j_min() ) * dIJKm1[0] + ( k - k_min() ) * dIJKm1[0] * dIJKm1[1];
}
unsigned long moab::ScdElementData::get_memory_use | ( | ) | const |
Definition at line 169 of file ScdElementData.cpp.
References vertexSeqRefs.
Referenced by moab::SweptElementSeq::get_const_memory_use(), and moab::StructuredElementSeq::get_const_memory_use().
{
return sizeof( *this ) + vertexSeqRefs.capacity() * sizeof( VertexDataRef );
}
ErrorCode moab::ScdElementData::get_params | ( | const EntityHandle | ehandle, |
int & | i, | ||
int & | j, | ||
int & | k | ||
) | const [inline] |
given a handle, get the corresponding parameters
Definition at line 232 of file ScdElementData.hpp.
References boxParams, dIJKm1, moab::HomCoord::i(), i_max(), i_min(), moab::HomCoord::j(), j_max(), j_min(), moab::HomCoord::k(), k_max(), k_min(), MB_SUCCESS, moab::SequenceData::size(), moab::SequenceData::start_handle(), and moab::TYPE_FROM_HANDLE().
Referenced by moab::SweptElementSeq::get_params(), and moab::StructuredElementSeq::get_params().
{
if( TYPE_FROM_HANDLE( ehandle ) != TYPE_FROM_HANDLE( start_handle() ) ) return MB_FAILURE;
int hdiff = ehandle - start_handle();
// use double ?: test below because on some platforms, both sides of the : are
// evaluated, and if dIJKm1[1] is zero, that'll generate a divide-by-zero
k = ( dIJKm1[1] > 0 ? hdiff / ( dIJKm1[1] > 0 ? dIJKm1[0] * dIJKm1[1] : 1 ) : 0 );
j = ( hdiff - ( k * dIJKm1[0] * dIJKm1[1] ) ) / dIJKm1[0];
i = hdiff % dIJKm1[0];
k += boxParams[0].k();
j += boxParams[0].j();
i += boxParams[0].i();
return ( ehandle >= start_handle() && ehandle < start_handle() + size() && i >= i_min() && i <= i_max() &&
j >= j_min() && j <= j_max() && k >= k_min() && k <= k_max() )
? MB_SUCCESS
: MB_FAILURE;
}
ErrorCode moab::ScdElementData::get_params_connectivity | ( | const int | i, |
const int | j, | ||
const int | k, | ||
std::vector< EntityHandle > & | connectivity | ||
) | const [inline] |
get connectivity of an entity given entity's parameters
Definition at line 362 of file ScdElementData.hpp.
References contains(), dIJKm1, moab::CN::Dimension(), get_vertex(), isPeriodic, MB_SUCCESS, moab::SequenceData::start_handle(), and moab::TYPE_FROM_HANDLE().
Referenced by moab::SweptElementSeq::get_params_connectivity(), and moab::StructuredElementSeq::get_params_connectivity().
{
if( contains( HomCoord( i, j, k ) ) == false ) return MB_FAILURE;
int ip1 = ( isPeriodic[0] ? ( i + 1 ) % dIJKm1[0] : i + 1 ),
jp1 = ( isPeriodic[1] ? ( j + 1 ) % dIJKm1[1] : j + 1 );
connectivity.push_back( get_vertex( i, j, k ) );
connectivity.push_back( get_vertex( ip1, j, k ) );
if( CN::Dimension( TYPE_FROM_HANDLE( start_handle() ) ) < 2 ) return MB_SUCCESS;
connectivity.push_back( get_vertex( ip1, jp1, k ) );
connectivity.push_back( get_vertex( i, jp1, k ) );
if( CN::Dimension( TYPE_FROM_HANDLE( start_handle() ) ) < 3 ) return MB_SUCCESS;
connectivity.push_back( get_vertex( i, j, k + 1 ) );
connectivity.push_back( get_vertex( ip1, j, k + 1 ) );
connectivity.push_back( get_vertex( ip1, jp1, k + 1 ) );
connectivity.push_back( get_vertex( i, jp1, k + 1 ) );
return MB_SUCCESS;
}
EntityHandle moab::ScdElementData::get_vertex | ( | const HomCoord & | coords | ) | const [inline] |
get handle of vertex at homogeneous coords
Definition at line 291 of file ScdElementData.hpp.
References boundary_complete(), and vertexSeqRefs.
Referenced by get_params_connectivity(), moab::SweptElementSeq::get_vertex(), moab::StructuredElementSeq::get_vertex(), and get_vertex().
{
assert( boundary_complete() );
for( std::vector< VertexDataRef >::const_iterator it = vertexSeqRefs.begin(); it != vertexSeqRefs.end(); ++it )
{
if( ( *it ).minmax[0] <= coords && ( *it ).minmax[1] >= coords )
{
// first get the vertex block-local parameters
HomCoord local_coords = coords / ( *it ).xform;
assert( ( *it ).srcSeq->contains( local_coords ) );
// now get the vertex handle for those coords
return ( *it ).srcSeq->get_vertex( local_coords );
}
}
// got here, it's an error
return 0;
}
EntityHandle moab::ScdElementData::get_vertex | ( | int | i, |
int | j, | ||
int | k | ||
) | const [inline] |
Definition at line 99 of file ScdElementData.hpp.
References get_vertex().
{
return get_vertex( HomCoord( i, j, k ) );
}
int moab::ScdElementData::i_max | ( | ) | const [inline] |
Definition at line 150 of file ScdElementData.hpp.
References boxParams.
Referenced by get_params().
{
return ( boxParams[1].hom_coord() )[0];
}
int moab::ScdElementData::i_min | ( | ) | const [inline] |
convenience functions for parameter extents
Definition at line 138 of file ScdElementData.hpp.
References boxParams.
Referenced by get_element(), and get_params().
{
return ( boxParams[0].hom_coord() )[0];
}
void moab::ScdElementData::is_periodic | ( | int | is_p[2] | ) | const [inline] |
Definition at line 131 of file ScdElementData.hpp.
References isPeriodic.
Referenced by moab::StructuredElementSeq::is_periodic().
{
is_p[0] = isPeriodic[0];
is_p[1] = isPeriodic[1];
}
int moab::ScdElementData::is_periodic_i | ( | ) | const [inline] |
return whether rectangle is periodic in i
Definition at line 120 of file ScdElementData.hpp.
References isPeriodic.
Referenced by moab::StructuredElementSeq::is_periodic_i().
{
return isPeriodic[0];
};
int moab::ScdElementData::is_periodic_j | ( | ) | const [inline] |
return whether rectangle is periodic in j
Definition at line 126 of file ScdElementData.hpp.
References isPeriodic.
Referenced by moab::StructuredElementSeq::is_periodic_j().
{
return isPeriodic[1];
};
int moab::ScdElementData::j_max | ( | ) | const [inline] |
Definition at line 154 of file ScdElementData.hpp.
References boxParams.
Referenced by get_params().
{
return ( boxParams[1].hom_coord() )[1];
}
int moab::ScdElementData::j_min | ( | ) | const [inline] |
Definition at line 142 of file ScdElementData.hpp.
References boxParams.
Referenced by get_element(), and get_params().
{
return ( boxParams[0].hom_coord() )[1];
}
int moab::ScdElementData::k_max | ( | ) | const [inline] |
Definition at line 158 of file ScdElementData.hpp.
References boxParams.
Referenced by get_params().
{
return ( boxParams[1].hom_coord() )[2];
}
int moab::ScdElementData::k_min | ( | ) | const [inline] |
Definition at line 146 of file ScdElementData.hpp.
References boxParams.
Referenced by get_element(), and get_params().
{
return ( boxParams[0].hom_coord() )[2];
}
const HomCoord & moab::ScdElementData::max_params | ( | ) | const [inline] |
get max params for this element
Definition at line 219 of file ScdElementData.hpp.
References boxParams.
Referenced by moab::SweptElementSeq::max_params(), moab::StructuredElementSeq::max_params(), and moab::ScdBox::ScdBox().
{
return boxParams[1];
}
const HomCoord & moab::ScdElementData::min_params | ( | ) | const [inline] |
get min params for this element
Definition at line 214 of file ScdElementData.hpp.
References boxParams.
Referenced by moab::SweptElementSeq::min_params(), moab::StructuredElementSeq::min_params(), and moab::ScdBox::ScdBox().
{
return boxParams[0];
}
void moab::ScdElementData::param_extents | ( | int & | di, |
int & | dj, | ||
int & | dk | ||
) | const [inline] |
get the number of vertices in each direction, inclusive
Definition at line 225 of file ScdElementData.hpp.
References dIJK.
Referenced by moab::SweptElementSeq::param_extents(), and moab::StructuredElementSeq::param_extents().
{
di = dIJK[0];
dj = dIJK[1];
dk = dIJK[2];
}
SequenceData * moab::ScdElementData::subset | ( | EntityHandle | start, |
EntityHandle | end, | ||
const int * | sequence_data_sizes, | ||
const int * | tag_data_sizes | ||
) | const |
Definition at line 161 of file ScdElementData.cpp.
{
return 0;
}
HomCoord moab::ScdElementData::boxParams[3] [private] |
parameter min/max/stride for vertices, in homogeneous coords ijkh; elements are one less than max
Definition at line 65 of file ScdElementData.hpp.
Referenced by boundary_complete(), contains(), contains_vertex(), get_params(), i_max(), i_min(), j_max(), j_min(), k_max(), k_min(), max_params(), min_params(), and ScdElementData().
int moab::ScdElementData::dIJK[3] [private] |
difference between max and min params plus one (i.e. # VERTICES in each parametric direction)
Definition at line 69 of file ScdElementData.hpp.
Referenced by contains_vertex(), param_extents(), and ScdElementData().
int moab::ScdElementData::dIJKm1[3] [private] |
difference between max and min params (i.e. # ELEMENTS in each parametric direction)
Definition at line 73 of file ScdElementData.hpp.
Referenced by contains(), get_element(), get_params(), get_params_connectivity(), and ScdElementData().
int moab::ScdElementData::isPeriodic[2] [private] |
whether scd element rectangle is periodic in i and possibly j
Definition at line 76 of file ScdElementData.hpp.
Referenced by get_params_connectivity(), is_periodic(), is_periodic_i(), is_periodic_j(), and ScdElementData().
std::vector< VertexDataRef > moab::ScdElementData::vertexSeqRefs [private] |
list of bounding vertex blocks
Definition at line 82 of file ScdElementData.hpp.
Referenced by add_vsequence(), boundary_complete(), get_memory_use(), and get_vertex().