cgma
RefEdge Class Reference

#include <RefEdge.hpp>

Inheritance diagram for RefEdge:
BasicTopologyEntity TopologyEntity RefEntity CubitEntity CubitObservable ToolDataUser CubitAttribUser

List of all members.

Public Types

typedef RefVertex ChildType
typedef RefFace ParentType

Public Member Functions

virtual ~RefEdge ()
virtual const char * class_name () const
DagType dag_type () const
 Returns the dag type of this enity.
const std::type_info & entity_type_info () const
 Returns the type info of this enity.
void get_parent_ref_entities (DLIList< RefEntity * > &entity_list)
virtual CubitStatus get_point_direction (CubitVector &origin, CubitVector &direction)
virtual CubitStatus get_center_radius (CubitVector &center, double &radius)
RefVertexstart_vertex ()
RefVertexend_vertex ()
virtual void reverse_topology ()
Chainget_chain_ptr ()
CubitStatus get_co_edges (DLIList< CoEdge * > &co_edges_found_list, RefFace *input_ref_face_ptr=NULL)
double angle_between (RefEdge *other_edge_ptr, RefFace *face_ptr)
CubitSense sense (RefFace *face)
virtual int dimension () const
int num_of_common_ref_face (RefEdge *other_edge)
RefFacecommon_ref_face (RefEdge *otherRefEdgePtr)
int common_ref_faces (RefEdge *otherRefEdgePtr, DLIList< RefFace * > &common_face_list)
RefVertexcommon_ref_vertex (RefEdge *otherRefEdgePtr)
RefVertexcommon_ref_vertex (RefEdge *next_ref_edge, RefFace *ref_face_ptr)
CubitBoolean common_vertices (RefEdge *otherRefEdgePtr, DLIList< RefVertex * > &common_verts)
RefEdgeget_other_curve (RefVertex *common_vertex, RefFace *ref_face_ptr)
CubitStatus get_two_co_edges (RefEdge *next_ref_edge, RefFace *ref_face_ptr, CoEdge *&co_edge_this, CoEdge *&co_edge_next)
RefFaceother_face (RefFace *not_face, RefVolume *ref_volume=NULL)
RefVertexother_vertex (RefVertex *refVertexPtr)
RefVertexclosest_vertex (const CubitVector &point3)
Curveget_curve_ptr ()
Curve const * get_curve_ptr () const
CubitVector start_coordinates ()
CubitVector end_coordinates ()
virtual void move_to_curve (CubitVector &vector)
CubitStatus get_interior_extrema (DLIList< CubitVector * > &interior_points, CubitSense &return_sense) const
CubitStatus closest_point_trimmed (CubitVector const &location, CubitVector &closest_location)
CubitStatus closest_point (CubitVector const &location, CubitVector &closest_location, CubitVector *tangent_ptr=NULL, CubitVector *curvature_ptr=NULL)
CubitPointContainment point_containment (const CubitVector &point)
void tangent (const CubitVector &point, CubitVector &tangent_vec)
CubitStatus tangent (const CubitVector &point, CubitVector &tangent_vec, RefFace *ref_face_ptr)
CubitStatus tangent (const CubitVector &point, CubitVector &tangent_vec, RefEdge *next_ref_edge, RefFace *ref_face_ptr)
virtual double measure ()
virtual CubitString measure_label ()
 Returns the type of measure: (volume, area, length, or N/A)
double get_arc_length ()
double get_arc_length (const CubitVector &point1, const CubitVector &point2)
double get_arc_length (const CubitVector &point1, int whichEnd)
double get_chord_length ()
CubitVector curve_center ()
virtual CubitVector center_point ()
 Return the approximate (spatial) center of this RefEntity.
CubitStatus mid_point (const CubitVector &point1, const CubitVector &point2, CubitVector &midPoint)
CubitStatus mid_point (CubitVector &mid_point)
CubitStatus position_from_fraction (double fraction_along_curve, CubitVector &ouput_position)
double start_param ()
double end_param ()
CubitBoolean get_param_range (double &start_param, double &end_param)
double u_from_position (const CubitVector &input_position)
CubitStatus position_from_u (double u_value, CubitVector &output_position)
double u_from_arc_length (double root_param, double arc_length)
double fraction_from_arc_length (RefVertex *root_vertex, double length)
CubitStatus point_from_arc_length (double root_param, double arc_length, CubitVector &new_point)
CubitStatus point_from_arc_length (const CubitVector &root_point, double arc_length, CubitVector &new_point)
double length_from_u (double parameter1, double parameter2)
CubitBoolean is_periodic ()
CubitBoolean is_periodic (double &period)
int get_mark ()
void set_mark (int set_value)
CubitStatus relative_sense (RefEdge *ref_edge_ptr_2, double tolerance_factor, CubitSense *sense, CubitBoolean &spatially_equal, CubitBoolean force_merge=CUBIT_FALSE)
CubitBoolean about_spatially_equal (RefEdge *ref_edge_ptr_2, double tolerance_factor=1.0, CubitSense *sensePtr=NULL, CubitBoolean notify_refEntity=CUBIT_FALSE)
virtual int validate ()
 Perform checks to see if entity valid.
CubitStatus evaluate_exterior_angle (double *exterior_angle)
CubitBoolean is_tolerant ()
CubitStatus get_graphics (GMem &polyline, double angle_tolerance=0, double distance_tolerance=0, double max_edge_length=0)
void reverse_tangent ()

Static Public Member Functions

static const char * get_class_name ()
static void suppress_edge_length_warning (bool flag)

Protected Member Functions

 RefEdge (Curve *curvePtr)

Private Member Functions

void initialize ()
double get_lower_param_bound ()
double get_upper_param_bound ()
 RefEdge (const RefEdge &)
void operator= (const RefEdge &)

Private Attributes

int refEdgeClone
cBit markedFlag: 2

Static Private Attributes

static bool mSuppressEdgeLengthWarning = false

Friends

class RefEntityFactory

Detailed Description

Definition at line 35 of file RefEdge.hpp.


Member Typedef Documentation

Definition at line 39 of file RefEdge.hpp.

Definition at line 40 of file RefEdge.hpp.


Constructor & Destructor Documentation

RefEdge::~RefEdge ( ) [virtual]

Definition at line 84 of file RefEdge.cpp.

{
}
RefEdge::RefEdge ( Curve curvePtr) [protected]

Definition at line 57 of file RefEdge.cpp.

{
    // Set the GeometryEntity pointer
  if (curvePtr != NULL)
  {
    set_geometry_entity_ptr(curvePtr) ;
  }
  else
  {
    PRINT_ERROR("In the RefEdge(Curve*) constructor\n");
    PRINT_ERROR("       Input Curve pointer is NULL\n");
    assert(CUBIT_FALSE);
  }
  
    // Initialize the member data
  initialize();
}
RefEdge::RefEdge ( const RefEdge ) [private]

Member Function Documentation

CubitBoolean RefEdge::about_spatially_equal ( RefEdge ref_edge_ptr_2,
double  tolerance_factor = 1.0,
CubitSense sensePtr = NULL,
CubitBoolean  notify_refEntity = CUBIT_FALSE 
)

Definition at line 1152 of file RefEdge.cpp.

{
  if( this == ref_edge_ptr_2)
  {
    if (sensePtr)
      *sensePtr = CUBIT_FORWARD;
    if (notify_refEntity)
      remove_compare_data();
    return CUBIT_TRUE;
  }
  
  CubitBoolean spatially_equal = CUBIT_FALSE;
  CubitSense rel_sense = CUBIT_FORWARD;
  CubitStatus stat = CUBIT_SUCCESS;
  stat = relative_sense( ref_edge_ptr_2, tol_factor,
                         &rel_sense, spatially_equal);

  if (stat != CUBIT_SUCCESS || !spatially_equal)
  {
    return CUBIT_FALSE;
  }

  if (sensePtr)
    *sensePtr = rel_sense;

    //compare the start and end vertices to be spatially equal.
  RefVertex* this_start = start_vertex();
  RefVertex* this_end = end_vertex();
  RefVertex* edge2_start = ref_edge_ptr_2->start_vertex();
  RefVertex* edge2_end = ref_edge_ptr_2->end_vertex();
  
    // Swap vertices to simplify things later.
  if (rel_sense == CUBIT_REVERSED)
    std::swap(edge2_start, edge2_end);

    //compare vertex locations unless force_merge is true.
      // closed curve case
    if (this_start == this_end || edge2_start == edge2_end)
    {
      if ((this_start != this_end)   ||
          (edge2_start != edge2_end) ||
          !this_start->about_spatially_equal(edge2_start, tol_factor, CUBIT_FALSE))
        return CUBIT_FALSE;
    }
    else
    {
      if ((this_start == edge2_end) ||
          (this_end == edge2_start) ||
          !this_start->about_spatially_equal(edge2_start, tol_factor, CUBIT_FALSE) ||
          !this_end->about_spatially_equal(edge2_end, tol_factor, CUBIT_FALSE))
        return CUBIT_FALSE;
    }

    //Now if they match report it.
    //Do vertices explicitly here rather than in
    //RefVertex::about_spatially_equal(..) because we don't call 
    //RefVertex::about_spatially_equal(..) if force_merge is true.
  if (notify_refEntity)
  {
    this->comparison_found(ref_edge_ptr_2);
    if (this_start != edge2_start)
      this_start->comparison_found(edge2_start);
    else
      this_start->remove_compare_data();
    if (this_end != edge2_end)
      this_end->comparison_found(edge2_end);
    else
      this_end->remove_compare_data();
  }
  
  return CUBIT_TRUE;
}
double RefEdge::angle_between ( RefEdge other_edge_ptr,
RefFace face_ptr 
)

Definition at line 687 of file RefEdge.cpp.

{
  DLIList<DLIList<RefEdge*> > ref_edge_loop;
  CubitVector vertex_vector, normal_vector;
  CubitVector left_vector, right_vector;

    // Loop through face's loops.
  face_ptr->ref_edge_loops( ref_edge_loop );
  double return_val = CUBIT_DBL_MAX;
  
  for( int i = 0; i<ref_edge_loop.size(); i++ )
  {
      // Look for this edge in the list.
    DLIList<RefEdge*>& ref_edge_list = ref_edge_loop[i];
    if( !ref_edge_list.move_to( this ) )
      continue;
      // Look for other_edge in the list.
    if( !ref_edge_list.move_to( other_edge_ptr ) )
    {
      PRINT_ERROR("Attempted to get angle between two "
                  "RefEdges not on the same Loop.\n");
      return 0.;
    }

      // If other_edge is before this edge in the loop . . .
    if( ref_edge_list.next() == this )
    {   // (other = left, this = right)
        // Get the vector of the the common vertex.
      vertex_vector =
        other_edge_ptr->common_ref_vertex(this, face_ptr)->coordinates();
        // Get the normal vector of the face.
      normal_vector = face_ptr->normal_at( vertex_vector );
        // Get the vectors of both edges.
      other_edge_ptr->tangent( vertex_vector, left_vector,
                               this, face_ptr );
      this->tangent( vertex_vector, right_vector,
                     ref_edge_list.next(2), face_ptr);
    }

      // If this edge is before other_edge in the loop . . .
    else if( ref_edge_list.prev() == this )
    {   // (this = left, other = right)
        // Get the vector of the the common vertex.
      vertex_vector =
          common_ref_vertex(other_edge_ptr,face_ptr)->coordinates();
        // Get the normal vector of the face.
      normal_vector = face_ptr->normal_at( vertex_vector );
        // Get the vectors of both edges.
      this->tangent( vertex_vector, left_vector,
                     other_edge_ptr, face_ptr);
      other_edge_ptr->tangent( vertex_vector, right_vector,
                               ref_edge_list.next(), face_ptr );
    }

      // Otherwise, error.
    else
    {
      PRINT_ERROR("Attempted to get angle between two "
                  "non-consecutive edges.\n");
      return 0.;
    }

      // Return the angle between the vectors.
    return_val = normal_vector.vector_angle( right_vector, -left_vector );
    break;
  }
    // If this far, then error.
  if (CUBIT_DBL_MAX == return_val) {
    PRINT_ERROR("Attempted to get angle between two edges "
                "where one is not found on the given face.\n");
    return_val = 0.0;
  }
  return return_val;
}  

Return the approximate (spatial) center of this RefEntity.

Reimplemented from RefEntity.

Definition at line 422 of file RefEdge.cpp.

{
    // Get the parameter range of the RefEdge
  Curve* curve_ptr = get_curve_ptr();
  
  assert (curve_ptr != NULL);
  
  return curve_ptr->center_point();
}
virtual const char* RefEdge::class_name ( ) const [inline, virtual]

Reimplemented from BasicTopologyEntity.

Definition at line 56 of file RefEdge.hpp.

     {
       return get_class_name();
     }
CubitStatus RefEdge::closest_point ( CubitVector const &  location,
CubitVector closest_location,
CubitVector tangent_ptr = NULL,
CubitVector curvature_ptr = NULL 
)

Definition at line 225 of file RefEdge.cpp.

{
  Curve *curve_ptr = get_curve_ptr();
  CubitStatus result = curve_ptr->closest_point(location, closest_location, 
                                           tangent_ptr, curvature_ptr);
  if (tangent_ptr && curve_ptr->bridge_sense() == CUBIT_REVERSED)
    *tangent_ptr = -(*tangent_ptr);
    
  return result;
}
CubitStatus RefEdge::closest_point_trimmed ( CubitVector const &  location,
CubitVector closest_location 
)

Definition at line 239 of file RefEdge.cpp.

{
  Curve *curve_ptr = get_curve_ptr();
  return curve_ptr->closest_point_trimmed(location, closest_location);
}
RefFace * RefEdge::common_ref_face ( RefEdge otherRefEdgePtr)

Definition at line 798 of file RefEdge.cpp.

{
  DLIList<RefFace*> ref_faces, ref_faces_for_other_edge;
  this->ref_faces(ref_faces);
  other_edge->ref_faces(ref_faces_for_other_edge);

  int i, j;
  for (i = 0; i < ref_faces.size(); i++)
  {
     RefFace * ref_face = ref_faces.get_and_step();
     for (j = 0; j < ref_faces_for_other_edge.size(); j++)
        if (ref_face == ref_faces_for_other_edge.get_and_step())
           return ref_face;
  }

  return NULL;
}
int RefEdge::common_ref_faces ( RefEdge otherRefEdgePtr,
DLIList< RefFace * > &  common_face_list 
)

Definition at line 816 of file RefEdge.cpp.

{
   int nedges = 0;
   DLIList<RefFace*> ref_faces, ref_faces_for_other_edge;
   this->ref_faces(ref_faces);
   input_edge->ref_faces(ref_faces_for_other_edge);                            

   int i, j;
   for (i = 0; i < ref_faces.size(); i++)
   {
      RefFace * ref_face = ref_faces.get_and_step();
      for (j = 0; j < ref_faces_for_other_edge.size(); j++)
         if (ref_face == ref_faces_for_other_edge.get_and_step())
         {
           common_face_list.append(ref_face);
           nedges++;
         }
   }

   return nedges;
}
RefVertex * RefEdge::common_ref_vertex ( RefEdge otherRefEdgePtr)

Definition at line 838 of file RefEdge.cpp.

{
  RefVertex *this_start = start_vertex();
  RefVertex *this_end   = end_vertex();
  RefVertex *other_start = other_edge->start_vertex();
  RefVertex *other_end   = other_edge->end_vertex();
  
  if ( this_start == other_start || this_start == other_end )
  {
    return this_start;
  }
  else if ( this_end == other_start || this_end == other_end )
  {
    return this_end;
  }
  else
  {
    return NULL;
  }
}
RefVertex * RefEdge::common_ref_vertex ( RefEdge next_ref_edge,
RefFace ref_face_ptr 
)

Definition at line 890 of file RefEdge.cpp.

{
  CoEdge *co_edge_this = NULL;
  CoEdge *co_edge_next = NULL;
  
    //First get the two coedges that corrispond to
    //this ref_edge an the next one, with reference to
    //the ref_face_ptr.
  CubitStatus status = get_two_co_edges( next_ref_edge,
                                         ref_face_ptr,
                                         co_edge_this,
                                         co_edge_next );
  if (status == CUBIT_FAILURE )
      return NULL;
  assert(co_edge_this != NULL );
  assert(co_edge_next != NULL );
  RefVertex *common_vertex;
  
    //Now according to the sense get the vertex at the
    //end of this edge (start if reversed).
  if ( co_edge_this->get_sense() == CUBIT_FORWARD )
      common_vertex = end_vertex();
  else
      common_vertex = start_vertex();
  
    //Lets just do a sanitiy check...
  if (common_vertex == NULL ||
      !next_ref_edge->is_directly_related( common_vertex )) {
    
      // let's check for bad sense, then print warning and return
      // correct vertex
    common_vertex = other_vertex(common_vertex);
    if (common_vertex != NULL &&
        next_ref_edge->is_directly_related( common_vertex )) {
      
      PRINT_ERROR(" bad sense between curve %d and surface %d; please"
                  " report this.\n", id(), ref_face_ptr->id());
    }
    else {
      PRINT_ERROR("unable to find common vertex (curves %d, %d)",
                  this->id(), next_ref_edge->id());
      assert(CUBIT_TRUE);
      return (RefVertex*)NULL;
    }
  }
    //Hurray, Success...
  return common_vertex;
}
CubitBoolean RefEdge::common_vertices ( RefEdge otherRefEdgePtr,
DLIList< RefVertex * > &  common_verts 
)

Definition at line 858 of file RefEdge.cpp.

{
  CubitBoolean result = CUBIT_FALSE;
  RefVertex *this_start = start_vertex();
  RefVertex *this_end   = end_vertex();
  RefVertex *other_start = other_edge->start_vertex();
  RefVertex *other_end   = other_edge->end_vertex();
  
  if ( this_start == other_start || this_start == other_end )
  {
    common_verts.append(this_start);
    result = CUBIT_TRUE;
  }
  if ( this_end == other_start || this_end == other_end )
  {
    common_verts.append(this_end);
    result = CUBIT_TRUE;
  }

  return result;
}

Definition at line 1427 of file RefEdge.cpp.

{
  CubitVector p1 = start_vertex()->coordinates();
  CubitVector p2 = end_vertex()->coordinates();
  if ( start_vertex() == end_vertex() )
    {
      mid_point(p2);
    }
  p1 += p2;
  p1 /= 2.0;
  return p1;
}
DagType RefEdge::dag_type ( ) const [inline, virtual]

Returns the dag type of this enity.

Implements BasicTopologyEntity.

Definition at line 61 of file RefEdge.hpp.

{ return DagType::ref_edge_type(); }
int RefEdge::dimension ( ) const [virtual]

Returns the geometric dimension of the entity. vertex == 0, edge == 1, etc.

Reimplemented from RefEntity.

Definition at line 763 of file RefEdge.cpp.

{
  return 1;
}

Definition at line 271 of file RefEdge.cpp.

{
  return this->end_vertex()->coordinates();
}
double RefEdge::end_param ( )

Definition at line 506 of file RefEdge.cpp.

{
    // Get the Curve associated with this RefEdge
  Curve* curve_ptr = get_curve_ptr();
  
  assert(curve_ptr != NULL);
  
  if (curve_ptr->bridge_sense() == CUBIT_REVERSED)
    return -(curve_ptr->start_param());
  else
    return curve_ptr->end_param();
}

Definition at line 1362 of file RefEdge.cpp.

{
    // Get the first (and only) Chain associated with this RefEdge.
  Chain* chain_ptr = get_chain_ptr();
  
    // Ask the Chain for its end (last) RefVertex
  if (chain_ptr == NULL)
  {
    return NULL;
  }
  
  else
  {
    return chain_ptr->end_vertex();
  }
}
const std::type_info& RefEdge::entity_type_info ( ) const [inline, virtual]

Returns the type info of this enity.

Implements RefEntity.

Definition at line 62 of file RefEdge.hpp.

{ return typeid(RefEdge); }
CubitStatus RefEdge::evaluate_exterior_angle ( double *  exterior_angle)

Definition at line 1475 of file RefEdge.cpp.

{
  //- Get the center point of this curve
  CubitVector center = this->center_point();
  CubitVector surf_norm[2];
  CubitVector tangent;
  RefFace* ref_face = NULL;
  RefEdge* next_curve = NULL;

  int i, j, k;

    //- Get surfaces on this curve
  DLIList<RefFace*> surf_list;
  this->ref_faces(surf_list);
  if(surf_list.size() == 2)
  {
    for( i = 0; i < surf_list.size(); i++)
    {
      RefFace* this_surf = surf_list.get_and_step();
      surf_norm[i] = this_surf->normal_at(center);
      if( i == 1)
      {
          //- Get the referance surface and curve to get exterior angle 
        ref_face = this_surf;
        DLIList<DLIList<RefEdge*> > loop_lists;
        ref_face->ref_edge_loops(loop_lists);
        for( j = 0; j < loop_lists.size(); j++ )
        {
          DLIList<RefEdge*>& current_list = loop_lists[j];
          for( k = 0; k < current_list.size(); k++ )
          {
            RefEdge* current_curve = current_list.get_and_step();
            if( current_curve == this )
            {
              next_curve = current_list.get();
              break;
            }
          }
          if(next_curve != NULL)
              break;
        }
      }
    }
  }
  else
  {
    *exterior_angle = 0.0;
    PRINT_ERROR("There aren't 2 surfaces on curve %d.  Can't compute angle.\n", this->id());
    return CUBIT_FAILURE;
  }

  this->tangent(center, tangent, next_curve, ref_face);

    //Find the angle from normal 1 to normal 0.
  double rad_angle = tangent.vector_angle_quick(surf_norm[1], surf_norm[0]); //angle in radians
  double angle = rad_angle * (180.0 / 3.14159); //angle in degrees

    //Now convert to the exterior angle between the two surfaces.
  angle += 180.0;
  if( angle > 360.0 )
      angle -= 360.0;
    
    //- Return the exterior angle for this curve in degrees.
  *exterior_angle = angle;
  return CUBIT_SUCCESS;
}
double RefEdge::fraction_from_arc_length ( RefVertex root_vertex,
double  length 
)

Definition at line 577 of file RefEdge.cpp.

{
  if (root_vertex != start_vertex() && root_vertex != end_vertex())
    return -1.0;
 
  if (geometry_type() == POINT_CURVE_TYPE || get_arc_length() < GEOMETRY_RESABS)
    return 0.0;
 
  if (length >= get_arc_length())
    return 1.0;

  if (root_vertex == start_vertex())
    return length/get_arc_length();

  else
    return 1-length/get_arc_length();
}

Definition at line 367 of file RefEdge.cpp.

{
    // Get the Curve associated with this RefEdge
  Curve* curve_ptr = get_curve_ptr();
  
  assert(curve_ptr != NULL);
  
    // Get the length
  return curve_ptr->get_arc_length();
}
double RefEdge::get_arc_length ( const CubitVector point1,
const CubitVector point2 
)

Definition at line 378 of file RefEdge.cpp.

{
    // Get the Curve associated with this RefEdge
  Curve* curve_ptr = get_curve_ptr();
  
  assert(curve_ptr != NULL);
  
    // Get the length between the 2 points
  return curve_ptr->get_arc_length(point1, point2);
}
double RefEdge::get_arc_length ( const CubitVector point1,
int  whichEnd 
)

Definition at line 390 of file RefEdge.cpp.

{
  Curve* curve_ptr = get_curve_ptr();
  
  assert (curve_ptr != NULL);
  
  if (curve_ptr->bridge_sense() == CUBIT_REVERSED)
    which_end = 1 - which_end;
  
  return curve_ptr->get_arc_length(point1, which_end);
}
CubitStatus RefEdge::get_center_radius ( CubitVector center,
double &  radius 
) [virtual]

Definition at line 161 of file RefEdge.cpp.

{
   Curve* curve_ptr = get_curve_ptr();

   if( curve_ptr != NULL )
   {

      if( curve_ptr->geometry_type() != ARC_CURVE_TYPE &&
          curve_ptr->geometry_type() != ELLIPSE_CURVE_TYPE )
         return CUBIT_FAILURE;

      if( curve_ptr->get_center_radius( center, radius ) == CUBIT_FAILURE )
         return CUBIT_FAILURE;
   }
   else 
   {
      PRINT_WARNING("In RefEdge::get_center_radius\n"
                    "         %s (curve %d) is not associated with a valid\n"
                    "         underlying geoemtric Curve\n",
                    entity_name().c_str(), id()) ;
      return CUBIT_FAILURE ;
   }
   
   return CUBIT_SUCCESS;
}

Definition at line 116 of file RefEdge.cpp.

{
  GroupingEntity* gpe_ptr = get_first_grouping_entity_ptr();
  if ( !gpe_ptr ||        // no chain
       gpe_ptr->next() )  // multiple chains
  {
    return 0;
  }
  
  return static_cast<Chain*>(gpe_ptr);
}

Definition at line 403 of file RefEdge.cpp.

{
  CubitVector start_pos(start_vertex()->coordinates());
  CubitVector end_pos(start_vertex()->coordinates());
  double distance = (start_pos - end_pos).length();
  return distance;
}
static const char* RefEdge::get_class_name ( ) [inline, static]

Reimplemented from BasicTopologyEntity.

Definition at line 51 of file RefEdge.hpp.

     {
       return "Curve";
     }
CubitStatus RefEdge::get_co_edges ( DLIList< CoEdge * > &  co_edges_found_list,
RefFace input_ref_face_ptr = NULL 
)

Definition at line 670 of file RefEdge.cpp.

{
  for (SenseEntity* coedge_ptr = get_first_sense_entity_ptr();
       coedge_ptr;
       coedge_ptr = coedge_ptr->next_on_bte())
  {
    if (!input_ref_face_ptr ||
      coedge_ptr->get_parent_basic_topology_entity_ptr() == input_ref_face_ptr)
    {
      co_edges_found_list.append (dynamic_cast<CoEdge*>(coedge_ptr));
    }
  }
  
  return CUBIT_SUCCESS;
}

Definition at line 97 of file RefEdge.cpp.

const Curve * RefEdge::get_curve_ptr ( ) const

Definition at line 102 of file RefEdge.cpp.

CubitStatus RefEdge::get_graphics ( GMem polyline,
double  angle_tolerance = 0,
double  distance_tolerance = 0,
double  max_edge_length = 0 
)

Definition at line 1440 of file RefEdge.cpp.

{
  Curve* curve_ptr = get_curve_ptr();
  if (!curve_ptr)
  {
    PRINT_ERROR("RefEdge %d is invalid -- no attached Curve.\n",id());
    return CUBIT_FAILURE;
  }
  
  return curve_ptr->get_geometry_query_engine()->get_graphics(curve_ptr,
    &polyline, angle_tolerance, distance_tolerance, max_edge_length );
}
CubitStatus RefEdge::get_interior_extrema ( DLIList< CubitVector * > &  interior_points,
CubitSense return_sense 
) const

Definition at line 212 of file RefEdge.cpp.

{
  CubitStatus result;
  
    // Cast away the constness of the Curve
  Curve* curve = const_cast<Curve*>(get_curve_ptr());
  result = curve->get_interior_extrema(interior_points, return_sense);
  if (curve->bridge_sense() == CUBIT_REVERSED)
    return_sense = CubitUtil::opposite_sense(return_sense);
  return result;
}
double RefEdge::get_lower_param_bound ( ) [private]
int RefEdge::get_mark ( ) [inline]

Definition at line 582 of file RefEdge.hpp.

{ return markedFlag; }
RefEdge * RefEdge::get_other_curve ( RefVertex common_vertex,
RefFace ref_face_ptr 
)

Definition at line 940 of file RefEdge.cpp.

{
  DLIList<RefEdge*> curves;
  ref_face_ptr->ref_edges(curves);
  for(int ii = curves.size(); ii>0; ii--)
  {
    RefEdge* temp_edge = curves.get_and_step();
    if((temp_edge->is_directly_related(common_vertex)) &&
       (this != temp_edge))
       return temp_edge;
  }
  return NULL;
}
CubitBoolean RefEdge::get_param_range ( double &  start_param,
double &  end_param 
)

Definition at line 519 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
    // Now get the parameter values of the start and end locations
    // of this RefEdge
  CubitBoolean result = curvePtr->get_param_range( start_param, end_param );
  
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
  {
    double tmp_start_param = start_param;
    start_param = -end_param;
    end_param = -tmp_start_param;
  }
  
  return result;
}
void RefEdge::get_parent_ref_entities ( DLIList< RefEntity * > &  entity_list) [virtual]

Appends all RefEntities that own this (parent RefEntities) to entity_list. (The query goes up just one dimension. For example, if this is a vertex, the resulting list contains only RefEdges).

Implements RefEntity.

Definition at line 1542 of file RefEdge.cpp.

{

  // First get the type of RefEntity that is a child of "this" one
  DagType parent_type = get_parent_ref_entity_type();;

  DLIList<TopologyEntity*> tempList ;

  CubitStatus result = ModelQueryEngine::instance()->
      query_model( *this, parent_type, tempList );
  if (result == CUBIT_FAILURE)
  {
    PRINT_ERROR("In RefEntity::get_parent_ref_entities\n");
    PRINT_ERROR("       Query failed for unknown reason.\n");
    return;
  }

  entity_list.clean_out();
  for(int i=0; i<tempList.size(); i++)
  {
    entity_list.append(static_cast<ParentType*>(tempList[i]));
  }
}
CubitStatus RefEdge::get_point_direction ( CubitVector origin,
CubitVector direction 
) [virtual]

Definition at line 128 of file RefEdge.cpp.

{
   Curve* curve_ptr = get_curve_ptr();

   if( curve_ptr != NULL )
   {

      if( curve_ptr->geometry_type() != STRAIGHT_CURVE_TYPE ) 
         return CUBIT_FAILURE;

      if( curve_ptr->get_point_direction( origin, direction ) == CUBIT_FAILURE )
      {
         origin = start_coordinates();
         direction = end_coordinates() - origin;
         direction.normalize();
      }
   }
   else 
   {
      PRINT_WARNING("In RefEdge::get_point_direction\n"
                    "         %s (curve %d) is not associated with a valid\n"
                    "         underlying geometric Curve\n",
                    entity_name().c_str(), id()) ;
      return CUBIT_FAILURE ;
   }
   
   if (curve_ptr->bridge_sense() == CUBIT_REVERSED)
     direction = -direction;
   
   return CUBIT_SUCCESS;

}
CubitStatus RefEdge::get_two_co_edges ( RefEdge next_ref_edge,
RefFace ref_face_ptr,
CoEdge *&  co_edge_this,
CoEdge *&  co_edge_next 
)

Definition at line 968 of file RefEdge.cpp.

{
  DLIList<Loop*> loop_list;
  Loop *loop_ptr;
  CubitStatus status;
  DLIList<CoEdge*> co_edge_list;
  
    //First get the loops for this ref_face;
  ref_face_ptr->loops( loop_list );
    //Now we want to find the coedge list that
    //has 'this' refedge followed by the next one
  assert( loop_list.size() != 0 );
  
  for ( int i = loop_list.size(); i--; )
  {
      // get the ordered co-edges of this loop
    loop_ptr = loop_list.get_and_step();
    co_edge_list.clean_out();   
    status = loop_ptr->ordered_co_edges( co_edge_list );

      //Now find the coedges corresponding to this ref_edge
      //  and the next ref_edge.
    if ( status == CUBIT_SUCCESS ) {
      for ( int j = co_edge_list.size(); j--; ) 
      {
          // candidates
        co_edge_this = co_edge_list.get_and_step();
        co_edge_next = co_edge_list.get();
        
          // really correspond to this and next edge?
        if ( co_edge_this->get_ref_edge_ptr() == this &&
             co_edge_next->get_ref_edge_ptr() == next_ref_edge ) 
        {
          return CUBIT_SUCCESS;          
        }        
      }
    }
  }  
  co_edge_this = co_edge_next = NULL;
  PRINT_ERROR("in RefEdge::get_two_co_edges.\n"
              "Couldn't find CoEdgeList with this edge\n"
              "and the next edge passed in.\n");
  return CUBIT_FAILURE;
}
double RefEdge::get_upper_param_bound ( ) [private]
void RefEdge::initialize ( ) [private]

Definition at line 1271 of file RefEdge.cpp.

{
    // Initialize some member data
  refEdgeClone = 0;
  markedFlag = CUBIT_FALSE;

    // Make sure the arc length is not zero if there are start and end 
    // RefVertex'es already assigned to this RefEdge
  if ( get_arc_length() < CUBIT_DBL_MIN )  
  {
    if ( start_vertex() && end_vertex() )
    {
      CubitVector start_pt = start_vertex()->coordinates();
      CubitVector end_pt   = end_vertex()->coordinates();
      PRINT_WARNING ( "(RefEdge::initialize): Edge has zero arclength.\n"
                      "  Start vertex location is (%9.2f, %9.2f, %9.2f ).\n"
                      "  End   vertex location is (%9.2f, %9.2f, %9.2f ).\n",
          start_pt.x(), start_pt.y(), start_pt.z(),
          end_pt.x(), end_pt.y(), end_pt.z() );
    }
    else if (!mSuppressEdgeLengthWarning)
    {
      PRINT_WARNING( "Edge found with zero arclength\n"
                     "  For cones, this may be normal.\n");
    }
    
  }
  
    // Set the Entity ID for this new RefEdge
   GeometryEntity* geom_ptr = get_geometry_entity_ptr();
   int saved_id = geom_ptr->get_saved_id();
   if ( !saved_id || RefEntityFactory::instance()->get_ref_edge(saved_id) )
   {
     saved_id =  RefEntityFactory::instance()->next_ref_edge_id();
     geom_ptr->set_saved_id(saved_id);
   }
   entityId = saved_id;
  
     // read and initialize attributes
   auto_read_cubit_attrib();
   auto_actuate_cubit_attrib();

     // Assign a default entity name
   assign_default_name();
}

Definition at line 640 of file RefEdge.cpp.

{
  double temp_val;
  return this->is_periodic(temp_val);
}
CubitBoolean RefEdge::is_periodic ( double &  period)

Definition at line 646 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
    // Now get the parameter values of the start and end locations
    // of this RefEdge
  CubitBoolean periodic = curvePtr->is_periodic(period);
  
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
    period = -period;
  
  return periodic;
}

Definition at line 1415 of file RefEdge.cpp.

{
   Curve* curve_ptr = get_curve_ptr();
   if (curve_ptr == NULL) {
      PRINT_WARNING("\tWARNING: Null underlying curve for %s, (%s %d)\n",
         entity_name().c_str(), class_name(), id());
      return CUBIT_FALSE;
   }

   return curve_ptr->is_tolerant();
}
double RefEdge::length_from_u ( double  parameter1,
double  parameter2 
)

Definition at line 628 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
    return -(curvePtr->length_from_u(-parameter1, -parameter2));
  else
    return curvePtr->length_from_u(parameter1, parameter2);
}
double RefEdge::measure ( ) [virtual]

A generic geometric extent function. Returns volume for RefVolumes, area for RefFaces, length for RefEdge, and 1.0 for RefVertices A RefGroup calculates the maximum dimension of its contained entities and returns the sum of the measures() of all entities of that dimension. Default return value is 0.0 for all other entities.

Reimplemented from BasicTopologyEntity.

Definition at line 768 of file RefEdge.cpp.

{
  return get_arc_length();
}

Returns the type of measure: (volume, area, length, or N/A)

Reimplemented from RefEntity.

Definition at line 773 of file RefEdge.cpp.

{
  return "length";
}
CubitStatus RefEdge::mid_point ( const CubitVector point1,
const CubitVector point2,
CubitVector midPoint 
)

Definition at line 432 of file RefEdge.cpp.

{
    // Get the Curve associated with this RefEdge
  Curve* curve_ptr = get_curve_ptr();
  
  assert(curve_ptr != NULL);
  
    // Get the global location of parameter3
  return curve_ptr->mid_point(point1, point2, mid_point);
  
}

Definition at line 446 of file RefEdge.cpp.

{
    // Get the Curve associated with this RefEdge
  Curve* curve_ptr = get_curve_ptr();
  
  assert(curve_ptr != NULL);
  
    // Get the global location of parameter3
  return curve_ptr->mid_point(mid_point);
  
}
void RefEdge::move_to_curve ( CubitVector vector) [virtual]

Definition at line 198 of file RefEdge.cpp.

{
    // Get the Curve associated with this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
    // Move the point to the Curve (the following call modifes the values
    // of the input "vector", if necessary)
  CubitVector closest_point;
  curvePtr->closest_point_trimmed(vector, closest_point);
  vector.set( closest_point.x(),
              closest_point.y(),
              closest_point.z() );
}

Definition at line 778 of file RefEdge.cpp.

{
  DLIList<RefFace*> ref_faces, ref_faces_for_other_edge;
  this->ref_faces(ref_faces);
  other_edge->ref_faces(ref_faces_for_other_edge);
                                                                                
  int i, j;
  int num = 0;
  for (i = 0; i < ref_faces.size(); i++)
  {
     RefFace * ref_face = ref_faces.get_and_step();
     for (j = 0; j < ref_faces_for_other_edge.size(); j++)
        if (ref_face == ref_faces_for_other_edge.get_and_step())
           num++;
  }
                                                                                
  return num;
}
void RefEdge::operator= ( const RefEdge ) [private]
RefFace * RefEdge::other_face ( RefFace not_face,
RefVolume ref_volume = NULL 
)

Definition at line 1033 of file RefEdge.cpp.

{
    //- return the (an) other face sharing this edge, which also borders
    //- ref_volume if non-NULL
  DLIList<RefFace*> temp_faces;
  ref_faces(temp_faces);
  int i;
  for (i = temp_faces.size(); i > 0; i--) {
    RefFace *other_face = temp_faces.get_and_step();
    if (other_face != not_face &&
        (!ref_volume || other_face->is_directly_related(ref_volume)))
      return other_face;
  }
  
  return NULL;
}
RefVertex * RefEdge::other_vertex ( RefVertex refVertexPtr)

Definition at line 1017 of file RefEdge.cpp.

{
  if ( vertex == start_vertex() )
  { 
    return end_vertex();
  }
  else if ( vertex == end_vertex() )
  {
    return start_vertex();
  }
  else
  {
    return NULL;
  }
}

Definition at line 246 of file RefEdge.cpp.

{
  Curve *curve_ptr = get_curve_ptr();
  return curve_ptr->point_containment(point);
}
CubitStatus RefEdge::point_from_arc_length ( double  root_param,
double  arc_length,
CubitVector new_point 
)

Definition at line 612 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
    arc_length = -arc_length;
    
    // Now get the parameter values of the start and end locations
    // of this RefEdge
  return curvePtr->point_from_arc_length (root_param, arc_length,
                                          new_point );
}
CubitStatus RefEdge::point_from_arc_length ( const CubitVector root_point,
double  arc_length,
CubitVector new_point 
)

Definition at line 596 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
    arc_length = -arc_length;
    
    // Now get the parameter values of the start and end locations
    // of this RefEdge
  return curvePtr->point_from_arc_length (root_point, arc_length,
                                          new_point );
}
CubitStatus RefEdge::position_from_fraction ( double  fraction_along_curve,
CubitVector ouput_position 
)

Definition at line 457 of file RefEdge.cpp.

{
  
    //Get the Curve of this RefEdge.
  Curve* curve_ptr = this->get_curve_ptr();
  
  assert( fraction_along_curve < 1.0000001 && fraction_along_curve > -0.0000001 );
    //Now get the postion from this fraction value.
  
  if (curve_ptr->bridge_sense() == CUBIT_REVERSED)
    fraction_along_curve = 1.0 - fraction_along_curve;  
    
  CubitStatus result = curve_ptr->position_from_fraction( fraction_along_curve,
                                                          output_position );
  return result;
}
CubitStatus RefEdge::position_from_u ( double  u_value,
CubitVector output_position 
)

Definition at line 552 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
    u_value = -u_value;
  
    // Now get the parameter values of the start and end locations
    // of this RefEdge
  return curvePtr->position_from_u(u_value, output_position);
}
CubitStatus RefEdge::relative_sense ( RefEdge ref_edge_ptr_2,
double  tolerance_factor,
CubitSense sense,
CubitBoolean spatially_equal,
CubitBoolean  force_merge = CUBIT_FALSE 
)

Definition at line 1051 of file RefEdge.cpp.

{
    // It is assumed that the endpoints have been determined to be
    // spatially equal before this function is called.
  
    // Algorithm: Get a point from 'this' (the 1/3 point).  See if it
    // lies on the second RefEdge also, within tolerance.  If so,
    // spatially_equal is set to true.  Next, compare the normal directions
    // at the common point.  If they point the same direction (dot product
    // of tangents is > 0), sense is forward.  Else, sense is reversed.  If
    // the tangent dot product is between -.5 and .5, a warning is printed.
  
  const double ONE_THIRD = 1.0/3.0;
  if (sense)
    *sense = CUBIT_FORWARD;
  CubitStatus result;
  
    // Find the point 1/3 along *this*
  CubitVector test_point_1, test_point_2;
  result = this->position_from_fraction( ONE_THIRD,
                                         test_point_1);
  if ( result != CUBIT_SUCCESS )
  {
    PRINT_ERROR("Error in RefEdge::compare(refedges).\n"
                "Can't find position 1/3 along curve.\n");
    return CUBIT_FAILURE;
  }
  
    // See if the 1/3 point on *this* lies on the other curve
  if ( ref_edge_ptr_2->closest_point_trimmed(test_point_1, test_point_2)
       != CUBIT_SUCCESS )
  {
    return CUBIT_FAILURE;
  }
  if ( GeometryQueryTool::instance()->
       about_spatially_equal(test_point_1, test_point_2,tolerance_factor ))
  {
    spatially_equal = CUBIT_TRUE;
  }
  else if ( !force_merge ) 
  {
    return CUBIT_FAILURE;
  }
  
    // Now find the sense
  if (sense)
  {
    CubitVector tangent_1, tangent_2;
    result = this->closest_point(test_point_2,
                                 test_point_1,
                                 &tangent_1);
    if (result == CUBIT_SUCCESS)
      result = ref_edge_ptr_2->closest_point(test_point_1,
                                             test_point_2,
                                             &tangent_2);
    if ( result != CUBIT_SUCCESS )
    {
      PRINT_ERROR("Error in RefEdge::relative_sense(refedges).\n"
                  "Can't find Curve tangents.\n");
      return CUBIT_FAILURE;
    }
      
      // Find the sense
            
    //If one of the curves is zero-length, it will have a zero
    //tangent vector.
    double len_product = tangent_1.length() * tangent_2.length();
    if( len_product > CUBIT_DBL_MIN )
    {
        
      double dot_product = (tangent_1 % tangent_2) / len_product;
      if (dot_product < 0)
      *sense = CUBIT_REVERSED;
        //if (dot_product > -.5 && dot_product < .5)
      //      tangent_warning = CUBIT_TRUE;
    }
    else
    {
      //If one of the tangents is zero-length, one of the curves had
      //better be as well.
      assert( (measure() * ref_edge_ptr_2->measure()) < CUBIT_RESABS );
    }
  }
  return CUBIT_SUCCESS;
}
void RefEdge::reverse_topology ( ) [virtual]

Reimplemented from TopologyEntity.

Definition at line 1379 of file RefEdge.cpp.

{
  Chain* chain_ptr = get_chain_ptr();
  chain_ptr->reverse_direction();

    // switch co_edge senses
  DLIList<SenseEntity*> co_edge_list;
  get_sense_entity_list(co_edge_list);
  for ( int i = co_edge_list.size(); i--; ) 
    co_edge_list.get_and_step()->reverse_sense();
}

Definition at line 1391 of file RefEdge.cpp.

{
  DLIList<CoEdge*> co_edge_list;
  get_co_edges(co_edge_list, face);
  CoEdge *coedge;
  CubitSense my_sense = CUBIT_UNKNOWN;
  for ( int i = co_edge_list.size(); i--; ) {
    coedge = co_edge_list.get_and_step();
    if(coedge->get_sense() == CUBIT_FORWARD)
      {
    if (my_sense == CUBIT_REVERSED)
      return CUBIT_UNKNOWN;
    my_sense = CUBIT_FORWARD;
      }
    else
      {
    if (my_sense == CUBIT_FORWARD)
      return CUBIT_UNKNOWN;
    my_sense = CUBIT_REVERSED;
      }
  }
  return my_sense;
}
void RefEdge::set_mark ( int  set_value) [inline]

Definition at line 586 of file RefEdge.hpp.

{ markedFlag = set_value; }

Definition at line 266 of file RefEdge.cpp.

{
  return this->start_vertex()->coordinates();
}

Definition at line 484 of file RefEdge.cpp.

{
    // Get the Curve associated with this RefEdge
  Curve* curve_ptr = get_curve_ptr();
  
  assert(curve_ptr != NULL);
  
  if (curve_ptr->bridge_sense() == CUBIT_REVERSED)
    return -(curve_ptr->end_param());
  else
    return curve_ptr->start_param();
}

Definition at line 1331 of file RefEdge.cpp.

{
    // Get the first (and only) Chain associated with this RefEdge.
  Chain* chain_ptr = this->get_chain_ptr();
  
    // Ask the Chain for its first RefVertex
  if (chain_ptr == NULL)
  {
    return NULL;
  }
  
  else
  {
    return chain_ptr->start_vertex();
  }
}
void RefEdge::suppress_edge_length_warning ( bool  flag) [static]

Definition at line 1470 of file RefEdge.cpp.

void RefEdge::tangent ( const CubitVector point,
CubitVector tangent_vec 
)

Definition at line 276 of file RefEdge.cpp.

{
    // Get the tangent at the point closest to "point" on the RefEdge.
    // The tangent is always pointing in the positive direction of the
    // RefEdge
  Curve* curve_ptr = this->get_curve_ptr();
  CubitVector closest_point;
  curve_ptr->closest_point( point, closest_point, &tangent_vec );
  if (curve_ptr->bridge_sense() == CUBIT_REVERSED)
    tangent_vec = -tangent_vec;
}
CubitStatus RefEdge::tangent ( const CubitVector point,
CubitVector tangent_vec,
RefFace ref_face_ptr 
)

Definition at line 300 of file RefEdge.cpp.

{
    //Get the tangent for this edge.
  tangent( point, tangent_vec );
  
    //Now allign the tangent with the face.
    //For this function we assume that there must only be one
    //co edge for this ref_edge.
  DLIList<CoEdge*> co_edge_list;
  get_co_edges( co_edge_list, ref_face_ptr );
 
  if( co_edge_list.size() != 1 )
  {
      PRINT_ERROR( "Ambiguous case in tangent calculation.\n" );
      return CUBIT_FAILURE;
  }
//  assert( co_edge_list.size() == 1 ); 
  
  if ( co_edge_list.get()->get_sense() == CUBIT_REVERSED )
      tangent_vec = -tangent_vec;
  return CUBIT_SUCCESS;
}
CubitStatus RefEdge::tangent ( const CubitVector point,
CubitVector tangent_vec,
RefEdge next_ref_edge,
RefFace ref_face_ptr 
)

Definition at line 337 of file RefEdge.cpp.

{
  CoEdge *co_edge_this = NULL;
  CoEdge *co_edge_next = NULL;
    //First get the two coedges that corrispond to
    //this ref_edge an the next one, with reference to
    //the ref_face_ptr.
  CubitStatus status = get_two_co_edges( next_ref_edge,
                                         ref_face_ptr,
                                         co_edge_this,
                                         co_edge_next );
  if (status == CUBIT_FAILURE )
      return status;
  
  assert(co_edge_this != NULL );
  assert(co_edge_next != NULL );
  
    //Now get the tangent from this curve.
  tangent ( point, tangent_vec );
  
    //with the go_edge data we have, we can get the tangent
    //going in the right direction...
  if ( co_edge_this->get_sense() == CUBIT_REVERSED )
      tangent_vec = -tangent_vec;
  return CUBIT_SUCCESS;
}
double RefEdge::u_from_arc_length ( double  root_param,
double  arc_length 
)

Definition at line 566 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
    return -(curvePtr->u_from_arc_length(-root_param, -arc_length));
  else
    return curvePtr->u_from_arc_length(root_param, arc_length);
}
double RefEdge::u_from_position ( const CubitVector input_position)

Definition at line 538 of file RefEdge.cpp.

{
    // Get the Curve of this RefEdge
  Curve* curvePtr = this->get_curve_ptr();
  
    // Now get the parameter values of the start and end locations
    // of this RefEdge
  double param = curvePtr->u_from_position(input_position);
  if (curvePtr->bridge_sense() == CUBIT_REVERSED)
    param = -param;
  
  return param;
}
int RefEdge::validate ( ) [virtual]

Perform checks to see if entity valid.

Reimplemented from RefEntity.

Definition at line 1230 of file RefEdge.cpp.

{
    //- This function determines whether the entity is valid.
    //- Several types of checks can be done, 
  int error = 0;
  
    // Perform general RefEntity checks (measure > 0)
  error += RefEntity::validate();
  
    // Pass through to curve and add in its validation
  Curve *curve = get_curve_ptr();
  
    // check curve ptr
  if (curve != NULL) {
      // Check underlying curve
    DLIList <TopologyEntity*> bad_entities;
    error += curve->validate(entity_name(),bad_entities);
  } else {
    PRINT_WARNING("\tWARNING: Null underlying curve for %s, (%s %d)\n",
                  entity_name().c_str(), class_name(), id());
    error++;
  }
  return error;
}

Friends And Related Function Documentation

friend class RefEntityFactory [friend]

Definition at line 45 of file RefEdge.hpp.


Member Data Documentation

Reimplemented from RefEntity.

Definition at line 569 of file RefEdge.hpp.

bool RefEdge::mSuppressEdgeLengthWarning = false [static, private]

Definition at line 573 of file RefEdge.hpp.

int RefEdge::refEdgeClone [private]

Definition at line 567 of file RefEdge.hpp.


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