cgma
FBPolyhedron Class Reference

#include <FBPolyhedron.hpp>

List of all members.

Public Member Functions

 FBPolyhedron ()
 ~FBPolyhedron ()
CubitStatus makepoly (const std::vector< double > &coords, const std::vector< int > &connections, std::vector< int > *f_c_indices)
bool boxintersection (FBPolyhedron *otherpoly)
CubitStatus retriangulate (std::vector< int > &newfacets, std::vector< int > &newfacetsindex)
CubitStatus retriangulate (std::vector< int > &newfacets)
bool edge_exists_in_tri (FB_Triangle &tri, int v0, int v1)

Private Member Functions

int makeahashvaluefrom_coord (double x, double y, double z)
int addavertex (double x, double y, double z)
void putnewpoints (std::vector< double > &newpoints)
void putedges (std::vector< int > &newedges)
bool edge_exists (int v0, int v1)
void add_new_triangle_data ()
void make_tri_plane_coeffs (FB_Triangle *tri)
void make_tri_boundingbox (FB_Triangle *tri)
void removeduddedtriangles ()
bool min_max_angles_in_fb_triangle (FB_Triangle *tri, double &min_angle, double &max_angle)
 Find the min and max angles in the given Triangle.
bool min_max_angles_in_polyhedron (double &min_angle, double &max_angle)
void debug_draw_fb_triangle (FB_Triangle *tri)
 Debug draw a triangle.
void debug_draw_boundary_edges (int color)
bool remove_small_angles_in_triangle_range (int lower_index, int upper_index)

Private Attributes

std::vector< FB_Coord * > verts
std::vector< FB_Triangle * > tris
std::vector< FB_Edge * > intersection_edges
std::vector< FB_Triangle * > dudded_tris
std::vector< FB_Edge * > orphaned_edges
std::multimap< unsigned int,
unsigned int > 
edgmmap
double polyxmin
double polyxmax
double polyymin
double polyymax
double polyzmin
double polyzmax
IntegerHashhashobj
unsigned int numpts
unsigned int numtris
KDTreekdtree
int original_numtris

Friends

class FBIntersect
class FBClassify
class FBImprint

Detailed Description

Definition at line 36 of file FBPolyhedron.hpp.


Constructor & Destructor Documentation

Definition at line 50 of file FBPolyhedron.cpp.

{
  unsigned int i;

  delete hashobj; 
  delete kdtree;
  for ( i = 0; i < verts.size(); i++ ) {
    delete verts[i];
  } 

  for ( i = 0; i < tris.size(); i++ ) {
    delete tris[i];  
  }

  for ( i = 0; i < dudded_tris.size(); i++ ) {
    delete dudded_tris[i];
  }

  for ( i = 0; i < orphaned_edges.size(); i++ )
  {
    delete orphaned_edges[i];
  }
}

Member Function Documentation

Definition at line 425 of file FBPolyhedron.cpp.

{
unsigned int i;
FB_Triangle *tri;

  for ( i = original_numtris; i < tris.size(); i++ ) {
  
    tri = tris[i];
//  make the bounding box
    make_tri_boundingbox(tri);
//  make the plane coefficients
    make_tri_plane_coeffs(tri);
  }
  
}
int FBPolyhedron::addavertex ( double  x,
double  y,
double  z 
) [private]

Definition at line 170 of file FBPolyhedron.cpp.

{
int hashvalue, i, ifoundit;
int *hasharrayptr, hasharraysize;
double xval, yval, zval;
FB_Coord *mycoord, *newcoord;

  hashvalue = makeahashvaluefrom_coord(x,y,z);
  
  hasharrayptr = hashobj->getHashBin(hashvalue,&hasharraysize);
  
  ifoundit = -1;
  for ( i = 0; i < hasharraysize; i++ ) {
    mycoord = verts[hasharrayptr[i]];
    xval = mycoord->coord[0];
    yval = mycoord->coord[1];
    zval = mycoord->coord[2];
    if ( ( fabs(xval-x) < EPSILON ) && 
    ( fabs(yval-y) < EPSILON ) &&
    ( fabs(zval-z) < EPSILON ) ) {
      ifoundit = hasharrayptr[i];
      break;
    }
  }
  if ( ifoundit == -1 ) {
    newcoord = new FB_Coord(x,y,z);
    ifoundit = verts.size();    
    verts.push_back(newcoord);
    hashobj->addtoHashList(hashvalue,ifoundit);
  }
  verts[ifoundit]->is_on_boundary = true;
  return ifoundit;
  
}

Definition at line 138 of file FBPolyhedron.cpp.

{
double xmin, ymin, zmin, xmax, ymax, zmax;

  xmin = otherpoly->polyxmin;
  ymin = otherpoly->polyymin;
  zmin = otherpoly->polyzmin;
  xmax = otherpoly->polyxmax;
  ymax = otherpoly->polyymax;
  zmax = otherpoly->polyzmax;
  
  if ( (polyxmin > xmax) || (polyxmax < xmin) ) return false;
  if ( (polyymin > ymax) || (polyymax < ymin) ) return false;
  if ( (polyzmin > zmax) || (polyzmax < zmin) ) return false;
   
  return true;
}
void FBPolyhedron::debug_draw_boundary_edges ( int  color) [private]

Debug draw the edges in the polyhedron that are marked boundary (ie, draw the edges that have a cubit_edge).

Definition at line 585 of file FBPolyhedron.cpp.

{
  unsigned int i;
  FB_Triangle* triangle=NULL;
  for(i = 0; i<tris.size(); i++){
    if(!tris[i]->dudded){
      triangle = tris[i];
      unsigned int counter = 0;
      CubitVector vert_0(verts[triangle->v0]->coord[0],
                           verts[triangle->v0]->coord[1],
                         verts[triangle->v0]->coord[2]);
      CubitVector vert_1(verts[triangle->v1]->coord[0],
                         verts[triangle->v1]->coord[1],
                         verts[triangle->v1]->coord[2]);
      CubitVector vert_2(verts[triangle->v2]->coord[0],
                         verts[triangle->v2]->coord[1],
                         verts[triangle->v2]->coord[2]);
      if(triangle->cubitedge0index){
        counter++;
        GfxDebug::draw_line(vert_0, vert_1, color);
      }
      if(triangle->cubitedge1index){
        counter++;
        GfxDebug::draw_line(vert_1, vert_2, color);
      }
      if(triangle->cubitedge2index){
        counter++;
        GfxDebug::draw_line(vert_2, vert_0, color);
      }
      if(counter != triangle->edge_list.size()){
        PRINT_WARNING("Possible debug problem.\n");
      }
    }
  }
}

Debug draw a triangle.

Definition at line 622 of file FBPolyhedron.cpp.

{
  CubitVector vert_0(verts[triangle->v0]->coord[0],
                     verts[triangle->v0]->coord[1],
                     verts[triangle->v0]->coord[2]);
  CubitVector vert_1(verts[triangle->v1]->coord[0],
                     verts[triangle->v1]->coord[1],
                     verts[triangle->v1]->coord[2]);
  CubitVector vert_2(verts[triangle->v2]->coord[0],
                     verts[triangle->v2]->coord[1],
                     verts[triangle->v2]->coord[2]);
  GfxDebug::draw_point(vert_0, CUBIT_RED_INDEX);
  GfxDebug::draw_point(vert_1, CUBIT_RED_INDEX);
  GfxDebug::draw_point(vert_2, CUBIT_RED_INDEX);
  GfxDebug::draw_line(vert_0, vert_1, CUBIT_BLUE_INDEX);
  GfxDebug::draw_line(vert_1, vert_2, CUBIT_BLUE_INDEX);
  GfxDebug::draw_line(vert_2, vert_0, CUBIT_BLUE_INDEX);
  int draw_normal = 1;
  if(draw_normal){
    CubitVector center_pos = (vert_0+vert_1+vert_2)/3.0;
    CubitVector opp_pos = center_pos + (vert_1-vert_0)*(vert_2-vert_0);
    GfxDebug::draw_line(center_pos, opp_pos, CUBIT_WHITE_INDEX);
  }
}
bool FBPolyhedron::edge_exists ( int  v0,
int  v1 
) [private]

Definition at line 348 of file FBPolyhedron.cpp.

{
bool exists = false;
unsigned int i;

  for ( i = 0; i < intersection_edges.size(); i++ ) {
    if ( ( (intersection_edges[i]->v0 == v0) && 
           (intersection_edges[i]->v1 == v1) ) ||
         ( (intersection_edges[i]->v0 == v1) && 
           (intersection_edges[i]->v1 == v0) ) ) {
      exists = true;
      break;
    }   
  }

  return exists;
}
bool FBPolyhedron::edge_exists_in_tri ( FB_Triangle tri,
int  v0,
int  v1 
)

Definition at line 408 of file FBPolyhedron.cpp.

{
FB_Edge *edge;
std::vector<FB_Edge*>::iterator dp;

  dp = tri.edge_list.begin();
  while ( dp != tri.edge_list.end() ) {
    edge = *dp;
    if ( ( ( (edge->v0) == v0 ) && ( (edge->v1) == v1 ) ) ||
         ( ( (edge->v0) == v1 ) && ( (edge->v1) == v0 ) ) )
      return true;
    dp++;  
  }
  
  return false;
}

Definition at line 477 of file FBPolyhedron.cpp.

{
double xmin, ymin, zmin, xmax, ymax, zmax;
int j;
int connections[3];
FB_Coord *mycoord;

     xmin = ymin = zmin = CUBIT_DBL_MAX;
     xmax = ymax = zmax = -xmin;
     connections[0] = tri->v0; connections[1] = tri->v1; connections[2] = tri->v2;
     for ( j = 0; j < 3; j++ ) { // make the bounding box
       mycoord = verts[connections[j]];
       xmin = ( xmin < mycoord->coord[0] ) ? xmin : mycoord->coord[0];
       xmax = ( xmax > mycoord->coord[0] ) ? xmax : mycoord->coord[0];
       ymin = ( ymin < mycoord->coord[1] ) ? ymin : mycoord->coord[1];
       ymax = ( ymax > mycoord->coord[1] ) ? ymax : mycoord->coord[1];
       zmin = ( zmin < mycoord->coord[2] ) ? zmin : mycoord->coord[2];
       zmax = ( zmax > mycoord->coord[2] ) ? zmax : mycoord->coord[2];
     }
     if ( (xmax - xmin) < BOX_CRACK ) {
       xmax += BOX_CRACK; xmin -= BOX_CRACK;
     }
     if ( (ymax - ymin) < BOX_CRACK ) {
       ymax += BOX_CRACK; ymin -= BOX_CRACK;
     }
     if ( (zmax - zmin) < BOX_CRACK ) {
       zmax += BOX_CRACK; zmin -= BOX_CRACK;
     }     
     tri->boundingbox.xmin = xmin; tri->boundingbox.xmax = xmax;
     tri->boundingbox.ymin = ymin; tri->boundingbox.ymax = ymax;
     tri->boundingbox.zmin = zmin; tri->boundingbox.zmax = zmax;

     //  Update the object's bounding box
     polyxmin = ( polyxmin < xmin ) ? polyxmin : xmin;
     polyymin = ( polyymin < ymin ) ? polyymin : ymin;
     polyzmin = ( polyzmin < zmin ) ? polyzmin : zmin;
     polyxmax = ( polyxmax > xmax ) ? polyxmax : xmax;
     polyymax = ( polyymax > ymax ) ? polyymax : ymax;
     polyzmax = ( polyzmax > zmax ) ? polyzmax : zmax;

}

Definition at line 441 of file FBPolyhedron.cpp.

{
FB_Coord *mycoord;
double x1, x2, x3, y1, y2, y3, z1, z2, z3, e1x, e1y, e1z, e2x, e2y, e2z;
double a, b, c, d, dtemp;

     mycoord = verts[tri->v0];
     x1 = mycoord->coord[0];
     y1 = mycoord->coord[1];
     z1 = mycoord->coord[2];
     mycoord = verts[tri->v1];
     x2 = mycoord->coord[0];
     y2 = mycoord->coord[1];
     z2 = mycoord->coord[2];
     mycoord = verts[tri->v2];
     x3 = mycoord->coord[0];
     y3 = mycoord->coord[1];
     z3 = mycoord->coord[2];
     e1x = x1 - x2; e1y = y1 - y2; e1z = z1 - z2;
     e2x = x3 - x2; e2y = y3 - y2; e2z = z3 - z2;
     a = e1z*e2y - e2z*e1y;
     b = e1x*e2z - e2x*e1z;
     c = e1y*e2x - e2y*e1x;
     dtemp = sqrt(a*a + b*b + c*c);
     if ( dtemp > EPSILON2 ) {
       a /= dtemp;
       b /= dtemp;
       c /= dtemp;
     } else {
       PRINT_WARNING("small-area triangle\n");
     }
     d = -(a*x1 + b*y1 + c*z1);
     tri->a = a; tri->b = b; tri->c = c; tri->d = d;

}
int FBPolyhedron::makeahashvaluefrom_coord ( double  x,
double  y,
double  z 
) [private]

Definition at line 156 of file FBPolyhedron.cpp.

{
double mantissasum;

      if ( fabs(x) < 1.e-3 ) x = 0.0;
      if ( fabs(y) < 1.e-3 ) y = 0.0;
      if ( fabs(z) < 1.e-3 ) z = 0.0;
      mantissasum = (int)(10000.0*fabs(x) + 0.5) + 
                    (int)(10000.0*fabs(y) + 0.5) + 
          (int)(10000.0*fabs(z) + 0.5);
      
      return (int)(mantissasum) % NUMHASHBINS;
}
CubitStatus FBPolyhedron::makepoly ( const std::vector< double > &  coords,
const std::vector< int > &  connections,
std::vector< int > *  f_c_indices 
)

Definition at line 77 of file FBPolyhedron.cpp.

{
  int hashvalue, parent, cubitfacetindex;
  int cubitedge0index, cubitedge1index, cubitedge2index;
  unsigned int i;
  FB_Coord *mycoord; 
  FB_Triangle *mytri;
  CubitStatus status;
  FSBOXVECTOR boxvector;
  std::vector<int>::iterator dpi;
  status = CUBIT_SUCCESS;
  
   for ( i = 0; i < coords.size(); i += 3 ) {

      mycoord = new FB_Coord(coords[i],coords[i+1],coords[i+2]); 
      hashvalue = makeahashvaluefrom_coord(coords[i],coords[i+1],coords[i+2]);
      hashobj->addtoHashList(hashvalue,verts.size());
      verts.push_back(mycoord);  
   }
   numpts = verts.size();
   parent = -1;
   if ( f_c_indices ){
     dpi = f_c_indices->begin();
   }
   for ( i = 0; i < connections.size(); i += 3 ) {
       //as a safety check make, not only do we see if we have the
       // f_c_indices vector, but make sure we don't go past the end of
       // that vector.
     if ( f_c_indices  && (i<f_c_indices->size())){
       cubitfacetindex = *dpi;
       cubitedge0index = *(dpi+1);
       cubitedge1index = *(dpi+2);
       cubitedge2index = *(dpi+3);
       dpi += 4;    
     }
     else {
       cubitfacetindex = 0;
       cubitedge0index = cubitedge1index = cubitedge2index = 0;
     }
     mytri = new FB_Triangle(connections[i],connections[i+1],connections[i+2],
                             parent,cubitfacetindex,cubitedge0index,
                             cubitedge1index,cubitedge2index);
     make_tri_boundingbox(mytri);

     boxvector.push_back(&mytri->boundingbox);  // for the KdTree

     make_tri_plane_coeffs(mytri);     
    
     tris.push_back(mytri);
     
   }
 
   numtris = original_numtris = tris.size();
   kdtree = new KDTree();
   kdtree->makeKDTree(numtris,boxvector);
   
   return status;        
}
bool FBPolyhedron::min_max_angles_in_fb_triangle ( FB_Triangle tri,
double &  min_angle,
double &  max_angle 
) [private]

Find the min and max angles in the given Triangle.

Definition at line 540 of file FBPolyhedron.cpp.

{
  CubitVector vert_0(verts[triangle->v0]->coord[0],
                     verts[triangle->v0]->coord[1],
                     verts[triangle->v0]->coord[2]);
  CubitVector vert_1(verts[triangle->v1]->coord[0],
                     verts[triangle->v1]->coord[1],
                     verts[triangle->v1]->coord[2]);
  CubitVector vert_2(verts[triangle->v2]->coord[0],
                     verts[triangle->v2]->coord[1],
                     verts[triangle->v2]->coord[2]);
  CubitVector sides[3];
  sides[0] = vert_1 - vert_0;
  sides[1] = vert_2 - vert_1;
  sides[2]= vert_0 - vert_2;
  if(sides[0].length_squared() < EPSILON ||
     sides[1].length_squared() < EPSILON ||
     sides[2].length_squared() < EPSILON ){
    min_angle =0.0;
    max_angle =180.0;
    return false;
  }
  double curr_angle;
  min_angle = sides[1].interior_angle(-sides[0]);
  max_angle = min_angle;
  curr_angle = sides[2].interior_angle(-sides[1]);
  if(curr_angle<min_angle){
    min_angle=curr_angle;
  }
  if(curr_angle>max_angle){
    min_angle=curr_angle;
  }
  curr_angle = sides[0].interior_angle(-sides[2]);
  if(curr_angle<min_angle){
    min_angle=curr_angle;
  }
  if(curr_angle>max_angle){
    min_angle=curr_angle;
  }
  return true;
  
}
bool FBPolyhedron::min_max_angles_in_polyhedron ( double &  min_angle,
double &  max_angle 
) [private]

Find the min and max angles in the polyhedron (ie, min and max angle in any triangle).

Definition at line 650 of file FBPolyhedron.cpp.

{
  double curr_min;
  double curr_max;
  if(tris.size() < 1){
    min_angle = 0.0;
    max_angle = 180.0;
    return false;
  }
  min_max_angles_in_fb_triangle(tris[0], min_angle, max_angle);
  unsigned int i;
  for(i = 1; i<tris.size(); i++){
    if(!min_max_angles_in_fb_triangle(tris[i], curr_min, curr_max)){
      return false;
    }
    if(curr_min<min_angle){
      min_angle=curr_min;
    }
    if(curr_max>max_angle){
      max_angle=curr_max;
    }
  }
  return true;
}
void FBPolyhedron::putedges ( std::vector< int > &  newedges) [private]

Definition at line 225 of file FBPolyhedron.cpp.

{
std::deque<unsigned int> edgedeque;
std::multimap<unsigned int,unsigned int>::iterator p3;
unsigned int startedge, startvert, edgnum;
std::multimap<unsigned int,unsigned int>::iterator pub;
unsigned int first, second;
unsigned int i, numedgesadded;
int v0, v1;
bool foundone;

  for ( i = 0; i < intersection_edges.size(); i++ ) {
    v0 = intersection_edges[i]->v0;
    v1 = intersection_edges[i]->v1;

    edgnum = i;
    edgmmap.insert(std::pair<const unsigned int, unsigned int>(v0,edgnum));        
    edgmmap.insert(std::pair<const unsigned int, unsigned int>(v1,edgnum));    
  }
    numedgesadded = 0;
    startedge = 0;

  while ( numedgesadded < intersection_edges.size() ) {
    if ( intersection_edges[startedge]->mark == true ) {
      startedge++;
      continue;
    }
    startvert = v0 = intersection_edges[startedge]->v0;
    v1 = intersection_edges[startedge]->v1;
    intersection_edges[startedge]->mark = true;
    numedgesadded += 1;
    edgedeque.push_back(v0);
    edgedeque.push_back(v1);
 
    foundone = true;

    while ( foundone == true ) {       
      p3 = edgmmap.find(v1);
      pub = edgmmap.upper_bound(v1);

      if ( p3 != edgmmap.end() ) {
        do {
          foundone = false;
     first = p3->first;
     second = p3->second;
     if ( (second == startedge) && (first == startvert) ) {
       break;
     }
     if ( (intersection_edges[second]->v0 != v0) &&
          (intersection_edges[second]->v1 != v0) ) {
       if ( intersection_edges[second]->mark == true ) {
         p3++;
         continue;
       }
       if (intersection_edges[second]->v0 == v1 ) {
         v1 = intersection_edges[second]->v1;
         v0 = intersection_edges[second]->v0;
       } else {
         v1 = intersection_edges[second]->v0;
         v0 = intersection_edges[second]->v1;
       }
            intersection_edges[second]->mark = true;
            numedgesadded += 1;
            edgedeque.push_back(v1);
            foundone = true;
     }
     p3++;
        } while ( (p3 != pub) && (foundone == false) );
        if ( foundone == false ) break;
      } 
    } //                     end while ( foundone == true )
    v1 = startvert;
    foundone = true;
    while ( foundone == true ) {       
        p3 = edgmmap.find(v1);
        pub = edgmmap.upper_bound(v1);

        if ( p3 != edgmmap.end() ) {
          do {
            foundone = false;
       first = p3->first;
       second = p3->second;
       if ( (intersection_edges[second]->v0 == v1) ||
            (intersection_edges[second]->v1 == v1) ) {
         if ( intersection_edges[second]->mark == true ) {
           p3++;
           continue;
         }
         if (intersection_edges[second]->v0 == v1 ) {
           v1 = intersection_edges[second]->v1;
           v0 = intersection_edges[second]->v0;
         } else {
           v1 = intersection_edges[second]->v0;
           v0 = intersection_edges[second]->v1;
         }
              intersection_edges[second]->mark = true;
         numedgesadded += 1;
              edgedeque.push_front(v1);
         foundone = true;
       }
       p3++;
          } while ( (p3 != pub) && (foundone == false) );
          if ( foundone == false ) break;
        } 
    }
  
    unsigned int size = edgedeque.size();

    if ( size > 0 ) {
      newedges.push_back(size);
      for ( i = 0; i < size; i++ ) {
        newedges.push_back(edgedeque[i]);
      }
    }
    edgedeque.clear();
    startedge++;
  } //                     end while ( numedgesadded < intersection_edges.size() )
 
  newedges.push_back(0);


}
void FBPolyhedron::putnewpoints ( std::vector< double > &  newpoints) [private]

Definition at line 205 of file FBPolyhedron.cpp.

{
double coordinate;
unsigned int i;

//  numpts is the original number of points.  Any new points were added  
//  onto the end of the verts vector.

  if ( verts.size() > numpts ) {
    for ( i = numpts; i < verts.size(); i++ ) {
      coordinate = verts[i]->coord[0];
      newpoints.push_back(coordinate);
      coordinate = verts[i]->coord[1];
      newpoints.push_back(coordinate);
      coordinate = verts[i]->coord[2];
      newpoints.push_back(coordinate);
    }  
  }
}
bool FBPolyhedron::remove_small_angles_in_triangle_range ( int  lower_index,
int  upper_index 
) [private]

This a function that attempts to remove small angles from the polyhedron.

Definition at line 519 of file FBPolyhedron.cpp.

{
unsigned int i;
int j;

  j = 0;
  for ( i = 0; i < tris.size(); i++ ) {
    if ( tris[i]->dudded == true )
    {
      dudded_tris.push_back(tris[i]);
    }
    else
    {
      tris[j] = tris[i];
      j++;
    }
  }
  tris.resize(j);
}
CubitStatus FBPolyhedron::retriangulate ( std::vector< int > &  newfacets,
std::vector< int > &  newfacetsindex 
)

Definition at line 366 of file FBPolyhedron.cpp.

{
CubitStatus status;
FBRetriangulate *retriangulater;
unsigned int i;

  status = CUBIT_SUCCESS;
  for ( i = 0; i < tris.size(); i++ ) {
    if ( tris[i]->dudded == true ) {
      tris[i]->parent = (int)i;
      retriangulater = new FBRetriangulate(verts, tris, newfacets, newfacetsindex); 
      status = retriangulater->retriangulate_this_tri(i, orphaned_edges);
      delete retriangulater;
    }
  
  }
  
  return status;
}
CubitStatus FBPolyhedron::retriangulate ( std::vector< int > &  newfacets)

Definition at line 387 of file FBPolyhedron.cpp.

{
CubitStatus status;
FBRetriangulate *retriangulater;
unsigned int i;

  status = CUBIT_SUCCESS;
  for ( i = 0; i < tris.size(); i++ ) {
    if ( tris[i]->dudded == true ) {
      tris[i]->parent = (int)i;
      retriangulater = new FBRetriangulate(verts, tris, newfacets); 
      status = retriangulater->retriangulate_this_tri(i, orphaned_edges);
      std::vector<FB_Edge*>::iterator iter;
      delete retriangulater;
    }
  
  }
  
  return status;
}

Friends And Related Function Documentation

friend class FBClassify [friend]

Definition at line 39 of file FBPolyhedron.hpp.

friend class FBImprint [friend]

Definition at line 40 of file FBPolyhedron.hpp.

friend class FBIntersect [friend]

Definition at line 38 of file FBPolyhedron.hpp.


Member Data Documentation

std::vector<FB_Triangle*> FBPolyhedron::dudded_tris [private]

Definition at line 58 of file FBPolyhedron.hpp.

std::multimap<unsigned int,unsigned int> FBPolyhedron::edgmmap [private]

Definition at line 60 of file FBPolyhedron.hpp.

Definition at line 62 of file FBPolyhedron.hpp.

std::vector<FB_Edge *> FBPolyhedron::intersection_edges [private]

Definition at line 57 of file FBPolyhedron.hpp.

Definition at line 69 of file FBPolyhedron.hpp.

unsigned int FBPolyhedron::numpts [private]

Definition at line 65 of file FBPolyhedron.hpp.

unsigned int FBPolyhedron::numtris [private]

Definition at line 65 of file FBPolyhedron.hpp.

Definition at line 70 of file FBPolyhedron.hpp.

std::vector<FB_Edge*> FBPolyhedron::orphaned_edges [private]

Definition at line 59 of file FBPolyhedron.hpp.

double FBPolyhedron::polyxmax [private]

Definition at line 61 of file FBPolyhedron.hpp.

double FBPolyhedron::polyxmin [private]

Definition at line 61 of file FBPolyhedron.hpp.

double FBPolyhedron::polyymax [private]

Definition at line 61 of file FBPolyhedron.hpp.

double FBPolyhedron::polyymin [private]

Definition at line 61 of file FBPolyhedron.hpp.

double FBPolyhedron::polyzmax [private]

Definition at line 61 of file FBPolyhedron.hpp.

double FBPolyhedron::polyzmin [private]

Definition at line 61 of file FBPolyhedron.hpp.

std::vector<FB_Triangle *> FBPolyhedron::tris [private]

Definition at line 56 of file FBPolyhedron.hpp.

std::vector<FB_Coord *> FBPolyhedron::verts [private]

Definition at line 55 of file FBPolyhedron.hpp.


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