1 // Copyright (C) 2007-2010 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // SMESH SMESH : implementaion of SMESH idl descriptions
24 // File : StdMeshers_FaceSide.hxx
25 // Created : Wed Jan 31 18:41:25 2007
26 // Author : Edward AGAPOV (eap)
29 #include "StdMeshers_FaceSide.hxx"
31 #include "SMDS_EdgePosition.hxx"
32 #include "SMDS_MeshNode.hxx"
33 #include "SMESHDS_Mesh.hxx"
34 #include "SMESHDS_SubMesh.hxx"
35 #include "SMESH_Algo.hxx"
36 #include "SMESH_Mesh.hxx"
37 #include "SMESH_MesherHelper.hxx"
38 #include "SMESH_ComputeError.hxx"
39 #include "SMESH_Block.hxx"
41 #include <Adaptor2d_Curve2d.hxx>
42 #include <BRepAdaptor_CompCurve.hxx>
43 #include <BRepAdaptor_Curve.hxx>
44 #include <BRep_Builder.hxx>
45 #include <BRep_Tool.hxx>
47 #include <TopExp_Explorer.hxx>
49 #include <TopoDS_Face.hxx>
50 #include <TopoDS_Vertex.hxx>
51 #include <TopoDS_Wire.hxx>
53 #include <GCPnts_AbscissaPoint.hxx>
54 #include <Geom2dAdaptor_Curve.hxx>
58 #include "utilities.h"
60 //================================================================================
62 * \brief Constructor of a side of one edge
63 * \param theFace - the face
64 * \param theEdge - the edge
66 //================================================================================
68 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
69 const TopoDS_Edge& theEdge,
71 const bool theIsForward,
72 const bool theIgnoreMediumNodes)
74 list<TopoDS_Edge> edges(1,theEdge);
75 *this = StdMeshers_FaceSide( theFace, edges, theMesh, theIsForward, theIgnoreMediumNodes );
78 //================================================================================
80 * \brief Constructor of a side of several edges
81 * \param theFace - the face
82 * \param theEdge - the edge
84 //================================================================================
86 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
87 list<TopoDS_Edge>& theEdges,
89 const bool theIsForward,
90 const bool theIgnoreMediumNodes)
92 int nbEdges = theEdges.size();
93 myEdge.resize( nbEdges );
94 myC2d.resize( nbEdges );
95 myC3dAdaptor.resize( nbEdges );
96 myFirst.resize( nbEdges );
97 myLast.resize( nbEdges );
98 myNormPar.resize( nbEdges );
99 myEdgeLength.resize( nbEdges );
100 myIsUniform.resize( nbEdges, true );
102 myNbPonits = myNbSegments = 0;
104 myMissingVertexNodes = false;
105 myIgnoreMediumNodes = theIgnoreMediumNodes;
106 myDefaultPnt2d = gp_Pnt2d( 1e+100, 1e+100 );
107 if ( nbEdges == 0 ) return;
109 SMESHDS_Mesh* meshDS = theMesh->GetMeshDS();
112 list<TopoDS_Edge>::iterator edge = theEdges.begin();
113 TopoDS_Iterator vExp;
114 for ( int index = 0; edge != theEdges.end(); ++index, ++edge )
116 int i = theIsForward ? index : nbEdges - index - 1;
117 myEdgeLength[i] = SMESH_Algo::EdgeLength( *edge );
118 if ( myEdgeLength[i] < DBL_MIN ) nbDegen++;
119 myLength += myEdgeLength[i];
121 if ( !theIsForward ) myEdge[i].Reverse();
123 if ( theFace.IsNull() )
124 BRep_Tool::Range( *edge, myFirst[i], myLast[i] );
126 myC2d[i] = BRep_Tool::CurveOnSurface( *edge, theFace, myFirst[i], myLast[i] );
127 if ( myEdge[i].Orientation() == TopAbs_REVERSED )
128 std::swap( myFirst[i], myLast[i] );
130 if ( SMESHDS_SubMesh* sm = meshDS->MeshElements( *edge )) {
131 int nbN = sm->NbNodes();
132 if ( theIgnoreMediumNodes ) {
133 SMDS_ElemIteratorPtr elemIt = sm->GetElements();
134 if ( elemIt->more() && elemIt->next()->IsQuadratic() )
135 nbN -= sm->NbElements();
138 myNbSegments += sm->NbElements();
140 // TopExp::FirstVertex() and TopExp::LastVertex() return NULL from INTERNAL edge
141 vExp.Initialize( *edge );
142 if ( vExp.Value().Orientation() == TopAbs_REVERSED ) vExp.Next();
143 if ( SMESH_Algo::VertexNode( TopoDS::Vertex( vExp.Value()), meshDS ))
144 myNbPonits += 1; // for the first end
146 myMissingVertexNodes = true;
148 // check if edge has non-uniform parametrization (issue 0020705)
149 if ( !myC2d[i].IsNull() && myEdgeLength[i] > DBL_MIN)
151 Geom2dAdaptor_Curve A2dC( myC2d[i] );
152 double p2 = myFirst[i]+(myLast[i]-myFirst[i])/2., p4 = myFirst[i]+(myLast[i]-myFirst[i])/4.;
153 double d2 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p2 );
154 double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
155 //cout<<"len = "<<len<<" d2 = "<<d2<<" fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
156 myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
157 if ( !myIsUniform[i] )
161 Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],L,fp,lp);
162 myC3dAdaptor[i].Load( C3d, fp,lp );
166 vExp.Initialize( theEdges.back() );
167 if ( vExp.Value().Orientation() != TopAbs_REVERSED ) vExp.Next();
170 if ( SMESH_Algo::VertexNode( TopoDS::Vertex( vExp.Value()), meshDS ))
171 myNbPonits++; // for the last end
173 myMissingVertexNodes = true;
175 if ( nbEdges > 1 && myLength > DBL_MIN ) {
176 const double degenNormLen = 1.e-5;
177 double totLength = myLength;
179 totLength += myLength * degenNormLen * nbDegen;
180 double prevNormPar = 0;
181 for ( int i = 0; i < nbEdges; ++i ) {
182 if ( myEdgeLength[ i ] < DBL_MIN )
183 myEdgeLength[ i ] = myLength * degenNormLen;
184 myNormPar[ i ] = prevNormPar + myEdgeLength[i]/totLength;
185 prevNormPar = myNormPar[ i ];
188 myNormPar[nbEdges-1] = 1.;
192 //================================================================================
194 * \brief Constructor of a side for vertex using data from other FaceSide
195 * \param theVertex - the vertex
196 * \param theSide - the side
198 //================================================================================
200 StdMeshers_FaceSide::StdMeshers_FaceSide(const SMDS_MeshNode* theNode,
201 const gp_Pnt2d thePnt2d,
202 const StdMeshers_FaceSide* theSide)
206 myMesh = theSide->GetMesh();
207 myDefaultPnt2d = thePnt2d;
209 myPoints = theSide->GetUVPtStruct();
210 myNbPonits = myNbSegments = myPoints.size();
211 std::vector<uvPtStruct>::iterator it = myPoints.begin();
212 for(; it!=myPoints.end(); it++) {
213 (*it).u = thePnt2d.X();
214 (*it).v = thePnt2d.Y();
216 (*it).node = theNode;
220 //================================================================================
222 * \brief Return info on nodes on the side
223 * \retval UVPtStruct* - array of data structures
225 //================================================================================
227 const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool isXConst,
228 double constValue) const
230 if ( myPoints.empty() ) {
232 if ( NbEdges() == 0 ) return myPoints;
234 SMESHDS_Mesh* meshDS = myMesh->GetMeshDS();
235 SMESH_MesherHelper helper(*myMesh);
238 // sort nodes of all edges putting them into a map
240 map< double, const SMDS_MeshNode*> u2node;
242 for ( int i = 0; i < myEdge.size(); ++i )
244 // Put 1st vertex node of a current edge
245 TopoDS_Vertex VV[2]; // TopExp::FirstVertex() returns NULL for INTERNAL edge
246 for ( TopoDS_Iterator vIt(myEdge[i]); vIt.More(); vIt.Next() )
247 VV[ VV[0].IsNull() ? 0 : 1 ] = TopoDS::Vertex(vIt.Value());
248 if ( VV[0].Orientation() == TopAbs_REVERSED ) std::swap ( VV[0], VV[1] );
249 const SMDS_MeshNode* node = SMESH_Algo::VertexNode( VV[0], meshDS );
250 double prevNormPar = ( i == 0 ? 0 : myNormPar[ i-1 ]); // normalized param
251 if ( node ) { // internal nodes may be missing
252 u2node.insert( make_pair( prevNormPar, node ));
255 MESSAGE(" NO NODE on VERTEX" );
259 // Put internal nodes
260 SMESHDS_SubMesh* sm = meshDS->MeshElements( myEdge[i] );
262 vector< pair< double, const SMDS_MeshNode*> > u2nodeVec;
263 u2nodeVec.reserve( sm->NbNodes() );
264 SMDS_NodeIteratorPtr nItr = sm->GetNodes();
265 double paramSize = myLast[i] - myFirst[i];
266 double r = myNormPar[i] - prevNormPar;
267 if ( !myIsUniform[i] )
268 while ( nItr->more() )
270 const SMDS_MeshNode* node = nItr->next();
271 if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
273 double u = helper.GetNodeU( myEdge[i], node, 0, ¶mOK );
274 double aLenU = GCPnts_AbscissaPoint::Length
275 ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), myFirst[i], u );
276 if ( myEdgeLength[i] < aLenU ) // nonregression test "3D_mesh_NETGEN/G6"
281 double normPar = prevNormPar + r*aLenU/myEdgeLength[i];
282 u2nodeVec.push_back( make_pair( normPar, node ));
284 nItr = sm->GetNodes();
285 if ( u2nodeVec.empty() )
286 while ( nItr->more() )
288 const SMDS_MeshNode* node = nItr->next();
289 if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
291 double u = helper.GetNodeU( myEdge[i], node, 0, ¶mOK );
293 // paramSize is signed so orientation is taken into account
294 double normPar = prevNormPar + r * ( u - myFirst[i] ) / paramSize;
295 u2nodeVec.push_back( make_pair( normPar, node ));
297 u2node.insert( u2nodeVec.begin(), u2nodeVec.end() );
299 // Put 2nd vertex node for a last edge
300 if ( i+1 == myEdge.size() ) {
301 node = SMESH_Algo::VertexNode( VV[1], meshDS );
303 MESSAGE(" NO NODE on VERTEX" );
306 u2node.insert( make_pair( 1., node ));
309 if ( u2node.size() != myNbPonits ) {
310 MESSAGE("Wrong node parameters on edges, u2node.size():"
311 <<u2node.size()<<" != myNbPonits:"<<myNbPonits);
315 // fill array of UVPtStruct
317 vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myPoints );
318 points->resize( myNbPonits );
321 double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
322 map< double, const SMDS_MeshNode*>::iterator u_node = u2node.begin();
323 for (int i = 0 ; u_node != u2node.end(); ++u_node, ++i ) {
324 UVPtStruct & uvPt = (*points)[i];
325 uvPt.node = u_node->second;
326 uvPt.x = uvPt.y = uvPt.normParam = u_node->first;
327 if ( isXConst ) uvPt.x = constValue;
328 else uvPt.y = constValue;
329 if ( myNormPar[ EdgeIndex ] < uvPt.normParam ) {
330 prevNormPar = myNormPar[ EdgeIndex ];
333 if ( EdgeIndex >= myEdge.size() ) {
335 MESSAGE ( "WRONg EdgeIndex " << 1+EdgeIndex
336 << " myNormPar.size()="<<myNormPar.size()
337 << " myNormPar["<< EdgeIndex<<"]="<< myNormPar[ EdgeIndex ]
338 << " uvPt.normParam="<<uvPt.normParam );
341 paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
343 const SMDS_EdgePosition* epos =
344 dynamic_cast<const SMDS_EdgePosition*>(uvPt.node->GetPosition().get());
346 uvPt.param = epos->GetUParameter();
349 double r = ( uvPt.normParam - prevNormPar )/ paramSize;
350 // uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
351 uvPt.param = ( r > 0.5 ? myLast[EdgeIndex] : myFirst[EdgeIndex] );
353 if ( !myC2d[ EdgeIndex ].IsNull() ) {
354 gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
359 uvPt.u = uvPt.v = 1e+100;
366 //================================================================================
368 * \brief Falsificate info on nodes
369 * \param nbSeg - nb of segments on the side
370 * \retval UVPtStruct* - array of data structures
372 //================================================================================
374 const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int nbSeg,
376 double constValue) const
378 if ( myFalsePoints.empty() ) {
380 if ( NbEdges() == 0 ) return myFalsePoints;
382 vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myFalsePoints );
383 points->resize( nbSeg+1 );
386 double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
387 for (int i = 0 ; i < myFalsePoints.size(); ++i ) {
388 double normPar = double(i) / double(nbSeg);
389 UVPtStruct & uvPt = (*points)[i];
391 uvPt.x = uvPt.y = uvPt.param = uvPt.normParam = normPar;
392 if ( isXConst ) uvPt.x = constValue;
393 else uvPt.y = constValue;
394 if ( myNormPar[ EdgeIndex ] < normPar ) {
395 prevNormPar = myNormPar[ EdgeIndex ];
397 paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
399 double r = ( normPar - prevNormPar )/ paramSize;
400 uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
401 if ( !myC2d[ EdgeIndex ].IsNull() ) {
402 gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
407 uvPt.u = uvPt.v = 1e+100;
411 return myFalsePoints;
413 // gp_Pnt StdMeshers_FaceSide::Value(double U) const
417 //================================================================================
419 * \brief reverse order of vector elements
420 * \param vec - vector to reverse
422 //================================================================================
424 template <typename T > void reverse(vector<T> & vec)
426 for ( int f=0, r=vec.size()-1; f < r; ++f, --r )
427 std::swap( vec[f], vec[r] );
430 //================================================================================
432 * \brief Change orientation of side geometry
434 //================================================================================
436 void StdMeshers_FaceSide::Reverse()
438 int nbEdges = myEdge.size();
439 for ( int i = nbEdges-1; i >= 0; --i ) {
440 std::swap( myFirst[i], myLast[i] );
442 if ( i > 0 ) // at the first loop 1. is overwritten
443 myNormPar[i] = 1 - myNormPar[i-1];
448 reverse( myC3dAdaptor );
451 reverse( myNormPar );
452 reverse( myEdgeLength );
453 reverse( myIsUniform );
455 myNormPar[nbEdges-1]=1.;
457 myFalsePoints.clear();
460 //================================================================================
462 * \brief Show side features
464 //================================================================================
466 void StdMeshers_FaceSide::dump(const char* msg) const
469 if (msg) MESSAGE ( std::endl << msg );
470 MESSAGE_BEGIN ("NB EDGES: "<< myEdge.size() );
471 MESSAGE_ADD ( "nbPoints: "<<myNbPonits<<" vecSize: " << myPoints.size()<<" "<<myFalsePoints.size() );
472 for ( int i=0; i<myEdge.size(); ++i)
474 MESSAGE_ADD ( "\t"<<i+1 );
475 MESSAGE_ADD ( "\tEDGE: " );
476 if (myEdge[i].IsNull()) {
477 MESSAGE_ADD ( "NULL" );
480 TopAbs::Print(myEdge[i].Orientation(),cout)<<" "<<myEdge[i].TShape().operator->()<<endl;
481 MESSAGE_ADD ( "\tV1: " << TopExp::FirstVertex( myEdge[i], 1).TShape().operator->()
482 << " V2: " << TopExp::LastVertex( myEdge[i], 1).TShape().operator->() );
484 MESSAGE_ADD ( "\tC2d: ");
486 if (myC2d[i].IsNull()) {
487 MESSAGE_ADD ( "NULL" );
490 MESSAGE_ADD ( myC2d[i].operator->() );
493 MESSAGE_ADD ( "\tF: "<<myFirst[i]<< " L: "<< myLast[i] );
494 MESSAGE_END ( "\tnormPar: "<<myNormPar[i]<<endl );
499 //================================================================================
501 * \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
502 * \retval Adaptor2d_Curve2d* -
504 //================================================================================
506 struct Adaptor2dCurve2d : public Adaptor2d_Curve2d
508 const StdMeshers_FaceSide* mySide;
509 Adaptor2dCurve2d(const StdMeshers_FaceSide* faceSide):mySide(faceSide) {}
510 gp_Pnt2d Value(const Standard_Real U) const { return mySide->Value2d( U ); }
511 Standard_Real FirstParameter() const { return 0; }
512 Standard_Real LastParameter() const { return 1; }
515 Adaptor2d_Curve2d* StdMeshers_FaceSide::GetCurve2d() const
517 return new Adaptor2dCurve2d( this );
520 //================================================================================
522 * \brief Creates a fully functional Adaptor_Curve
524 //================================================================================
526 BRepAdaptor_CompCurve* StdMeshers_FaceSide::GetCurve3d() const
528 if ( myEdge.empty() )
531 // if ( myEdge.size() == 1 )
532 // return new BRepAdaptor_Curve( myEdge[0] );
535 BRep_Builder aBuilder;
536 aBuilder.MakeWire(aWire);
537 for ( int i=0; i<myEdge.size(); ++i )
538 aBuilder.Add( aWire, myEdge[i] );
539 return new BRepAdaptor_CompCurve( aWire );
543 //================================================================================
545 * \brief Return 2D point by normalized parameter
546 * \param U - normalized parameter value
547 * \retval gp_Pnt2d - point
549 //================================================================================
551 gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
553 if ( !myC2d[0].IsNull() ) {
554 int i = EdgeIndex( U );
555 double prevU = i ? myNormPar[ i-1 ] : 0;
556 double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
558 double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
560 // check parametrization of curve
561 if( !myIsUniform[i] )
563 double aLen3dU = r * myEdgeLength[i] * ( myFirst[i]>myLast[i] ? -1. : 1.);
564 GCPnts_AbscissaPoint AbPnt
565 ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
566 if( AbPnt.IsDone() ) {
567 par = AbPnt.Parameter();
570 return myC2d[ i ]->Value(par);
573 return myDefaultPnt2d;
576 //================================================================================
578 * \brief Return wires of a face as StdMeshers_FaceSide's
580 //================================================================================
582 TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face& theFace,
583 SMESH_Mesh & theMesh,
584 const bool theIgnoreMediumNodes,
588 list< TopoDS_Edge > edges, internalEdges;
589 list< int > nbEdgesInWires;
590 int nbWires = SMESH_Block::GetOrderedEdges (theFace, V1, edges, nbEdgesInWires);
592 // split list of all edges into separate wires
593 TSideVector wires( nbWires );
594 list< int >::iterator nbE = nbEdgesInWires.begin();
595 list< TopoDS_Edge >::iterator from = edges.begin(), to = from;
596 for ( int iW = 0; iW < nbWires; ++iW, ++nbE )
598 std::advance( to, *nbE );
599 if ( *nbE == 0 ) // Issue 0020676
603 wires.resize( nbWires );
606 list< TopoDS_Edge > wireEdges( from, to );
607 // assure that there is a node on the first vertex
608 // as StdMeshers_FaceSide::GetUVPtStruct() requires
609 if ( wireEdges.front().Orientation() != TopAbs_INTERNAL ) // Issue 0020676
611 while ( !SMESH_Algo::VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
612 theMesh.GetMeshDS()))
614 wireEdges.splice(wireEdges.end(), wireEdges,
615 wireEdges.begin(), ++wireEdges.begin());
616 if ( from->IsSame( wireEdges.front() )) {
618 ( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"No nodes on vertices"));
619 return TSideVector(0);
623 else if ( *nbE > 1 ) // Issue 0020676 (Face_pb_netgen.brep) - several internal edges in a wire
625 internalEdges.splice( internalEdges.end(), wireEdges, ++wireEdges.begin(), wireEdges.end());
628 StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, wireEdges, &theMesh,
629 /*isForward=*/true, theIgnoreMediumNodes);
630 wires[ iW ] = StdMeshers_FaceSidePtr( wire );
633 while ( !internalEdges.empty() )
635 StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, internalEdges.back(), &theMesh,
636 /*isForward=*/true, theIgnoreMediumNodes);
637 wires.push_back( StdMeshers_FaceSidePtr( wire ));
638 internalEdges.pop_back();
643 //================================================================================
645 * \brief Return 1st vertex of the i-the edge
647 //================================================================================
649 TopoDS_Vertex StdMeshers_FaceSide::FirstVertex(int i) const
654 v = myEdge[i].Orientation() <= TopAbs_REVERSED ? // FORWARD || REVERSED
655 TopExp::FirstVertex( myEdge[i], 1 ) :
656 TopoDS::Vertex( TopoDS_Iterator( myEdge[i] ).Value() );
661 //================================================================================
663 * \brief Return last vertex of the i-the edge
665 //================================================================================
667 TopoDS_Vertex StdMeshers_FaceSide::LastVertex(int i) const
672 const TopoDS_Edge& edge = i<0 ? myEdge[ NbEdges() + i ] : myEdge[i];
673 if ( edge.Orientation() <= TopAbs_REVERSED ) // FORWARD || REVERSED
674 v = TopExp::LastVertex( edge, 1 );
676 for ( TopoDS_Iterator vIt( edge ); vIt.More(); vIt.Next() )
677 v = TopoDS::Vertex( vIt.Value() );