cgma
CubitPoint Class Reference

#include <CubitPoint.hpp>

Inheritance diagram for CubitPoint:
FacetEntity ToolDataUser CubitPointData FaceterPointData ImprintPointData

List of all members.

Public Member Functions

 CubitPoint ()
virtual ~CubitPoint ()
virtual int id ()=0
virtual void set_id (int)
virtual double x ()=0
virtual double y ()=0
virtual double z ()=0
virtual void set (const CubitVector &pos)=0
virtual void marked (int marked)
virtual int marked ()
virtual CubitVector coordinates () const =0
virtual void coordinates (double point_array[3])=0
virtual void add_facet (CubitFacet *facet)=0
virtual void remove_facet (CubitFacet *facet)=0
virtual int num_adj_facets ()=0
virtual void facets (DLIList< CubitFacet * > &facet_list)=0
virtual void edges (DLIList< CubitFacetEdge * > &edge_list)=0
virtual void points (DLIList< CubitPoint * > &point_list)=0
virtual void tris (DLIList< CubitFacet * > &facet_list)
virtual void normal (CubitVector &surf_norm)
virtual CubitVector normal ()
virtual CubitVectornormal_ptr ()
virtual void set_normal (CubitVector &surf_norm)
virtual void reset_normal ()
virtual void d_coef (const double d_coefficient)
virtual double d_coef ()
virtual double u ()
virtual double v ()
virtual double size ()
virtual void set_uv (double u, double v)
virtual void set_uvs (double u, double v, double s)
virtual CubitStatus get_uv (CubitFacet *facet, double &u, double &v)
virtual CubitStatus get_uvs (CubitFacet *facet, double &u, double &v, double &s)
virtual void du (CubitVector &duvec)
virtual CubitVector du ()
virtual void dv (CubitVector &dvvec)
virtual CubitVector dv ()
virtual double * coef_vector ()
virtual void coef_vector (const CubitMatrix &coef)
virtual CubitStatus merge_points (CubitPoint *cp, CubitBoolean keep_point=CUBIT_FALSE)
void shared_facets (CubitPoint *other_point, CubitFacet *&facet1, CubitFacet *&facet2)
void shared_facets (CubitPoint *other_point, DLIList< CubitFacet * > &result_list)
CubitFacetEdgeshared_edge (CubitPoint *other_point)
void adjacent_points (CubitPoint **adj_points, int &num_adj_points)
void adjacent_points (DLIList< CubitPoint * > &result_list)
CubitBox bounding_box ()
void facets_on_surf (int surf_id, DLIList< CubitFacet * > &facet_list, CubitBoolean &on_internal_boundary)
CubitVector normal (CubitFacet *facet_ptr)
CubitVector normal (CubitQuadFacet *qfacet_ptr)
CubitVector normal (CubitFacetEdge *edge_ptr)
CubitVector tangent (CubitFacetEdge *edge_ptr, double mindot)
void next_feature_edges (CubitFacetEdge *this_edge_ptr, DLIList< CubitFacetEdge * > feature_edge_list)
CubitVector project_to_tangent_plane (CubitVector &pt)
void transform_to_local (CubitVector &glob_vec, CubitVector &loc_vec)
void transform_to_global (CubitVector &loc_vec, CubitVector &glob_vec)
void define_tangent_vectors ()
void get_parents (DLIList< FacetEntity * > &facet_list)
void debug_draw (int color=-1, int flush=1, int draw_uv=0)
void compute_avg_normal ()
CubitFacetEdgeget_edge (CubitPoint *other_point)
void transform (CubitTransformMatrix &tfmat)
void rotate_normal (CubitTransformMatrix &rotmat)
CubitStatus check_inverted_facets (const CubitVector &new_position)
void set_as_feature ()
CubitBoolean is_feature ()

Static Public Member Functions

static void set_box_tol (double tol)

Protected Attributes

int markedFlag
CubitVectorsurfNormal
double dCoef
double uVal
double vVal
double sizeVal
CubitVectorsurfU
CubitVectorsurfV
double * coefVector
IttyBit isFeature

Static Private Attributes

static double boxTol = GEOMETRY_RESABS

Detailed Description

Definition at line 45 of file CubitPoint.hpp.


Constructor & Destructor Documentation

Definition at line 20 of file CubitPoint.cpp.

    : markedFlag(0), surfNormal(NULL), dCoef(0.0), 
      uVal(0.0), vVal(0.0), sizeVal(0.0),
      surfU(NULL), surfV(NULL),
      coefVector(NULL), isFeature(0)
{
}
CubitPoint::~CubitPoint ( ) [virtual]

Definition at line 36 of file CubitPoint.cpp.

{
  if (surfNormal) {
    delete surfNormal;
  }
  if (surfU) {
    delete surfU;
  }
  if (surfV) {
    delete surfV;
  }
  if (coefVector) {
    delete coefVector;
  }
}

Member Function Documentation

virtual void CubitPoint::add_facet ( CubitFacet facet) [pure virtual]
void CubitPoint::adjacent_points ( CubitPoint **  adj_points,
int &  num_adj_points 
)

Definition at line 349 of file CubitPoint.cpp.

{
  int i, j, k, index = -1, nextindex = -1;
  CubitBoolean found;
  CubitFacet *facet;
  num_adj_points = 0;
  DLIList <CubitFacet *> attached_facets;
  facets( attached_facets );
  for(i=0; i<attached_facets.size(); i++) {
    facet = attached_facets.get_and_step();
    found = CUBIT_FALSE;
    for (j=0; j<3 && !found; j++) {
      if (facet->point(j) == this) {
        index = (j+1)%3;
        nextindex = (j+2)%3;
        found = CUBIT_TRUE;
      }
    }
    if (found) {
      found = CUBIT_FALSE;
      adj_points[num_adj_points++] = facet->point(index);
      for (k=0; k<num_adj_points-1 && !found; k++) {
        if(adj_points[k] == facet->point(nextindex)){
          found = CUBIT_TRUE;
        }
      }
      if (!found) {
        adj_points[num_adj_points++] = facet->point(nextindex);
      }
    }
  }
}
void CubitPoint::adjacent_points ( DLIList< CubitPoint * > &  result_list)

Definition at line 382 of file CubitPoint.cpp.

{
  DLIList <CubitFacet *> attached_facets;
  facets( attached_facets );
  for( int i = attached_facets.size(); i--; )
  {
    CubitPoint *pt1, *pt2;
    attached_facets.get_and_step()->opposite_edge( this, pt1, pt2 );
    result.append_unique(pt1);
    result.append_unique(pt2);  
  }
}

Definition at line 593 of file CubitPoint.cpp.

{
  CubitVector ptmin( coordinates().x() - boxTol,
                     coordinates().y() - boxTol,
                     coordinates().z() - boxTol );
  CubitVector ptmax( coordinates().x() + boxTol,
                     coordinates().y() + boxTol,
                     coordinates().z() + boxTol );
  CubitBox ptbox( ptmin, ptmax );
  return ptbox;
}

Definition at line 781 of file CubitPoint.cpp.

{
  DLIList<CubitFacet*> facets;
  this->facets(facets);
  while( facets.size() )
  {
    CubitFacet* facet = facets.pop();
    int index = facet->point_index(this);
    CubitVector corner = facet->point((index+1)%3)->coordinates();
    CubitVector opposite_edge = facet->point((index+2)%3)->coordinates();
    opposite_edge -= corner;
    CubitVector old_edge = corner - coordinates();
    CubitVector new_edge = corner - pos;
    old_edge *= opposite_edge;
    new_edge *= opposite_edge;
    if ( (old_edge % new_edge) <= 0.0 )
      return CUBIT_FAILURE;
  }
  return CUBIT_SUCCESS;
}
double * CubitPoint::coef_vector ( ) [inline, virtual]

Definition at line 269 of file CubitPoint.hpp.

{
  assert (coefVector != NULL);
  return coefVector;
}
void CubitPoint::coef_vector ( const CubitMatrix coef) [inline, virtual]

Definition at line 210 of file CubitPoint.hpp.

{
  if (!coefVector) coefVector = new double[5];
  for (int i=0; i<5; i++) {
    coefVector[i] = input_matrix.get(i,0);
  }
}

Reimplemented in ImprintPointData, FaceterPointData, and CubitPointData.

Definition at line 546 of file CubitPoint.cpp.

{
  int j;
  DLIList<CubitFacet*> adj_facet_list;
  facets(adj_facet_list);
  if (adj_facet_list.size() > 0) {
    CubitVector avg_normal(0.0e0, 0.0e0, 0.0e0);
    double totangle = 0.0e0;

      // weight the normal by the spanning angle at the point

    for (j = 0; j < adj_facet_list.size(); j++)
    {
      CubitFacet* facet = adj_facet_list.get_and_step();
      double angle = facet->angle( this );
      facet->weight( angle );
      totangle += angle;
    }
    for (j = 0; j < adj_facet_list.size(); j++)
    {
      CubitFacet* facet = adj_facet_list.get_and_step();
      CubitVector normal = facet->normal();
      normal.normalize();
      avg_normal += (facet->weight() / totangle) * normal;
    }
    avg_normal.normalize();
    if(!surfNormal) {
      surfNormal = new CubitVector ( avg_normal );
    }
    else
    {
      *surfNormal = avg_normal;
    }
    dCoef = -(this->coordinates()%avg_normal);
  }
}
virtual CubitVector CubitPoint::coordinates ( ) const [pure virtual]
virtual void CubitPoint::coordinates ( double  point_array[3]) [pure virtual]
virtual void CubitPoint::d_coef ( const double  d_coefficient) [inline, virtual]

Definition at line 111 of file CubitPoint.hpp.

{dCoef = d_coefficient;};
virtual double CubitPoint::d_coef ( ) [inline, virtual]

Definition at line 112 of file CubitPoint.hpp.

{return dCoef;};
void CubitPoint::debug_draw ( int  color = -1,
int  flush = 1,
int  draw_uv = 0 
) [virtual]

Implements FacetEntity.

Definition at line 527 of file CubitPoint.cpp.

{
  if ( color == -1 )
    color = CUBIT_YELLOW_INDEX;
  CubitVector vec = this->coordinates();
  GfxDebug::draw_point(vec, color);
  if (flush)
    GfxDebug::flush();
}

Definition at line 403 of file CubitPoint.cpp.

{
  // define orthogonal vectors to the normal that are tangent to the
  // the surface.  Note that du and dv are not defined in any global
  // parametric space - they are only defined locally.  Their directions
  // are defined arbitrarily in the tangent plane by taking the smallest
  // components of the normal vector and setting them to "1" for the
  // du vector and then solving for the other component so that the
  // dot product of du and the normal will be zero.  dv is just the cross
  // product of the normal and du

  CubitVector absnorm, duvec, dvvec;
  CubitVector surf_normal = normal();
  absnorm.x( fabs(surf_normal.x()) );
  absnorm.y( fabs(surf_normal.y()) );
  absnorm.z( fabs(surf_normal.z()) );
  if (absnorm.x() >= absnorm.y() && absnorm.x() >= absnorm.z()) {
    duvec.x( (-surf_normal.y() - surf_normal.z()) / surf_normal.x() );
    duvec.y( 1.0e0 );
    duvec.z( 1.0e0 );
  }
  else if (absnorm.y() >= absnorm.z() ) {
    duvec.x( 1.0e0 );
    duvec.y( (-surf_normal.x() - surf_normal.z()) / surf_normal.y() );
    duvec.z( 1.0e0 );
  }
  else {
    duvec.x( 1.0e0 );
    duvec.y( 1.0e0 );
    duvec.z( (-surf_normal.x() - surf_normal.y()) / surf_normal.z() );
  }
  duvec.normalize();
  dvvec = surf_normal * duvec;

//  //sjowen debug
//  CubitVector test = dvvec * surf_normal;
//  double dot = test % duvec;
//  if (dot < 0.999999999999) {
//    PRINT_ERROR("Error in define_tangent_vectors");
//  }

  du( duvec );
  dv( dvvec );
}
void CubitPoint::du ( CubitVector duvec) [inline, virtual]

Definition at line 245 of file CubitPoint.hpp.

{
  if(!surfU) surfU = new CubitVector (duvec);
  else *surfU = duvec;
}
CubitVector CubitPoint::du ( ) [inline, virtual]

Definition at line 257 of file CubitPoint.hpp.

{
  assert(surfU != NULL);
  return *surfU;
}
void CubitPoint::dv ( CubitVector dvvec) [inline, virtual]

Definition at line 251 of file CubitPoint.hpp.

{
  if(!surfV) surfV = new CubitVector (dvvec);
  else *surfV = dvvec;
}
CubitVector CubitPoint::dv ( ) [inline, virtual]

Definition at line 263 of file CubitPoint.hpp.

{
  assert(surfV != NULL);
  return *surfV;
}
virtual void CubitPoint::edges ( DLIList< CubitFacetEdge * > &  edge_list) [pure virtual]

Implements FacetEntity.

Implemented in ImprintPointData, FaceterPointData, and CubitPointData.

virtual void CubitPoint::facets ( DLIList< CubitFacet * > &  facet_list) [pure virtual]

Implements FacetEntity.

Implemented in ImprintPointData, FaceterPointData, and CubitPointData.

void CubitPoint::facets_on_surf ( int  surf_id,
DLIList< CubitFacet * > &  facet_list,
CubitBoolean on_internal_boundary 
)

Definition at line 615 of file CubitPoint.cpp.

{
  DLIList<CubitFacet *>all_facets;
  facets( all_facets );
  int ii;
  CubitFacet *facet;
  for (ii=0; ii<all_facets.size(); ii++)
  {
    facet = all_facets.get_and_step();
    if (facet->tool_id() == surf_id)
    {
      facet_list.append( facet );
    }
  } 
  if (facet_list.size() == 0)
  {
    on_internal_boundary = CUBIT_FALSE;
  }
  else if (all_facets.size() == facet_list.size())
  {
    on_internal_boundary = CUBIT_FALSE;
  }
  else
  {
    on_internal_boundary = CUBIT_TRUE;
  }
}

Definition at line 722 of file CubitPoint.cpp.

{
  DLIList<CubitFacetEdge *>edge_list;
  edges(edge_list);
  int ii;
  CubitFacetEdge *edge_ptr;
  for (ii=0; ii<edge_list.size(); ii++)
  {
    edge_ptr = edge_list.get_and_step();
    if (edge_ptr->other_point( this ) == other_point)
      return edge_ptr;
  }
  return (CubitFacetEdge *)NULL;
}
void CubitPoint::get_parents ( DLIList< FacetEntity * > &  facet_list) [virtual]

Implements FacetEntity.

Definition at line 510 of file CubitPoint.cpp.

{
  DLIList<CubitFacetEdge*> edge_list;
  edges( edge_list );
  int ii;
  for (ii=0; ii<edge_list.size(); ii++)
    facet_list.append( edge_list.get_and_step() );
}
CubitStatus CubitPoint::get_uv ( CubitFacet facet,
double &  u,
double &  v 
) [virtual]

Definition at line 654 of file CubitPoint.cpp.

{
  CubitStatus stat;
  TDFacetBoundaryPoint *td_bfp = 
    TDFacetBoundaryPoint::get_facet_boundary_point( this );
  if (!td_bfp)
  {
    u = uVal;
    v = vVal;
    stat = CUBIT_SUCCESS;
  }
  else
  {
    stat = td_bfp->get_uv( facet, u, v );
  }
  return stat;
}
CubitStatus CubitPoint::get_uvs ( CubitFacet facet,
double &  u,
double &  v,
double &  s 
) [virtual]

Definition at line 681 of file CubitPoint.cpp.

{
  CubitStatus stat;
  TDFacetBoundaryPoint *td_bfp = 
    TDFacetBoundaryPoint::get_facet_boundary_point( this );
  if (!td_bfp)
  {
    u = uVal;
    v = vVal;
    s = sizeVal;
    stat = CUBIT_SUCCESS;
  }
  else
  {
    stat = td_bfp->get_uvs( facet, u, v, s );
  }
  return stat;
}
virtual int CubitPoint::id ( ) [pure virtual]

Definition at line 204 of file CubitPoint.hpp.

{return (isFeature ? CUBIT_TRUE : CUBIT_FALSE); }
virtual void CubitPoint::marked ( int  marked) [inline, virtual]

Reimplemented in ImprintPointData, FaceterPointData, and CubitPointData.

Definition at line 88 of file CubitPoint.hpp.

virtual int CubitPoint::marked ( ) [inline, virtual]

Reimplemented in ImprintPointData, FaceterPointData, and CubitPointData.

Definition at line 89 of file CubitPoint.hpp.

{return markedFlag;};

Reimplemented in CubitPointData.

Definition at line 707 of file CubitPoint.cpp.

{
  // this virtual function must be defined in the inheriting class if you get here
  assert(0);
  return CUBIT_FAILURE;
}
void CubitPoint::next_feature_edges ( CubitFacetEdge this_edge_ptr,
DLIList< CubitFacetEdge * >  feature_edge_list 
)

Definition at line 299 of file CubitPoint.cpp.

{
  //CubitFacetEdge *next_edge_ptr = NULL;

  DLIList<CubitFacetEdge*> edge_list;
  edges( edge_list );
  int ii;

  CubitFacetEdge *edge_ptr = NULL;
  for (ii=0; ii<edge_list.size(); ii++)
  {
    edge_ptr = edge_list.get_and_step();
    if (edge_ptr != this_edge_ptr)
    {
      if (edge_ptr->is_feature())
      {
        feature_edge_list.append(edge_ptr);
      }
    }
  }
}
void CubitPoint::normal ( CubitVector surf_norm) [inline, virtual]

Definition at line 218 of file CubitPoint.hpp.

{
  if(!surfNormal) surfNormal = new CubitVector (surf_norm);
  else *surfNormal = surf_norm;
}
CubitVector CubitPoint::normal ( ) [inline, virtual]

Definition at line 224 of file CubitPoint.hpp.

{
  if (surfNormal==NULL) compute_avg_normal();
  return *surfNormal;
}

Definition at line 144 of file CubitPoint.cpp.

{
  TDFacetBoundaryPoint *td_bfp =
    TDFacetBoundaryPoint::get_facet_boundary_point( this );
  if (td_bfp == NULL)
  {
    return normal();
  }
  else
  {
    CubitVector norm;
    td_bfp->get_normal( facet_ptr, norm );
    return norm;
  }
}

Definition at line 169 of file CubitPoint.cpp.

{
  CubitFacet *facet_ptr = qfacet_ptr->get_tri_facet_at_point( this );
  return normal( facet_ptr );
}

Definition at line 185 of file CubitPoint.cpp.

{
  TDFacetBoundaryPoint *td_bfp =
    TDFacetBoundaryPoint::get_facet_boundary_point( this );
  if (td_bfp == NULL)
  {
    return normal();
  }
  else
  {
    CubitVector norm;
    td_bfp->get_normal( edge_ptr, norm );
    return norm;
  }
}
CubitVector * CubitPoint::normal_ptr ( ) [inline, virtual]

Definition at line 230 of file CubitPoint.hpp.

{
  return surfNormal;
}
virtual int CubitPoint::num_adj_facets ( ) [pure virtual]
virtual void CubitPoint::points ( DLIList< CubitPoint * > &  point_list) [pure virtual]

Implements FacetEntity.

Implemented in ImprintPointData, FaceterPointData, and CubitPointData.

Definition at line 331 of file CubitPoint.cpp.

{
  CubitVector surf_normal = normal();
  double dist = (surf_normal)%pt + d_coef();
  CubitVector point_on_plane( pt.x() - surf_normal.x() * dist,
                              pt.y() - surf_normal.y() * dist,
                              pt.z() - surf_normal.z() * dist );
  return point_on_plane;
}
virtual void CubitPoint::remove_facet ( CubitFacet facet) [pure virtual]
void CubitPoint::reset_normal ( ) [inline, virtual]

Definition at line 240 of file CubitPoint.hpp.

Definition at line 760 of file CubitPoint.cpp.

{
  if (surfNormal)
  {
    *surfNormal = rotmat * (*surfNormal); 
  }
  TDFacetBoundaryPoint *td = TDFacetBoundaryPoint::get_facet_boundary_point( this );
  if (td)
  {
    td->rotate_normal(rotmat);
  }
}
virtual void CubitPoint::set ( const CubitVector pos) [pure virtual]
void CubitPoint::set_as_feature ( ) [inline]

Definition at line 203 of file CubitPoint.hpp.

{ isFeature = 1; }
static void CubitPoint::set_box_tol ( double  tol) [inline, static]

Definition at line 207 of file CubitPoint.hpp.

{boxTol = tol;}
virtual void CubitPoint::set_id ( int  ) [inline, virtual]

Reimplemented in CubitPointData.

Definition at line 81 of file CubitPoint.hpp.

{};
void CubitPoint::set_normal ( CubitVector surf_norm) [inline, virtual]

Definition at line 235 of file CubitPoint.hpp.

{
  *surfNormal = surf_norm;
}
virtual void CubitPoint::set_uv ( double  u,
double  v 
) [inline, virtual]

Definition at line 117 of file CubitPoint.hpp.

{ uVal = u; vVal = v; };
virtual void CubitPoint::set_uvs ( double  u,
double  v,
double  s 
) [inline, virtual]

Definition at line 118 of file CubitPoint.hpp.

{uVal = u; vVal = v; sizeVal = s; }

Definition at line 112 of file CubitPoint.cpp.

{
  
//   CubitFacetEdge *edge = NULL;
//   DLIList <CubitFacetEdge *> attached_edges;
//   edges( attached_edges );
//   if( attached_edges.size() > 0 && other_pt )
//   {
//     for( int i = attached_edges.size(); i > 0; i-- )
//     {
//       edge = attached_edges.get_and_step();
//       if( edge->contains( other_pt ) )
//       {
//         return edge;
//       }
//     }
//   }
//   return edge;
  return get_edge(other_pt);

}
void CubitPoint::shared_facets ( CubitPoint other_point,
CubitFacet *&  facet1,
CubitFacet *&  facet2 
)

Definition at line 60 of file CubitPoint.cpp.

{
  f1 = f2 = 0;
  DLIList <CubitFacet *> attached_facets;
  facets( attached_facets );
  if( attached_facets.size() > 0 && other_pt )
  {
    for( int i = attached_facets.size(); i > 0; i-- )
    {
      CubitFacet* facet = attached_facets.get_and_step();
      if( facet->contains( other_pt ) )
      {
        //three facets??
        assert( !f2 );
        
        if( f1 ) f2 = facet;
        else f1 = facet;
      }
    }
  }
}
void CubitPoint::shared_facets ( CubitPoint other_point,
DLIList< CubitFacet * > &  result_list 
)

Definition at line 84 of file CubitPoint.cpp.

{
  DLIList <CubitFacet *> attached_facets;
  facets( attached_facets );
  if( attached_facets.size() > 0 && other_pt )
  {
    for( int i = attached_facets.size(); i > 0; i-- )
    {
      CubitFacet* facet = attached_facets.get_and_step();
      if( facet->contains( other_pt ) )
      {
        result_set.append( facet );
      }
    }
  }
}
virtual double CubitPoint::size ( ) [inline, virtual]

Definition at line 116 of file CubitPoint.hpp.

{return sizeVal; };
CubitVector CubitPoint::tangent ( CubitFacetEdge edge_ptr,
double  mindot 
)

Definition at line 211 of file CubitPoint.cpp.

{
  CubitPoint *p0 = edge_ptr->point( 0 );
  CubitPoint *p1 = edge_ptr->point( 1 );
  int ii;

  assert( p0 == this || p1 == this ); // the point isn't on the edge

    // if this isn't a feature edge, just return the tangent vector of 
    // the edge.  Otherwise compute the tangent based on neighboring 
    // feature edges

  CubitVector pt_tangent;
  if (!edge_ptr->is_feature())
  {
    CubitVector tmp_vec = coordinates();
    
    edge_ptr->edge_tangent( tmp_vec, pt_tangent );
  }
  else
  {
    // compute tangent for feature edge at previous
    if (p0 == this)
    {
      CubitFacetEdge *prev_edge;
      DLIList <CubitFacetEdge *>feature_edge_list;
      next_feature_edges( edge_ptr, feature_edge_list );

      // average the edges that meet the min_dot criteria
      CubitVector e1 = p1->coordinates() - p0->coordinates();         
      pt_tangent = e1;
      e1.normalize();
      for (ii=0; ii<feature_edge_list.size(); ii++)
      {
        prev_edge = feature_edge_list.get_and_step();
        CubitPoint *p2 = prev_edge->other_point( p0 );          
        CubitVector e0 = p0->coordinates() - p2->coordinates();     
        e0.normalize();
        if (e0 % e1 >= min_dot)
        {
          pt_tangent += (p0->coordinates() - p2->coordinates());
        }
      }
      if (feature_edge_list.size() == 0)
        pt_tangent = e1;
      else
        pt_tangent.normalize();
    }

    // compute tangent for feature edge at next
    else if (p1 == this)
    {
      CubitFacetEdge *next_edge;
      DLIList <CubitFacetEdge *>feature_edge_list;
      next_feature_edges( edge_ptr, feature_edge_list );

      // average the edges that meet the min_dot criteria
      CubitVector e1 = p1->coordinates() - p0->coordinates();         
      pt_tangent = e1;
      e1.normalize();
      for (ii=0; ii<feature_edge_list.size(); ii++)
      {
        next_edge = feature_edge_list.get_and_step();
        CubitPoint *p2 = next_edge->other_point( p1 );          
        CubitVector e0 = p2->coordinates() - p1->coordinates();     
        e0.normalize();
        if (e0 % e1 >= min_dot)
        {
          pt_tangent += (p2->coordinates() - p1->coordinates());
        }
      }
      if (feature_edge_list.size() == 0)
        pt_tangent = e1;
      else
        pt_tangent.normalize();
    }
  } 
  return pt_tangent;
}

Definition at line 745 of file CubitPoint.cpp.

{
  CubitVector loc;
  loc = tfmat * coordinates();
  set( loc );
}
void CubitPoint::transform_to_global ( CubitVector loc_vec,
CubitVector glob_vec 
)

Definition at line 478 of file CubitPoint.cpp.

{
  // Multiply by transformation matrix

  CubitVector vect;
  CubitVector surf_u = du();
  CubitVector surf_v = dv();
  CubitVector surf_normal = normal();
  vect.x( loc_vec.x() * surf_u.x ()+
          loc_vec.y() * surf_v.x() +
          loc_vec.z() * surf_normal.x() );
  vect.y( loc_vec.x() * surf_u.y ()+
          loc_vec.y() * surf_v.y() +
          loc_vec.z() * surf_normal.y() );
  vect.z( loc_vec.x() * surf_u.z ()+
          loc_vec.y() * surf_v.z() +
          loc_vec.z() * surf_normal.z() );

  // Translate from origin

  glob_vec = vect + this->coordinates();
}
void CubitPoint::transform_to_local ( CubitVector glob_vec,
CubitVector loc_vec 
)

Definition at line 456 of file CubitPoint.cpp.

{
  // Translate to local origin at point

  CubitVector vect = glob_vec - this->coordinates();

  // Multiply by transpose (inverse) of transformation vector */

  loc_vec.x( vect % du() );
  loc_vec.y( vect % dv() );
  loc_vec.z( vect % normal() );
}
virtual void CubitPoint::tris ( DLIList< CubitFacet * > &  facet_list) [inline, virtual]

Reimplemented in CubitPointData.

Definition at line 102 of file CubitPoint.hpp.

{ facets(facet_list); }
virtual double CubitPoint::u ( ) [inline, virtual]

Definition at line 114 of file CubitPoint.hpp.

{ return uVal; };
virtual double CubitPoint::v ( ) [inline, virtual]

Definition at line 115 of file CubitPoint.hpp.

{ return vVal; };
virtual double CubitPoint::x ( ) [pure virtual]
virtual double CubitPoint::y ( ) [pure virtual]
virtual double CubitPoint::z ( ) [pure virtual]

Member Data Documentation

double CubitPoint::boxTol = GEOMETRY_RESABS [static, private]

Definition at line 49 of file CubitPoint.hpp.

double* CubitPoint::coefVector [protected]

Definition at line 69 of file CubitPoint.hpp.

double CubitPoint::dCoef [protected]

Definition at line 59 of file CubitPoint.hpp.

Definition at line 72 of file CubitPoint.hpp.

int CubitPoint::markedFlag [protected]

Definition at line 53 of file CubitPoint.hpp.

double CubitPoint::sizeVal [protected]

Definition at line 63 of file CubitPoint.hpp.

Definition at line 56 of file CubitPoint.hpp.

Definition at line 66 of file CubitPoint.hpp.

Definition at line 66 of file CubitPoint.hpp.

double CubitPoint::uVal [protected]

Definition at line 63 of file CubitPoint.hpp.

double CubitPoint::vVal [protected]

Definition at line 63 of file CubitPoint.hpp.


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