moab
moab::SweptElementData Class Reference

#include <SweptElementData.hpp>

Inheritance diagram for moab::SweptElementData:
moab::SequenceData

List of all members.

Classes

class  VertexDataRef
 structure to hold references to bounding vertex blocks More...

Public Member Functions

 SweptElementData (EntityHandle start_handle, const int imin, const int jmin, const int kmin, const int imax, const int jmax, const int kmax, const int *Cq)
 constructor
virtual ~SweptElementData ()
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 HomCoordmin_params () const
 get min params for this element
const HomCoordmax_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 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
bool contains (const HomCoord &coords) const
 test whether 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 (SweptVertexData *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])
SequenceDatasubset (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)

Private Member Functions

 SweptElementData ()
 bare constructor, so compiler doesn't create one for me

Private Attributes

HomCoord elementParams [3]
 parameter min/max/stride, in homogeneous coords ijkh
int dIJK [3]
int dIJKm1 [3]
std::vector< VertexDataRefvertexSeqRefs
 list of bounding vertex blocks

Detailed Description

Definition at line 42 of file SweptElementData.hpp.


Constructor & Destructor Documentation

bare constructor, so compiler doesn't create one for me

moab::SweptElementData::SweptElementData ( EntityHandle  start_handle,
const int  imin,
const int  jmin,
const int  kmin,
const int  imax,
const int  jmax,
const int  kmax,
const int *  Cq 
)

constructor

Definition at line 39 of file SweptElementData.cpp.

    : SequenceData(0, shandle,
                   shandle + 
                   calc_num_entities( shandle, imax-imin, jmax-jmin, kmax-kmin )
                   - 1)
{
    // need to have meaningful parameters
  assert(imax >= imin && jmax >= jmin && kmax >= kmin);
    
  elementParams[0] = HomCoord(imin, jmin, kmin);
  elementParams[1] = HomCoord(imax, jmax, kmax);
  elementParams[2] = HomCoord(1, 1, 1);
  
    // assign and compute parameter stuff
  dIJK[0] = elementParams[1][0] - elementParams[0][0] + 1;
  dIJK[1] = elementParams[1][1] - elementParams[0][1] + 1;
  dIJK[2] = elementParams[1][2] - elementParams[0][2] + 1;
  dIJKm1[0] = dIJK[0] - 1;
  dIJKm1[1] = dIJK[1] - 1;
  dIJKm1[2] = dIJK[2] - 1;
}

Definition at line 65 of file SweptElementData.cpp.

{
}

Member Function Documentation

ErrorCode moab::SweptElementData::add_vsequence ( SweptVertexData 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 242 of file SweptElementData.hpp.

{
    // 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;
}

test the bounding vertex sequences and determine whether they fully define the vertices covering this element block's parameter space

Definition at line 69 of file SweptElementData.cpp.

{
    // 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] == elementParams[0] && 
      maxlist[0].minmax[1] == elementParams[1])
      //   complete
    return true;
    // else

  return false;
}
EntityID moab::SweptElementData::calc_num_entities ( EntityHandle  start_handle,
int  irange,
int  jrange,
int  krange 
) [static]

Definition at line 26 of file SweptElementData.cpp.

{
  size_t result = 1;
  switch (CN::Dimension(TYPE_FROM_HANDLE(start_handle))) {
    default: result = 0; assert( false ); 
    case 3: result *= krange;
    case 2: result *= jrange;
    case 1: result *= irange;
  }
  return result;
}
bool moab::SweptElementData::contains ( const HomCoord coords) const [inline]

test whether this sequence contains these parameters

Definition at line 204 of file SweptElementData.hpp.

{
    // upper bound is < instead of <= because element params max is one less
    // than vertex params max
  return (temp >= elementParams[0] && temp < elementParams[1]);
}
EntityHandle moab::SweptElementData::get_element ( const int  i,
const int  j,
const int  k 
) const [inline]

get handle of element at i, j, k

Definition at line 157 of file SweptElementData.hpp.

{
  return start_handle() + (i-i_min()) + (j-j_min())*dIJKm1[0] + (k-k_min())*dIJKm1[0]*dIJKm1[1];
}
unsigned long moab::SweptElementData::get_memory_use ( ) const

Definition at line 149 of file SweptElementData.cpp.

{
  return sizeof(*this) + vertexSeqRefs.capacity() * sizeof(VertexDataRef);
}
ErrorCode moab::SweptElementData::get_params ( const EntityHandle  ehandle,
int &  i,
int &  j,
int &  k 
) const [inline]

given a handle, get the corresponding parameters

Definition at line 180 of file SweptElementData.hpp.

{
  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 += elementParams[0].k();
  j += elementParams[0].j();
  i += elementParams[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::SweptElementData::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 291 of file SweptElementData.hpp.

{
  if (contains(HomCoord(i, j, k)) == false) return MB_FAILURE;
  
  connectivity.push_back(get_vertex(i, j, k));
  connectivity.push_back(get_vertex(i+1, j, k));
  if (CN::Dimension(TYPE_FROM_HANDLE(start_handle())) < 2) return MB_SUCCESS;
  connectivity.push_back(get_vertex(i+1, j+1, k));
  connectivity.push_back(get_vertex(i, j+1, 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(i+1, j, k+1));
  connectivity.push_back(get_vertex(i+1, j+1, k+1));
  connectivity.push_back(get_vertex(i, j+1, k+1));
  return MB_SUCCESS;
}
EntityHandle moab::SweptElementData::get_vertex ( const HomCoord coords) const [inline]

get handle of vertex at homogeneous coords

Definition at line 224 of file SweptElementData.hpp.

{
  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;
    
      // 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::SweptElementData::get_vertex ( int  i,
int  j,
int  k 
) const [inline]

Definition at line 92 of file SweptElementData.hpp.

    { return get_vertex(HomCoord(i,j,k)); }
int moab::SweptElementData::i_max ( ) const [inline]

Definition at line 115 of file SweptElementData.hpp.

{return (elementParams[1].hom_coord())[0];}
int moab::SweptElementData::i_min ( ) const [inline]

convenience functions for parameter extents

Definition at line 112 of file SweptElementData.hpp.

{return (elementParams[0].hom_coord())[0];}
int moab::SweptElementData::j_max ( ) const [inline]

Definition at line 116 of file SweptElementData.hpp.

{return (elementParams[1].hom_coord())[1];}
int moab::SweptElementData::j_min ( ) const [inline]

Definition at line 113 of file SweptElementData.hpp.

{return (elementParams[0].hom_coord())[1];}
int moab::SweptElementData::k_max ( ) const [inline]

Definition at line 117 of file SweptElementData.hpp.

{return (elementParams[1].hom_coord())[2];}
int moab::SweptElementData::k_min ( ) const [inline]

Definition at line 114 of file SweptElementData.hpp.

{return (elementParams[0].hom_coord())[2];}
const HomCoord & moab::SweptElementData::max_params ( ) const [inline]

get max params for this element

Definition at line 167 of file SweptElementData.hpp.

{
  return elementParams[1];
}
const HomCoord & moab::SweptElementData::min_params ( ) const [inline]

get min params for this element

Definition at line 162 of file SweptElementData.hpp.

{
  return elementParams[0];
}
void moab::SweptElementData::param_extents ( int &  di,
int &  dj,
int &  dk 
) const [inline]

get the number of vertices in each direction, inclusive

Definition at line 173 of file SweptElementData.hpp.

{
  di = dIJK[0];
  dj = dIJK[1];
  dk = dIJK[2];
}
SequenceData * moab::SweptElementData::subset ( EntityHandle  start,
EntityHandle  end,
const int *  sequence_data_sizes,
const int *  tag_data_sizes 
) const

Definition at line 141 of file SweptElementData.cpp.

{
  return 0;
}

Member Data Documentation

int moab::SweptElementData::dIJK[3] [private]

difference between max and min params plus one (i.e. # VERTICES in each parametric direction)

Definition at line 68 of file SweptElementData.hpp.

difference between max and min params (i.e. # ELEMENTS in each parametric direction)

Definition at line 72 of file SweptElementData.hpp.

parameter min/max/stride, in homogeneous coords ijkh

Definition at line 64 of file SweptElementData.hpp.

list of bounding vertex blocks

Definition at line 78 of file SweptElementData.hpp.


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