1 // Copyright (C) 2007-2015 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, or (at your option) any later version.
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
22 // File : SMESH_MAT2d.cxx
23 // Created : Thu May 28 17:49:53 2015
24 // Author : Edward AGAPOV (eap)
26 #include "SMESH_MAT2d.hxx"
30 #include <BRepAdaptor_CompCurve.hxx>
31 #include <BRepBuilderAPI_MakeEdge.hxx>
32 #include <BRepBuilderAPI_MakeVertex.hxx>
33 #include <BRep_Builder.hxx>
34 #include <BRep_Tool.hxx>
35 #include <Bnd_B2d.hxx>
36 //#include <GCPnts_AbscissaPoint.hxx>
37 #include <GCPnts_TangentialDeflection.hxx>
38 // #include <GCPnts_UniformAbscissa.hxx>
39 // #include <GCPnts_UniformDeflection.hxx>
40 #include <Geom2d_Curve.hxx>
41 //#include <GeomAdaptor_Curve.hxx>
42 #include <Geom2dAdaptor_Curve.hxx>
43 #include <Geom_Curve.hxx>
44 #include <Geom_Surface.hxx>
46 #include <TopoDS_Vertex.hxx>
47 #include <TopoDS_Wire.hxx>
51 #include "SMESH_File.hxx"
52 #include "SMESH_Comment.hxx"
56 using boost::polygon::x;
57 using boost::polygon::y;
58 using SMESH_MAT2d::TVD;
59 using SMESH_MAT2d::TVDEdge;
60 using SMESH_MAT2d::TVDCell;
61 using SMESH_MAT2d::TVDVertex;
65 // Input data for construct_voronoi()
66 // -------------------------------------------------------------------------------------
70 int _a, _b; // coordinates
71 double _param; // param on EDGE
72 InPoint(int x, int y, double param) : _a(x), _b(y), _param(param) {}
73 InPoint() : _a(0), _b(0), _param(0) {}
76 list< const TVDEdge* > _edges; // MA edges of a concave InPoint in CCW order
78 size_t index( const vector< InPoint >& inPoints ) const { return this - &inPoints[0]; }
79 bool operator==( const InPoint& other ) const { return _a == other._a && _b == other._b; }
80 bool operator==( const TVDVertex* v ) const { return ( Abs( _a - v->x() ) < 1. &&
81 Abs( _b - v->y() ) < 1. ); }
83 // -------------------------------------------------------------------------------------
91 size_t _geomEdgeInd; // EDGE index within the FACE
93 list< const TVDEdge* > _edges; // MA edges in CCW order within _cell
95 InSegment( InPoint * p0, InPoint * p1, size_t iE)
96 : _p0(p0), _p1(p1), _geomEdgeInd(iE) {}
97 InSegment() : _p0(0), _p1(0), _geomEdgeInd(0) {}
99 const InPoint& point0() const { return *_p0; }
100 const InPoint& point1() const { return *_p1; }
102 inline bool isConnected( const TVDEdge* edge );
104 inline bool isExternal( const TVDEdge* edge );
106 static void setGeomEdgeToCell( const TVDCell* cell, size_t eID ) { cell->color( eID ); }
108 static size_t getGeomEdge( const TVDCell* cell ) { return cell->color(); }
111 // check if a TVDEdge begins at my end or ends at my start
112 inline bool InSegment::isConnected( const TVDEdge* edge )
114 return (( edge->vertex0() && edge->vertex1() )
116 ((Abs( edge->vertex0()->x() - _p1->_a ) < 1.&&
117 Abs( edge->vertex0()->y() - _p1->_b ) < 1. ) ||
118 (Abs( edge->vertex1()->x() - _p0->_a ) < 1.&&
119 Abs( edge->vertex1()->y() - _p0->_b ) < 1. )));
122 // check if a MA TVDEdge is outside of a domain
123 inline bool InSegment::isExternal( const TVDEdge* edge )
125 double dot = // x1*x2 + y1*y2; (x1,y1) - internal normal of InSegment
126 ( _p0->_b - _p1->_b ) * ( 0.5 * ( edge->vertex0()->x() + edge->vertex1()->x() ) - _p0->_a ) +
127 ( _p1->_a - _p0->_a ) * ( 0.5 * ( edge->vertex0()->y() + edge->vertex1()->y() ) - _p0->_b );
131 // // -------------------------------------------------------------------------------------
132 // const size_t theExternMA = 111; // to mark external MA edges
134 // bool isExternal( const TVDEdge* edge )
136 // return ( SMESH_MAT2d::Branch::getBndSegment( edge ) == theExternMA );
139 // // mark external MA edges
140 // void markExternalEdges( const TVDEdge* edge )
142 // if ( isExternal( edge ))
144 // SMESH_MAT2d::Branch::setBndSegment( theExternMA, edge );
145 // SMESH_MAT2d::Branch::setBndSegment( theExternMA, edge->twin() );
146 // if ( edge->is_primary() && edge->vertex1() )
148 // const TVDVertex * v = edge->vertex1();
149 // edge = v->incident_edge();
151 // markExternalEdges( edge );
152 // edge = edge->rot_next();
153 // } while ( edge != v->incident_edge() );
157 // -------------------------------------------------------------------------------------
159 // writes segments into a txt file readable by voronoi_visualizer
160 void inSegmentsToFile( vector< InSegment>& inSegments)
162 if ( inSegments.size() > 1000 )
164 const char* fileName = "/misc/dn25/salome/eap/salome/misc/Code/C++/MAdebug.txt";
165 SMESH_File file(fileName, false );
167 file.openForWriting();
169 text << "0\n"; // nb points
170 text << inSegments.size() << "\n"; // nb segments
171 for ( size_t i = 0; i < inSegments.size(); ++i )
173 text << inSegments[i]._p0->_a << " "
174 << inSegments[i]._p0->_b << " "
175 << inSegments[i]._p1->_a << " "
176 << inSegments[i]._p1->_b << "\n";
179 file.write( text.c_str(), text.size() );
180 cout << "Write " << fileName << endl;
182 void dumpEdge( const TVDEdge* edge )
184 cout << "*Edge_" << edge;
185 if ( !edge->vertex0() )
186 cout << " ( INF, INF";
188 cout << " ( " << edge->vertex0()->x() << ", " << edge->vertex0()->y();
189 if ( !edge->vertex1() )
190 cout << ") -> ( INF, INF";
192 cout << ") -> ( " << edge->vertex1()->x() << ", " << edge->vertex1()->y();
193 cout << ")\t cell=" << edge->cell()
194 << " iBnd=" << edge->color()
195 << " twin=" << edge->twin()
196 << " twin_cell=" << edge->twin()->cell()
197 << " prev=" << edge->prev() << " next=" << edge->next()
198 << ( edge->is_primary() ? " MA " : " SCND" )
199 << ( edge->is_linear() ? " LIN " : " CURV" )
202 void dumpCell( const TVDCell* cell )
204 cout << "**Cell_" << cell << " GEOM=" << cell->color() << " ";
205 cout << ( cell->contains_segment() ? " SEG " : " PNT " );
206 if ( cell-> is_degenerate() )
211 const TVDEdge* edge = cell->incident_edge();
215 cout << " - " << ++i << " ";
217 } while (edge != cell->incident_edge());
221 void inSegmentsToFile( vector< InSegment>& inSegments) {}
222 void dumpEdge( const TVDEdge* edge ) {}
223 void dumpCell( const TVDCell* cell ) {}
226 // -------------------------------------------------------------------------------------
232 struct geometry_concept<InPoint> {
233 typedef point_concept type;
236 struct point_traits<InPoint> {
237 typedef int coordinate_type;
239 static inline coordinate_type get(const InPoint& point, orientation_2d orient) {
240 return (orient == HORIZONTAL) ? point._a : point._b;
245 struct geometry_concept<InSegment> {
246 typedef segment_concept type;
250 struct segment_traits<InSegment> {
251 typedef int coordinate_type;
252 typedef InPoint point_type;
254 static inline point_type get(const InSegment& segment, direction_1d dir) {
255 return *(dir.to_int() ? segment._p1 : segment._p0);
258 } // namespace polygon
260 // -------------------------------------------------------------------------------------
264 const int theNoBrachID = 0;
265 double theScale[2]; // scale used in bndSegsToMesh()
266 const size_t theNoEdgeID = std::numeric_limits<size_t>::max() / 1000;
268 // -------------------------------------------------------------------------------------
270 * \brief Intermediate DS to create InPoint's
276 UVU( gp_Pnt2d uv, double u ): _uv(uv), _u(u) {}
277 InPoint getInPoint( double scale[2] )
279 return InPoint( int( _uv.X() * scale[0]), int( _uv.Y() * scale[1]), _u );
282 // -------------------------------------------------------------------------------------
284 * \brief Segment of EDGE, used to create BndPoints
289 const TVDEdge* _edge;
291 BndSeg* _prev; // previous BndSeg in FACE boundary
292 int _branchID; // negative ID means reverse direction
294 BndSeg( InSegment* seg, const TVDEdge* edge, double u ):
295 _inSeg(seg), _edge(edge), _uLast(u), _prev(0), _branchID( theNoBrachID ) {}
297 void setIndexToEdge( size_t id )
299 SMESH_MAT2d::Branch::setBndSegment( id, _edge );
302 int branchID() const { return Abs( _branchID ); }
304 size_t geomEdge() const { return _inSeg->_geomEdgeInd; }
306 static BndSeg* getBndSegOfEdge( const TVDEdge* edge,
307 vector< vector< BndSeg > >& bndSegsPerEdge )
312 size_t oppSegIndex = SMESH_MAT2d::Branch::getBndSegment( edge );
313 size_t oppEdgeIndex = SMESH_MAT2d::Branch::getGeomEdge ( edge );
314 if ( oppEdgeIndex < bndSegsPerEdge.size() &&
315 oppSegIndex < bndSegsPerEdge[ oppEdgeIndex ].size() )
317 seg = & bndSegsPerEdge[ oppEdgeIndex ][ oppSegIndex ];
323 void setBranch( int branchID, vector< vector< BndSeg > >& bndSegsPerEdge )
325 _branchID = branchID;
327 // pass branch to an opposite BndSeg
329 if ( BndSeg* oppSeg = getBndSegOfEdge( _edge->twin(), bndSegsPerEdge ))
331 if ( oppSeg->_branchID == theNoBrachID )
332 oppSeg->_branchID = -branchID;
335 bool hasOppositeEdge()
337 if ( !_edge ) return false;
338 return ( _inSeg->getGeomEdge( _edge->twin()->cell() ) != theNoEdgeID );
341 // check a next segment in CW order
342 bool isSameBranch( const BndSeg& seg2 )
344 if ( !_edge || !seg2._edge )
347 if ( _edge->twin() == seg2._edge )
350 const TVDCell* cell1 = this->_edge->twin()->cell();
351 const TVDCell* cell2 = seg2. _edge->twin()->cell();
352 if ( cell1 == cell2 )
355 const TVDEdge* edgeMedium1 = this->_edge->twin()->next();
356 const TVDEdge* edgeMedium2 = seg2. _edge->twin()->prev();
358 if ( edgeMedium1->is_secondary() && edgeMedium2->is_secondary() )
360 if ( edgeMedium1->twin() == edgeMedium2 )
362 // edgeMedium's are edges whose twin()->cell is built on an end point of inSegment
363 // and is located between cell1 and cell2
364 if ( edgeMedium1->twin() == edgeMedium2->twin() ) // is this possible???
366 if ( edgeMedium1->twin() == edgeMedium2->twin()->next() &&
367 edgeMedium1->twin()->cell()->contains_point() )
370 else if ( edgeMedium1->is_primary() && edgeMedium2->is_primary() )
372 if ( edgeMedium1->twin() == edgeMedium2 &&
373 SMESH_MAT2d::Branch::getGeomEdge( edgeMedium1 ) ==
374 SMESH_MAT2d::Branch::getGeomEdge( edgeMedium2 ))
375 // this is an ignored MA edge between inSegment's on one EDGE forming a convex corner
383 // -------------------------------------------------------------------------------------
387 struct BranchIterator
390 const std::vector<const TVDEdge*> & _edges;
393 BranchIterator(const std::vector<const TVDEdge*> & edges, int i )
394 :_i( i ), _size( edges.size() ), _edges( edges )
396 _closed = ( edges[0]->vertex1() == edges.back()->vertex0() || // closed branch
397 edges[0]->vertex0() == edges.back()->vertex1() );
399 const TVDEdge* operator++() { ++_i; return edge(); }
400 const TVDEdge* operator--() { --_i; return edge(); }
401 bool operator<( const BranchIterator& other ) { return _i < other._i; }
402 BranchIterator& operator=( const BranchIterator& other ) { _i = other._i; return *this; }
403 void set(int i) { _i = i; }
404 int index() const { return _i; }
405 int indexMod() const { return ( _i + _size ) % _size; }
406 const TVDEdge* edge() const {
407 return _closed ? _edges[ indexMod() ] : ( _i < 0 || _i >= _size ) ? 0 : _edges[ _i ];
409 const TVDEdge* edgePrev() { --_i; const TVDEdge* e = edge(); ++_i; return e; }
410 const TVDEdge* edgeNext() { ++_i; const TVDEdge* e = edge(); --_i; return e; }
413 //================================================================================
415 * \brief debug: to visually check found MA edges
417 //================================================================================
419 void bndSegsToMesh( const vector< vector< BndSeg > >& bndSegsPerEdge )
422 if ( !getenv("bndSegsToMesh")) return;
423 map< const TVDVertex *, int > v2Node;
424 map< const TVDVertex *, int >::iterator v2n;
425 set< const TVDEdge* > addedEdges;
427 const char* fileName = "/misc/dn25/salome/eap/salome/misc/Code/C++/MAedges.py";
428 SMESH_File file(fileName, false );
430 file.openForWriting();
432 text << "import salome, SMESH\n";
433 text << "salome.salome_init()\n";
434 text << "from salome.smesh import smeshBuilder\n";
435 text << "smesh = smeshBuilder.New(salome.myStudy)\n";
436 text << "m=smesh.Mesh()\n";
437 for ( size_t iE = 0; iE < bndSegsPerEdge.size(); ++iE )
439 const vector< BndSeg >& bndSegs = bndSegsPerEdge[ iE ];
440 for ( size_t i = 0; i < bndSegs.size(); ++i )
442 if ( !bndSegs[i]._edge )
443 text << "# E=" << iE << " i=" << i << " NULL edge\n";
444 else if ( !bndSegs[i]._edge->vertex0() ||
445 !bndSegs[i]._edge->vertex1() )
446 text << "# E=" << iE << " i=" << i << " INFINITE edge\n";
447 else if ( addedEdges.insert( bndSegs[i]._edge ).second &&
448 addedEdges.insert( bndSegs[i]._edge->twin() ).second )
450 v2n = v2Node.insert( make_pair( bndSegs[i]._edge->vertex0(), v2Node.size() + 1 )).first;
451 size_t n0 = v2n->second;
452 if ( n0 == v2Node.size() )
453 text << "n" << n0 << " = m.AddNode( "
454 << bndSegs[i]._edge->vertex0()->x() / theScale[0] << ", "
455 << bndSegs[i]._edge->vertex0()->y() / theScale[1] << ", 0 )\n";
457 v2n = v2Node.insert( make_pair( bndSegs[i]._edge->vertex1(), v2Node.size() + 1 )).first;
458 size_t n1 = v2n->second;
459 if ( n1 == v2Node.size() )
460 text << "n" << n1 << " = m.AddNode( "
461 << bndSegs[i]._edge->vertex1()->x() / theScale[0] << ", "
462 << bndSegs[i]._edge->vertex1()->y() / theScale[1] << ", 0 )\n";
464 text << "e" << i << " = m.AddEdge([ n" << n0 << ", n" << n1 << " ])\n";
469 file.write( text.c_str(), text.size() );
470 cout << "execfile( '" << fileName << "')" << endl;
474 //================================================================================
476 * \brief Computes length of a TVDEdge
478 //================================================================================
480 double length( const TVDEdge* edge )
482 gp_XY d( edge->vertex0()->x() - edge->vertex1()->x(),
483 edge->vertex0()->y() - edge->vertex1()->y() );
487 //================================================================================
489 * \brief Compute scale to have the same 2d proportions as in 3d
491 //================================================================================
493 void computeProportionScale( const TopoDS_Face& face,
494 const Bnd_B2d& uvBox,
497 scale[0] = scale[1] = 1.;
498 if ( uvBox.IsVoid() ) return;
501 Handle(Geom_Surface) surface = BRep_Tool::Surface( face, loc );
503 const int nbDiv = 30;
504 gp_XY uvMin = uvBox.CornerMin(), uvMax = uvBox.CornerMax();
505 gp_XY uvMid = 0.5 * ( uvMin + uvMax );
506 double du = ( uvMax.X() - uvMin.X() ) / nbDiv;
507 double dv = ( uvMax.Y() - uvMin.Y() ) / nbDiv;
509 double uLen3d = 0, vLen3d = 0;
510 gp_Pnt uPrevP = surface->Value( uvMin.X(), uvMid.Y() );
511 gp_Pnt vPrevP = surface->Value( uvMid.X(), uvMin.Y() );
512 for (int i = 1; i <= nbDiv; i++)
514 double u = uvMin.X() + du * i;
515 double v = uvMin.Y() + dv * i;
516 gp_Pnt uP = surface->Value( u, uvMid.Y() );
517 gp_Pnt vP = surface->Value( uvMid.X(), v );
518 uLen3d += uP.Distance( uPrevP );
519 vLen3d += vP.Distance( vPrevP );
523 scale[0] = uLen3d / ( uvMax.X() - uvMin.X() );
524 scale[1] = vLen3d / ( uvMax.Y() - uvMin.Y() );
527 //================================================================================
529 * \brief Fill input data for construct_voronoi()
531 //================================================================================
533 bool makeInputData(const TopoDS_Face& face,
534 const std::vector< TopoDS_Edge >& edges,
535 const double minSegLen,
536 vector< InPoint >& inPoints,
537 vector< InSegment>& inSegments,
540 const double theDiscrCoef = 0.5; // to decrease minSegLen for discretization
543 // discretize the EDGEs to get 2d points and segments
545 vector< vector< UVU > > uvuVec( edges.size() );
547 for ( size_t iE = 0; iE < edges.size(); ++iE )
549 vector< UVU > & points = uvuVec[ iE ];
552 Handle(Geom_Curve) c3d = BRep_Tool::Curve ( edges[ iE ], loc, f, l );
553 Handle(Geom2d_Curve) c2d = BRep_Tool::CurveOnSurface( edges[ iE ], face, f, l );
554 if ( c2d.IsNull() ) return false;
556 points.push_back( UVU( c2d->Value( f ), f ));
557 uvBox.Add( points.back()._uv );
559 Geom2dAdaptor_Curve c2dAdaptor (c2d, f,l );
560 double curDeflect = 0.3; //0.01; //Curvature deflection
561 double angDeflect = 0.2; // 0.09; //Angular deflection
563 GCPnts_TangentialDeflection discret(c2dAdaptor, angDeflect, curDeflect);
564 // if ( discret.NbPoints() > 2 )
569 // discret.Initialize( c2dAdaptor, 100, curDeflect );
570 // cout << "C " << curDeflect << " " << discret.NbPoints() << endl;
571 // curDeflect *= 1.5;
573 // while ( discret.NbPoints() > 5 );
577 // discret.Initialize( c2dAdaptor, angDeflect, 100 );
578 // cout << "A " << angDeflect << " " << discret.NbPoints() << endl;
579 // angDeflect *= 1.5;
581 // while ( discret.NbPoints() > 5 );
585 pPrev = c3d->Value( f );
586 if ( discret.NbPoints() > 2 )
587 for ( int i = 2; i <= discret.NbPoints(); i++ ) // skip the 1st point
589 double u = discret.Parameter(i);
593 int nbDiv = int( p.Distance( pPrev ) / minSegLen / theDiscrCoef );
594 double dU = ( u - points.back()._u ) / nbDiv;
595 for ( int iD = 1; iD < nbDiv; ++iD )
597 double uD = points.back()._u + dU;
598 points.push_back( UVU( c2d->Value( uD ), uD ));
602 points.push_back( UVU( c2d->Value( u ), u ));
603 uvBox.Add( points.back()._uv );
605 // if ( !c3d.IsNull() )
607 // vector<double> params;
608 // GeomAdaptor_Curve c3dAdaptor( c3d,f,l );
611 // const double deflection = minSegLen * 0.1;
612 // GCPnts_UniformDeflection discret( c3dAdaptor, deflection, f, l, true );
613 // if ( !discret.IsDone() )
615 // int nbP = discret.NbPoints();
616 // for ( int i = 2; i < nbP; i++ ) // skip 1st and last points
617 // params.push_back( discret.Parameter(i) );
621 // double eLen = GCPnts_AbscissaPoint::Length( c3dAdaptor );
622 // int nbSeg = Max( 1, int( eLen / minSegLen / theDiscrCoef ));
623 // double segLen = eLen / nbSeg;
624 // GCPnts_UniformAbscissa discret( c3dAdaptor, segLen, f, l );
625 // int nbP = Min( discret.NbPoints(), nbSeg + 1 );
626 // for ( int i = 2; i < nbP; i++ ) // skip 1st and last points
627 // params.push_back( discret.Parameter(i) );
629 // for ( size_t i = 0; i < params.size(); ++i )
631 // points.push_back( UVU( c2d->Value( params[i] ), params[i] ));
632 // uvBox.Add( points.back()._uv );
635 if ( points.size() < 2 )
637 points.push_back( UVU( c2d->Value( l ), l ));
638 uvBox.Add( points.back()._uv );
640 if ( edges[ iE ].Orientation() == TopAbs_REVERSED )
641 std::reverse( points.begin(), points.end() );
644 // make connected EDGEs have same UV at shared VERTEX
645 TopoDS_Vertex vShared;
646 for ( size_t iE = 0; iE < edges.size(); ++iE )
648 size_t iE2 = (iE+1) % edges.size();
649 if ( !TopExp::CommonVertex( edges[iE], edges[iE2], vShared ))
651 if ( !vShared.IsSame( TopExp::LastVertex( edges[iE], true )))
653 vector< UVU > & points1 = uvuVec[ iE ];
654 vector< UVU > & points2 = uvuVec[ iE2 ];
655 gp_Pnt2d & uv1 = points1.back() ._uv;
656 gp_Pnt2d & uv2 = points2.front()._uv;
657 uv1 = uv2 = 0.5 * ( uv1.XY() + uv2.XY() );
660 // get scale to have the same 2d proportions as in 3d
661 computeProportionScale( face, uvBox, scale );
663 // make scale to have coordinates precise enough when converted to int
665 gp_XY uvMin = uvBox.CornerMin(), uvMax = uvBox.CornerMax();
666 uvMin.ChangeCoord(1) = uvMin.X() * scale[0];
667 uvMin.ChangeCoord(2) = uvMin.Y() * scale[1];
668 uvMax.ChangeCoord(1) = uvMax.X() * scale[0];
669 uvMax.ChangeCoord(2) = uvMax.Y() * scale[1];
670 double vMax[2] = { Max( Abs( uvMin.X() ), Abs( uvMax.X() )),
671 Max( Abs( uvMin.Y() ), Abs( uvMax.Y() )) };
672 int iMax = ( vMax[0] > vMax[1] ) ? 0 : 1;
673 const double precision = 1e-5;
674 double preciScale = Min( vMax[iMax] / precision,
675 std::numeric_limits<int>::max() / vMax[iMax] );
676 preciScale /= scale[iMax];
677 double roundedScale = 10; // to ease debug
678 while ( roundedScale * 10 < preciScale )
680 scale[0] *= roundedScale;
681 scale[1] *= roundedScale;
683 // create input points and segments
688 for ( size_t iE = 0; iE < uvuVec.size(); ++iE )
689 nbPnt += uvuVec[ iE ].size();
690 inPoints.resize( nbPnt );
691 inSegments.reserve( nbPnt );
694 if ( face.Orientation() == TopAbs_REVERSED )
696 for ( int iE = uvuVec.size()-1; iE >= 0; --iE )
698 vector< UVU > & points = uvuVec[ iE ];
699 inPoints[ iP++ ] = points.back().getInPoint( scale );
700 for ( size_t i = points.size()-1; i >= 1; --i )
702 inPoints[ iP++ ] = points[i-1].getInPoint( scale );
703 inSegments.push_back( InSegment( & inPoints[ iP-2 ], & inPoints[ iP-1 ], iE ));
709 for ( size_t iE = 0; iE < uvuVec.size(); ++iE )
711 vector< UVU > & points = uvuVec[ iE ];
712 inPoints[ iP++ ] = points[0].getInPoint( scale );
713 for ( size_t i = 1; i < points.size(); ++i )
715 inPoints[ iP++ ] = points[i].getInPoint( scale );
716 inSegments.push_back( InSegment( & inPoints[ iP-2 ], & inPoints[ iP-1 ], iE ));
721 theScale[0] = scale[0];
722 theScale[1] = scale[1];
727 //================================================================================
729 * \brief Update a branch joined to another one
731 //================================================================================
733 void updateJoinedBranch( vector< const TVDEdge* > & branchEdges,
735 vector< vector< BndSeg > > & bndSegs,
741 for ( size_t i = 0; i < branchEdges.size(); ++i )
743 if (( seg1 = BndSeg::getBndSegOfEdge( branchEdges[i], bndSegs )) &&
744 ( seg2 = BndSeg::getBndSegOfEdge( branchEdges[i]->twin(), bndSegs )))
746 seg1->_branchID /= seg1->branchID();
747 seg2->_branchID /= seg2->branchID();
748 seg1->_branchID *= -newID;
749 seg2->_branchID *= -newID;
750 branchEdges[i] = branchEdges[i]->twin();
753 std::reverse( branchEdges.begin(), branchEdges.end() );
757 for ( size_t i = 0; i < branchEdges.size(); ++i )
759 if (( seg1 = BndSeg::getBndSegOfEdge( branchEdges[i], bndSegs )) &&
760 ( seg2 = BndSeg::getBndSegOfEdge( branchEdges[i]->twin(), bndSegs )))
762 seg1->_branchID /= seg1->branchID();
763 seg2->_branchID /= seg2->branchID();
764 seg1->_branchID *= newID;
765 seg2->_branchID *= newID;
771 //================================================================================
773 * \brief Create MA branches and FACE boundary data
774 * \param [in] vd - voronoi diagram of \a inSegments
775 * \param [in] inPoints - FACE boundary points
776 * \param [in,out] inSegments - FACE boundary segments
777 * \param [out] branch - MA branches to fill
778 * \param [out] branchEnd - ends of MA branches to fill
779 * \param [out] boundary - FACE boundary to fill
781 //================================================================================
783 void makeMA( const TVD& vd,
784 const bool ignoreCorners,
785 vector< InPoint >& inPoints,
786 vector< InSegment > & inSegments,
787 vector< SMESH_MAT2d::Branch >& branch,
788 vector< const SMESH_MAT2d::BranchEnd* >& branchPnt,
789 SMESH_MAT2d::Boundary& boundary )
791 // Associate MA cells with geom EDGEs
792 for (TVD::const_cell_iterator it = vd.cells().begin(); it != vd.cells().end(); ++it)
794 const TVDCell* cell = &(*it);
795 if ( cell->contains_segment() )
797 InSegment& seg = inSegments[ cell->source_index() ];
799 seg.setGeomEdgeToCell( cell, seg._geomEdgeInd );
803 InSegment::setGeomEdgeToCell( cell, theNoEdgeID );
807 vector< bool > inPntChecked( inPoints.size(), false );
809 // Find MA edges of each inSegment
811 for ( size_t i = 0; i < inSegments.size(); ++i )
813 InSegment& inSeg = inSegments[i];
815 // get edges around the cell lying on MA
816 bool hasSecondary = false;
817 const TVDEdge* edge = inSeg._cell->incident_edge();
819 edge = edge->next(); // Returns the CCW next edge within the cell.
820 if ( edge->is_primary() && !inSeg.isExternal( edge ) )
821 inSeg._edges.push_back( edge ); // edge equidistant from two InSegments
824 } while (edge != inSeg._cell->incident_edge());
826 // there can be several continuous MA edges but maEdges can begin in the middle of
827 // a chain of continuous MA edges. Make the chain continuous.
828 list< const TVDEdge* >& maEdges = inSeg._edges;
829 if ( maEdges.empty() )
832 while ( maEdges.back()->next() == maEdges.front() )
833 maEdges.splice( maEdges.end(), maEdges, maEdges.begin() );
835 // remove maEdges equidistant from two neighbor InSegments of the same geom EDGE
836 list< const TVDEdge* >::iterator e = maEdges.begin();
837 while ( e != maEdges.end() )
839 const TVDCell* cell2 = (*e)->twin()->cell(); // cell on the other side of a MA edge
840 size_t geoE2 = InSegment::getGeomEdge( cell2 );
841 bool toRemove = ( inSeg._geomEdgeInd == geoE2 && inSeg.isConnected( *e ));
843 e = maEdges.erase( e );
847 if ( maEdges.empty() )
850 // add MA edges corresponding to concave InPoints
851 for ( int is2nd = 0; is2nd < 2; ++is2nd ) // loop on two ends of inSeg
853 InPoint& inPnt = *( is2nd ? inSeg._p1 : inSeg._p0 );
854 size_t pInd = inPnt.index( inPoints );
855 if ( inPntChecked[ pInd ] )
858 inPntChecked[ pInd-1 ] &&
859 inPoints[ pInd-1 ] == inPnt )
861 inPntChecked[ pInd ] = true;
863 const TVDEdge* maE = is2nd ? maEdges.front() : maEdges.back();
864 if ( inPnt == ( is2nd ? maE->vertex0() : maE->vertex1() ))
866 const TVDEdge* edge = // a secondary TVDEdge connecting inPnt and maE
867 is2nd ? maE->prev() : maE->next();
868 while ( inSeg.isConnected( edge ))
870 if ( edge->is_primary() ) break; // this should not happen
871 const TVDEdge* edge2 = edge->twin(); // we are in a neighbor cell, add MA edges to inPnt
872 if ( inSeg.getGeomEdge( edge2->cell() ) != theNoEdgeID )
873 break; // cell of an InSegment
874 bool hasInfinite = false;
875 list< const TVDEdge* > pointEdges;
879 edge = edge->next(); // Returns the CCW next edge within the cell.
880 if ( edge->is_infinite() )
882 else if ( edge->is_primary() && !inSeg.isExternal( edge ))
883 pointEdges.push_back( edge );
885 while ( edge != edge2 && !hasInfinite );
887 if ( hasInfinite || pointEdges.empty() )
889 inPnt._edges.splice( inPnt._edges.end(), pointEdges );
890 inSeg.setGeomEdgeToCell( edge->cell(), inSeg._geomEdgeInd );
892 edge = is2nd ? inPnt._edges.front()->prev() : inPnt._edges.back()->next();
894 } // add MA edges corresponding to concave InPoints
896 } // loop on inSegments to find corresponding MA edges
899 // -------------------------------------------
900 // Create Branches and BndPoints for each EDGE
901 // -------------------------------------------
903 if ( inPoints.front() == inPoints.back() /*&& !inPoints[0]._edges.empty()*/ )
905 inPntChecked[0] = false; // do not use the 1st point twice
906 //InSegment::setGeomEdgeToCell( inPoints[0]._edges.back()->cell(), theNoEdgeID );
907 inPoints[0]._edges.clear();
910 // Divide InSegment's into BndSeg's (so that each BndSeg corresponds to one MA edge)
912 vector< vector< BndSeg > > bndSegsPerEdge( boundary.nbEdges() ); // all BndSeg's
914 vector< BndSeg > bndSegs; // bndSeg's of a current EDGE
915 size_t prevGeomEdge = theNoEdgeID;
917 list< const TVDEdge* >::reverse_iterator e;
918 for ( size_t i = 0; i < inSegments.size(); ++i )
920 InSegment& inSeg = inSegments[i];
922 if ( inSeg._geomEdgeInd != prevGeomEdge )
924 if ( !bndSegs.empty() )
925 bndSegsPerEdge[ prevGeomEdge ].swap( bndSegs );
926 prevGeomEdge = inSeg._geomEdgeInd;
929 // segments around 1st concave point
930 size_t ip0 = inSeg._p0->index( inPoints );
931 if ( inPntChecked[ ip0 ] )
932 for ( e = inSeg._p0->_edges.rbegin(); e != inSeg._p0->_edges.rend(); ++e )
933 bndSegs.push_back( BndSeg( &inSeg, *e, inSeg._p0->_param ));
934 inPntChecked[ ip0 ] = false;
936 // segments of InSegment's
937 const size_t nbMaEdges = inSeg._edges.size();
938 switch ( nbMaEdges ) {
939 case 0: // "around" circle center
940 bndSegs.push_back( BndSeg( &inSeg, 0, inSeg._p1->_param )); break;
942 bndSegs.push_back( BndSeg( &inSeg, inSeg._edges.back(), inSeg._p1->_param )); break;
944 gp_XY inSegDir( inSeg._p1->_a - inSeg._p0->_a,
945 inSeg._p1->_b - inSeg._p0->_b );
946 const double inSegLen2 = inSegDir.SquareModulus();
947 e = inSeg._edges.rbegin();
948 for ( size_t iE = 1; iE < nbMaEdges; ++e, ++iE )
950 gp_XY toMA( (*e)->vertex0()->x() - inSeg._p0->_a,
951 (*e)->vertex0()->y() - inSeg._p0->_b );
952 double r = toMA * inSegDir / inSegLen2;
953 double u = r * inSeg._p1->_param + ( 1. - r ) * inSeg._p0->_param;
954 bndSegs.push_back( BndSeg( &inSeg, *e, u ));
956 bndSegs.push_back( BndSeg( &inSeg, *e, inSeg._p1->_param ));
958 // segments around 2nd concave point
959 size_t ip1 = inSeg._p1->index( inPoints );
960 if ( inPntChecked[ ip1 ] )
961 for ( e = inSeg._p1->_edges.rbegin(); e != inSeg._p1->_edges.rend(); ++e )
962 bndSegs.push_back( BndSeg( &inSeg, *e, inSeg._p1->_param ));
963 inPntChecked[ ip1 ] = false;
965 if ( !bndSegs.empty() )
966 bndSegsPerEdge[ prevGeomEdge ].swap( bndSegs );
969 // prepare to MA branch search
970 for ( size_t iE = 0; iE < bndSegsPerEdge.size(); ++iE )
972 // 1) make TVDEdge's know it's BndSeg to enable passing branchID to
973 // an opposite BndSeg in BndSeg::setBranch(); geom EDGE ID is known from TVDCell
974 // 2) connect bndSegs via BndSeg::_prev
976 vector< BndSeg >& bndSegs = bndSegsPerEdge[ iE ];
977 if ( bndSegs.empty() ) continue;
979 for ( size_t i = 1; i < bndSegs.size(); ++i )
981 bndSegs[i]._prev = & bndSegs[i-1];
982 bndSegs[i].setIndexToEdge( i );
984 // look for the last bndSeg of previous EDGE to set bndSegs[0]._prev
985 const InPoint& p0 = bndSegs[0]._inSeg->point0();
986 for ( size_t iE2 = 0; iE2 < bndSegsPerEdge.size(); ++iE2 )
987 if ( p0 == bndSegsPerEdge[ iE2 ].back()._inSeg->point1() )
989 bndSegs[0]._prev = & bndSegsPerEdge[ iE2 ].back();
992 bndSegs[0].setIndexToEdge( 0 );
995 bndSegsToMesh( bndSegsPerEdge ); // debug: visually check found MA edges
998 // Find TVDEdge's of Branches and associate them with bndSegs
1000 vector< vector<const TVDEdge*> > branchEdges;
1001 branchEdges.reserve( boundary.nbEdges() * 4 );
1003 map< const TVDVertex*, SMESH_MAT2d::BranchEndType > endType;
1005 int branchID = 1; // we code orientation as branchID sign
1006 branchEdges.resize( branchID );
1008 for ( size_t iE = 0; iE < bndSegsPerEdge.size(); ++iE )
1010 vector< BndSeg >& bndSegs = bndSegsPerEdge[ iE ];
1011 for ( size_t i = 0; i < bndSegs.size(); ++i )
1013 if ( bndSegs[i].branchID() )
1015 if ( bndSegs[i]._prev &&
1016 bndSegs[i]._branchID == -bndSegs[i]._prev->_branchID &&
1019 SMESH_MAT2d::BranchEndType type =
1020 ( bndSegs[i]._inSeg->isConnected( bndSegs[i]._edge ) ?
1021 SMESH_MAT2d::BE_ON_VERTEX :
1022 SMESH_MAT2d::BE_END );
1023 endType.insert( make_pair( bndSegs[i]._edge->vertex1(), type ));
1027 if ( !bndSegs[i]._prev &&
1028 !bndSegs[i].hasOppositeEdge() )
1031 if ( !bndSegs[i]._prev ||
1032 !bndSegs[i]._prev->isSameBranch( bndSegs[i] ))
1034 branchEdges.resize(( branchID = branchEdges.size()) + 1 );
1035 if ( bndSegs[i]._edge && bndSegs[i]._prev )
1036 endType.insert( make_pair( bndSegs[i]._edge->vertex1(), SMESH_MAT2d::BE_BRANCH_POINT ));
1038 else if ( bndSegs[i]._prev->_branchID )
1040 branchID = bndSegs[i]._prev->_branchID; // with sign
1042 else if ( bndSegs[i]._edge ) // 1st bndSeg of a WIRE
1044 branchEdges.resize(( branchID = branchEdges.size()) + 1 );
1045 if ( bndSegs[i]._inSeg->isConnected( bndSegs[i]._edge ))
1047 if ( bndSegs[i]._inSeg->point0() == bndSegs[i]._edge->vertex1() )
1048 endType.insert( make_pair( bndSegs[i]._edge->vertex1(), SMESH_MAT2d::BE_ON_VERTEX ));
1050 endType.insert( make_pair( bndSegs[i]._edge->vertex0(), SMESH_MAT2d::BE_ON_VERTEX ));
1054 bndSegs[i].setBranch( branchID, bndSegsPerEdge ); // set to i-th and to the opposite bndSeg
1055 if ( bndSegs[i].hasOppositeEdge() )
1056 branchEdges[ bndSegs[i].branchID() ].push_back( bndSegs[i]._edge );
1060 // join the 1st and the last branch edges if it is the same branch
1061 // if ( bndSegs.back().branchID() != bndSegs.front().branchID() &&
1062 // bndSegs.back().isSameBranch( bndSegs.front() ))
1064 // vector<const TVDEdge*> & br1 = branchEdges[ bndSegs.front().branchID() ];
1065 // vector<const TVDEdge*> & br2 = branchEdges[ bndSegs.back().branchID() ];
1066 // br1.insert( br1.begin(), br2.begin(), br2.end() );
1070 // remove branches ending at BE_ON_VERTEX
1072 vector<bool> isBranchRemoved( branchEdges.size(), false );
1074 if ( ignoreCorners && branchEdges.size() > 2 && !branchEdges[2].empty() )
1076 // find branches to remove
1077 map< const TVDVertex*, SMESH_MAT2d::BranchEndType >::iterator v2et;
1078 for ( size_t iB = 1; iB < branchEdges.size(); ++iB )
1080 if ( branchEdges[iB].empty() )
1082 const TVDVertex* v0 = branchEdges[iB][0]->vertex1();
1083 const TVDVertex* v1 = branchEdges[iB].back()->vertex0();
1084 v2et = endType.find( v0 );
1085 if ( v2et != endType.end() && v2et->second == SMESH_MAT2d::BE_ON_VERTEX )
1086 isBranchRemoved[ iB ] = true;
1087 v2et = endType.find( v1 );
1088 if ( v2et != endType.end() && v2et->second == SMESH_MAT2d::BE_ON_VERTEX )
1089 isBranchRemoved[ iB ] = true;
1091 // try to join not removed branches into one
1092 for ( size_t iB = 1; iB < branchEdges.size(); ++iB )
1094 if ( branchEdges[iB].empty() || isBranchRemoved[iB] )
1096 const TVDVertex* v0 = branchEdges[iB][0]->vertex1();
1097 const TVDVertex* v1 = branchEdges[iB].back()->vertex0();
1098 v2et = endType.find( v0 );
1099 if ( v2et == endType.end() || v2et->second != SMESH_MAT2d::BE_BRANCH_POINT )
1101 v2et = endType.find( v1 );
1102 if ( v2et == endType.end() || v2et->second != SMESH_MAT2d::BE_BRANCH_POINT )
1107 for ( int isV0 = 0; isV0 < 2; ++isV0 )
1109 const TVDVertex* v = isV0 ? v0 : v1;
1110 size_t iBrToJoin = 0;
1111 for ( size_t iB2 = 1; iB2 < branchEdges.size(); ++iB2 )
1113 if ( branchEdges[iB2].empty() || isBranchRemoved[iB2] || iB == iB2 )
1115 const TVDVertex* v02 = branchEdges[iB2][0]->vertex1();
1116 const TVDVertex* v12 = branchEdges[iB2].back()->vertex0();
1117 if ( v == v02 || v == v12 )
1119 if ( iBrToJoin > 0 )
1122 break; // more than 2 not removed branches meat at a TVDVertex
1127 if ( iBrToJoin > 0 )
1129 vector<const TVDEdge*>& branch = branchEdges[ iBrToJoin ];
1130 const TVDVertex* v02 = branch[0]->vertex1();
1131 const TVDVertex* v12 = branch.back()->vertex0();
1132 updateJoinedBranch( branch, iB, bndSegsPerEdge, /*reverse=*/(v0 == v02 || v1 == v12 ));
1133 if ( v0 == v02 || v0 == v12 )
1134 branchEdges[iB].insert( branchEdges[iB].begin(), branch.begin(), branch.end() );
1136 branchEdges[iB].insert( branchEdges[iB].end(), branch.begin(), branch.end() );
1140 } // loop on branchEdges
1141 } // if ( ignoreCorners )
1143 // associate branchIDs and the input branch vector (arg)
1144 vector< SMESH_MAT2d::Branch* > branchByID( branchEdges.size(), 0 );
1146 for ( size_t i = 0; i < branchEdges.size(); ++i )
1148 nbBranches += ( !branchEdges[i].empty() );
1150 branch.resize( nbBranches );
1152 for ( size_t brID = 1; brID < branchEdges.size(); ++brID ) // 1st - not removed
1154 if ( !branchEdges[ brID ].empty() && !isBranchRemoved[ brID ])
1155 branchByID[ brID ] = & branch[ iBr++ ];
1157 for ( size_t brID = 1; brID < branchEdges.size(); ++brID ) // then - removed
1159 if ( !branchEdges[ brID ].empty() && isBranchRemoved[ brID ])
1160 branchByID[ brID ] = & branch[ iBr++ ];
1163 // Fill in BndPoints of each EDGE of the boundary
1166 int edgeInd = -1, dInd = 0;
1167 for ( size_t iE = 0; iE < bndSegsPerEdge.size(); ++iE )
1169 vector< BndSeg >& bndSegs = bndSegsPerEdge[ iE ];
1170 SMESH_MAT2d::BndPoints & bndPoints = boundary.getPoints( iE );
1172 // make TVDEdge know an index of bndSegs within BndPoints
1173 for ( size_t i = 0; i < bndSegs.size(); ++i )
1174 if ( bndSegs[i]._edge )
1175 SMESH_MAT2d::Branch::setBndSegment( i, bndSegs[i]._edge );
1177 // parameters on EDGE
1179 bndPoints._params.reserve( bndSegs.size() + 1 );
1180 bndPoints._params.push_back( bndSegs[ 0 ]._inSeg->_p0->_param );
1182 for ( size_t i = 0; i < bndSegs.size(); ++i )
1183 bndPoints._params.push_back( bndSegs[ i ]._uLast );
1187 bndPoints._maEdges.reserve( bndSegs.size() );
1189 for ( size_t i = 0; i < bndSegs.size(); ++i )
1191 const size_t brID = bndSegs[ i ].branchID();
1192 const SMESH_MAT2d::Branch* br = branchByID[ brID ];
1194 if ( bndSegs[ i ]._edge && !branchEdges[ brID ].empty() )
1198 if (( edgeInd < 0 ||
1199 edgeInd >= (int) branchEdges[ brID ].size() ) ||
1200 ( branchEdges[ brID ][ edgeInd ] != bndSegs[ i ]._edge &&
1201 branchEdges[ brID ][ edgeInd ]->twin() != bndSegs[ i ]._edge ))
1203 if ( bndSegs[ i ]._branchID < 0 )
1206 for ( edgeInd = branchEdges[ brID ].size() - 1; edgeInd > 0; --edgeInd )
1207 if ( branchEdges[ brID ][ edgeInd ]->twin() == bndSegs[ i ]._edge )
1210 else // bndSegs[ i ]._branchID > 0
1213 for ( edgeInd = 0; edgeInd < (int)branchEdges[ brID ].size(); ++edgeInd )
1214 if ( branchEdges[ brID ][ edgeInd ] == bndSegs[ i ]._edge )
1221 // no MA edge, bndSeg corresponds to an end point of a branch
1222 if ( bndPoints._maEdges.empty() )
1225 edgeInd = branchEdges[ brID ].size();
1226 dInd = bndSegs[ i ]._branchID > 0 ? +1 : -1;
1229 bndPoints._maEdges.push_back( make_pair( br, ( 1 + edgeInd ) * dInd ));
1231 } // loop on bndSegs of an EDGE
1232 } // loop on all bndSegs to construct Boundary
1234 // Initialize branches
1236 // find a not removed branch
1237 size_t iBrNorRemoved = 0;
1238 for ( size_t brID = 1; brID < branchEdges.size(); ++brID )
1239 if ( !branchEdges[brID].empty() && !isBranchRemoved[brID] )
1241 iBrNorRemoved = brID;
1244 // fill the branches with MA edges
1245 for ( size_t brID = 1; brID < branchEdges.size(); ++brID )
1246 if ( !branchEdges[brID].empty() )
1248 branchByID[ brID ]->init( branchEdges[brID], & boundary, endType );
1250 // mark removed branches
1251 for ( size_t brID = 1; brID < branchEdges.size(); ++brID )
1252 if ( isBranchRemoved[brID] && iBrNorRemoved > 0 )
1254 SMESH_MAT2d::Branch* branch = branchByID[ brID ];
1255 SMESH_MAT2d::Branch* mainBranch = branchByID[ iBrNorRemoved ];
1256 bool is1stBrPnt = ( branch->getEnd(0)->_type == SMESH_MAT2d::BE_BRANCH_POINT );
1257 const TVDVertex* branchVextex =
1258 is1stBrPnt ? branch->getEnd(0)->_vertex : branch->getEnd(1)->_vertex;
1259 SMESH_MAT2d::BranchPoint bp = mainBranch->getPoint( branchVextex );
1260 branch->setRemoved( bp );
1262 // set branches to branch ends
1263 for ( size_t i = 0; i < branch.size(); ++i )
1264 if ( !branch[i].isRemoved() )
1265 branch[i].setBranchesToEnds( branch );
1267 // fill branchPnt arg
1268 map< const TVDVertex*, const SMESH_MAT2d::BranchEnd* > v2end;
1269 for ( size_t i = 0; i < branch.size(); ++i )
1271 if ( branch[i].getEnd(0)->_branches.size() > 2 )
1272 v2end.insert( make_pair( branch[i].getEnd(0)->_vertex, branch[i].getEnd(0) ));
1273 if ( branch[i].getEnd(1)->_branches.size() > 2 )
1274 v2end.insert( make_pair( branch[i].getEnd(1)->_vertex, branch[i].getEnd(1) ));
1276 branchPnt.resize( v2end.size() );
1277 map< const TVDVertex*, const SMESH_MAT2d::BranchEnd* >::iterator v2e = v2end.begin();
1278 for ( size_t i = 0; v2e != v2end.end(); ++v2e, ++i )
1279 branchPnt[ i ] = v2e->second;
1285 //================================================================================
1287 * \brief MedialAxis constructor
1288 * \param [in] face - a face to create MA for
1289 * \param [in] edges - edges of the face (possibly not all) on the order they
1290 * encounter in the face boundary.
1291 * \param [in] minSegLen - minimal length of a mesh segment used to discretize
1292 * the edges. It is used to define precision of MA approximation
1294 //================================================================================
1296 SMESH_MAT2d::MedialAxis::MedialAxis(const TopoDS_Face& face,
1297 const std::vector< TopoDS_Edge >& edges,
1298 const double minSegLen,
1299 const bool ignoreCorners):
1300 _face( face ), _boundary( edges.size() )
1302 // input to construct_voronoi()
1303 vector< InPoint > inPoints;
1304 vector< InSegment> inSegments;
1305 if ( !makeInputData( face, edges, minSegLen, inPoints, inSegments, _scale ))
1308 inSegmentsToFile( inSegments );
1310 // build voronoi diagram
1311 construct_voronoi( inSegments.begin(), inSegments.end(), &_vd );
1314 makeMA( _vd, ignoreCorners, inPoints, inSegments, _branch, _branchPnt, _boundary );
1316 // count valid branches
1317 _nbBranches = _branch.size();
1318 for ( size_t i = 0; i < _branch.size(); ++i )
1319 if ( _branch[i].isRemoved() )
1323 //================================================================================
1325 * \brief Returns the i-th branch
1327 //================================================================================
1329 const SMESH_MAT2d::Branch* SMESH_MAT2d::MedialAxis::getBranch(size_t i) const
1331 return i < _nbBranches ? &_branch[i] : 0;
1334 //================================================================================
1336 * \brief Return UVs of ends of MA edges of a branch
1338 //================================================================================
1340 void SMESH_MAT2d::MedialAxis::getPoints( const Branch* branch,
1341 std::vector< gp_XY >& points) const
1343 branch->getPoints( points, _scale );
1346 //================================================================================
1348 * \brief Returns a BranchPoint corresponding to a given point on a geom EDGE
1349 * \param [in] iEdge - index of geom EDGE within a vector passed at MA construction
1350 * \param [in] u - parameter of the point on EDGE curve
1351 * \param [out] p - the found BranchPoint
1352 * \return bool - is OK
1354 //================================================================================
1356 bool SMESH_MAT2d::Boundary::getBranchPoint( const std::size_t iEdge,
1358 BranchPoint& p ) const
1360 if ( iEdge >= _pointsPerEdge.size() || _pointsPerEdge[iEdge]._params.empty() )
1363 const BndPoints& points = _pointsPerEdge[ iEdge ];
1364 const bool edgeReverse = ( points._params[0] > points._params.back() );
1366 if ( u < ( edgeReverse ? points._params.back() : points._params[0] ))
1367 u = edgeReverse ? points._params.back() : points._params[0];
1368 else if ( u > ( edgeReverse ? points._params[0] : points._params.back()) )
1369 u = edgeReverse ? points._params[0] : points._params.back();
1371 double r = ( u - points._params[0] ) / ( points._params.back() - points._params[0] );
1372 int i = int( r * double( points._maEdges.size()-1 ));
1375 while ( points._params[i ] < u ) --i;
1376 while ( points._params[i+1] > u ) ++i;
1380 while ( points._params[i ] > u ) --i;
1381 while ( points._params[i+1] < u ) ++i;
1384 if ( points._params[i] == points._params[i+1] ) // coincident points at some end
1386 int di = ( points._params[0] == points._params[i] ) ? +1 : -1;
1387 while ( points._params[i] == points._params[i+1] )
1389 if ( i < 0 || i+1 >= (int)points._params.size() )
1393 double edgeParam = ( u - points._params[i] ) / ( points._params[i+1] - points._params[i] );
1395 if ( !points._maEdges[ i ].second ) // no branch at the EDGE end, look for a closest branch
1397 if ( i < (int)points._maEdges.size() / 2 ) // near 1st point
1399 while ( i < (int)points._maEdges.size()-1 && !points._maEdges[ i ].second )
1401 edgeParam = edgeReverse;
1403 else // near last point
1405 while ( i > 0 && !points._maEdges[ i ].second )
1407 edgeParam = !edgeReverse;
1410 const std::pair< const Branch*, int >& maE = points._maEdges[ i ];
1411 bool maReverse = ( maE.second < 0 );
1413 p._branch = maE.first;
1414 p._iEdge = ( maReverse ? -maE.second : maE.second ) - 1; // countered from 1 to store sign
1415 p._edgeParam = ( maE.first && maReverse ) ? ( 1. - edgeParam ) : edgeParam;
1420 //================================================================================
1422 * \brief Returns a BranchPoint corresponding to a given BoundaryPoint on a geom EDGE
1423 * \param [in] bp - the BoundaryPoint
1424 * \param [out] p - the found BranchPoint
1425 * \return bool - is OK
1427 //================================================================================
1429 bool SMESH_MAT2d::Boundary::getBranchPoint( const BoundaryPoint& bp,
1430 BranchPoint& p ) const
1432 return getBranchPoint( bp._edgeIndex, bp._param, p );
1435 //================================================================================
1437 * \brief Check if a given boundary segment is a null-length segment on a concave
1439 * \param [in] iEdge - index of a geom EDGE
1440 * \param [in] iSeg - index of a boundary segment
1441 * \return bool - true if the segment is on concave corner
1443 //================================================================================
1445 bool SMESH_MAT2d::Boundary::isConcaveSegment( std::size_t iEdge, std::size_t iSeg ) const
1447 if ( iEdge >= _pointsPerEdge.size() || _pointsPerEdge[iEdge]._params.empty() )
1450 const BndPoints& points = _pointsPerEdge[ iEdge ];
1451 if ( points._params.size() <= iSeg+1 )
1454 return Abs( points._params[ iSeg ] - points._params[ iSeg+1 ]) < 1e-20;
1457 //================================================================================
1459 * \brief Moves (changes _param) a given BoundaryPoint to a closest EDGE end
1461 //================================================================================
1463 bool SMESH_MAT2d::Boundary::moveToClosestEdgeEnd( BoundaryPoint& bp ) const
1465 if ( bp._edgeIndex >= _pointsPerEdge.size() )
1468 const BndPoints& points = _pointsPerEdge[ bp._edgeIndex ];
1469 if ( Abs( bp._param - points._params[0]) < Abs( points._params.back() - bp._param ))
1470 bp._param = points._params[0];
1472 bp._param = points._params.back();
1477 //================================================================================
1479 * \brief Creates a 3d curve corresponding to a Branch
1480 * \param [in] branch - the Branch
1481 * \return Adaptor3d_Curve* - the new curve the caller is to delete
1483 //================================================================================
1485 Adaptor3d_Curve* SMESH_MAT2d::MedialAxis::make3DCurve(const Branch& branch) const
1487 Handle(Geom_Surface) surface = BRep_Tool::Surface( _face );
1488 if ( surface.IsNull() )
1492 branch.getPoints( uv, _scale );
1493 if ( uv.size() < 2 )
1496 vector< TopoDS_Vertex > vertex( uv.size() );
1497 for ( size_t i = 0; i < uv.size(); ++i )
1498 vertex[i] = BRepBuilderAPI_MakeVertex( surface->Value( uv[i].X(), uv[i].Y() ));
1501 BRep_Builder aBuilder;
1502 aBuilder.MakeWire(aWire);
1503 for ( size_t i = 1; i < vertex.size(); ++i )
1505 TopoDS_Edge edge = BRepBuilderAPI_MakeEdge( vertex[i-1], vertex[i] );
1506 aBuilder.Add( aWire, edge );
1509 // if ( myEdge.size() == 2 && FirstVertex().IsSame( LastVertex() ))
1510 // aWire.Closed(true); // issue 0021141
1512 return new BRepAdaptor_CompCurve( aWire );
1515 //================================================================================
1517 * \brief Copy points of an EDGE
1519 //================================================================================
1521 void SMESH_MAT2d::Branch::init( vector<const TVDEdge*>& maEdges,
1522 const Boundary* boundary,
1523 map< const TVDVertex*, BranchEndType > endType )
1525 if ( maEdges.empty() ) return;
1527 _boundary = boundary;
1528 _maEdges.swap( maEdges );
1531 _params.reserve( _maEdges.size() + 1 );
1532 _params.push_back( 0. );
1533 for ( size_t i = 0; i < _maEdges.size(); ++i )
1534 _params.push_back( _params.back() + length( _maEdges[i] ));
1536 for ( size_t i = 1; i < _params.size(); ++i )
1537 _params[i] /= _params.back();
1540 _endPoint1._vertex = _maEdges.front()->vertex1();
1541 _endPoint2._vertex = _maEdges.back ()->vertex0();
1543 if ( endType.count( _endPoint1._vertex ))
1544 _endPoint1._type = endType[ _endPoint1._vertex ];
1545 if ( endType.count( _endPoint2._vertex ))
1546 _endPoint2._type = endType[ _endPoint2._vertex ];
1549 //================================================================================
1551 * \brief fills BranchEnd::_branches of its ends
1553 //================================================================================
1555 void SMESH_MAT2d::Branch::setBranchesToEnds( const vector< Branch >& branches )
1557 for ( size_t i = 0; i < branches.size(); ++i )
1559 if ( this->_endPoint1._vertex == branches[i]._endPoint1._vertex ||
1560 this->_endPoint1._vertex == branches[i]._endPoint2._vertex )
1561 this->_endPoint1._branches.push_back( &branches[i] );
1563 if ( this->_endPoint2._vertex == branches[i]._endPoint1._vertex ||
1564 this->_endPoint2._vertex == branches[i]._endPoint2._vertex )
1565 this->_endPoint2._branches.push_back( &branches[i] );
1569 //================================================================================
1571 * \brief returns a BranchPoint corresponding to a TVDVertex
1573 //================================================================================
1575 SMESH_MAT2d::BranchPoint SMESH_MAT2d::Branch::getPoint( const TVDVertex* vertex ) const
1581 if ( vertex == _maEdges[0]->vertex1() )
1587 for ( ; p._iEdge < _maEdges.size(); ++p._iEdge )
1588 if ( vertex == _maEdges[ p._iEdge ]->vertex0() )
1590 p._edgeParam = _params[ p._iEdge ];
1597 //================================================================================
1599 * \brief Sets a proxy point for a removed branch
1600 * \param [in] proxyPoint - a point of another branch to which all points of this
1603 //================================================================================
1605 void SMESH_MAT2d::Branch::setRemoved( const BranchPoint& proxyPoint )
1607 _proxyPoint = proxyPoint;
1610 //================================================================================
1612 * \brief Returns points on two EDGEs, equidistant from a given point of this Branch
1613 * \param [in] param - [0;1] normalized param on the Branch
1614 * \param [out] bp1 - BoundaryPoint on EDGE with a lower index
1615 * \param [out] bp2 - BoundaryPoint on EDGE with a higher index
1616 * \return bool - true if the BoundaryPoint's found
1618 //================================================================================
1620 bool SMESH_MAT2d::Branch::getBoundaryPoints(double param,
1622 BoundaryPoint& bp2 ) const
1624 if ( param < _params[0] || param > _params.back() )
1627 // look for an index of a MA edge by param
1628 double ip = param * _params.size();
1629 size_t i = size_t( Min( int( _maEdges.size()-1), int( ip )));
1631 while ( param < _params[i ] ) --i;
1632 while ( param > _params[i+1] ) ++i;
1634 double r = ( param - _params[i] ) / ( _params[i+1] - _params[i] );
1636 return getBoundaryPoints( i, r, bp1, bp2 );
1639 //================================================================================
1641 * \brief Returns points on two EDGEs, equidistant from a given point of this Branch
1642 * \param [in] iMAEdge - index of a MA edge within this Branch
1643 * \param [in] maEdgeParam - [0;1] normalized param on the \a iMAEdge
1644 * \param [out] bp1 - BoundaryPoint on EDGE with a lower index
1645 * \param [out] bp2 - BoundaryPoint on EDGE with a higher index
1646 * \return bool - true if the BoundaryPoint's found
1648 //================================================================================
1650 bool SMESH_MAT2d::Branch::getBoundaryPoints(std::size_t iMAEdge,
1653 BoundaryPoint& bp2 ) const
1656 return _proxyPoint._branch->getBoundaryPoints( _proxyPoint, bp1, bp2 );
1658 if ( iMAEdge > _maEdges.size() )
1660 if ( iMAEdge == _maEdges.size() )
1661 iMAEdge = _maEdges.size() - 1;
1663 size_t iGeom1 = getGeomEdge( _maEdges[ iMAEdge ] );
1664 size_t iGeom2 = getGeomEdge( _maEdges[ iMAEdge ]->twin() );
1665 size_t iSeg1 = getBndSegment( _maEdges[ iMAEdge ] );
1666 size_t iSeg2 = getBndSegment( _maEdges[ iMAEdge ]->twin() );
1668 return ( _boundary->getPoint( iGeom1, iSeg1, maEdgeParam, bp1 ) &&
1669 _boundary->getPoint( iGeom2, iSeg2, maEdgeParam, bp2 ));
1672 //================================================================================
1674 * \brief Returns points on two EDGEs, equidistant from a given point of this Branch
1676 //================================================================================
1678 bool SMESH_MAT2d::Branch::getBoundaryPoints(const BranchPoint& p,
1680 BoundaryPoint& bp2 ) const
1682 return ( p._branch ? p._branch : this )->getBoundaryPoints( p._iEdge, p._edgeParam, bp1, bp2 );
1685 //================================================================================
1687 * \brief Return a parameter of a BranchPoint normalized within this Branch
1689 //================================================================================
1691 bool SMESH_MAT2d::Branch::getParameter(const BranchPoint & p, double & u ) const
1693 if ( this != p._branch && p._branch )
1694 return p._branch->getParameter( p, u );
1697 return _proxyPoint._branch->getParameter( _proxyPoint, u );
1699 if ( p._iEdge > _params.size()-1 )
1701 if ( p._iEdge == _params.size()-1 )
1704 u = ( _params[ p._iEdge ] * ( 1 - p._edgeParam ) +
1705 _params[ p._iEdge+1 ] * p._edgeParam );
1710 //================================================================================
1712 * \brief Check type of both ends
1714 //================================================================================
1716 bool SMESH_MAT2d::Branch::hasEndOfType(BranchEndType type) const
1718 return ( _endPoint1._type == type || _endPoint2._type == type );
1721 //================================================================================
1723 * \brief Returns MA points
1724 * \param [out] points - the 2d points
1725 * \param [in] scale - the scale that was used to scale the 2d space of MA
1727 //================================================================================
1729 void SMESH_MAT2d::Branch::getPoints( std::vector< gp_XY >& points,
1730 const double scale[2]) const
1732 points.resize( _maEdges.size() + 1 );
1734 points[0].SetCoord( _maEdges[0]->vertex1()->x() / scale[0], // CCW order! -> vertex1 not vertex0
1735 _maEdges[0]->vertex1()->y() / scale[1] );
1737 for ( size_t i = 0; i < _maEdges.size(); ++i )
1738 points[i+1].SetCoord( _maEdges[i]->vertex0()->x() / scale[0],
1739 _maEdges[i]->vertex0()->y() / scale[1] );
1742 //================================================================================
1744 * \brief Return indices of EDGEs equidistant from this branch
1746 //================================================================================
1748 void SMESH_MAT2d::Branch::getGeomEdges( std::vector< std::size_t >& edgeIDs1,
1749 std::vector< std::size_t >& edgeIDs2 ) const
1751 edgeIDs1.push_back( getGeomEdge( _maEdges[0] ));
1752 edgeIDs2.push_back( getGeomEdge( _maEdges[0]->twin() ));
1754 for ( size_t i = 1; i < _maEdges.size(); ++i )
1756 size_t ie1 = getGeomEdge( _maEdges[i] );
1757 size_t ie2 = getGeomEdge( _maEdges[i]->twin() );
1759 if ( edgeIDs1.back() != ie1 ) edgeIDs1.push_back( ie1 );
1760 if ( edgeIDs2.back() != ie2 ) edgeIDs2.push_back( ie2 );
1764 //================================================================================
1766 * \brief Looks for a BranchPoint position around a concave VERTEX
1768 //================================================================================
1770 bool SMESH_MAT2d::Branch::addDivPntForConcaVertex( std::vector< std::size_t >& edgeIDs1,
1771 std::vector< std::size_t >& edgeIDs2,
1772 std::vector< BranchPoint >& divPoints,
1773 const vector<const TVDEdge*>& maEdges,
1774 const vector<const TVDEdge*>& maEdgesTwin,
1777 // if there is a concave vertex between EDGEs
1778 // then position of a dividing BranchPoint is undefined, it is somewhere
1779 // on an arc-shaped part of the Branch around the concave vertex.
1780 // Chose this position by a VERTEX of the opposite EDGE, or put it in the middle
1781 // of the arc if there is no opposite VERTEX.
1782 // All null-length segments around a VERTEX belong to one of EDGEs.
1784 BranchPoint divisionPnt;
1785 divisionPnt._branch = this;
1787 BranchIterator iCur( maEdges, i );
1789 size_t ie1 = getGeomEdge( maEdges [i] );
1790 size_t ie2 = getGeomEdge( maEdgesTwin[i] );
1792 size_t iSeg1 = getBndSegment( iCur.edgePrev() );
1793 size_t iSeg2 = getBndSegment( iCur.edge() );
1794 bool isConcaPrev = _boundary->isConcaveSegment( edgeIDs1.back(), iSeg1 );
1795 bool isConcaNext = _boundary->isConcaveSegment( ie1, iSeg2 );
1796 if ( !isConcaNext && !isConcaPrev )
1799 bool isConcaveV = false;
1802 BranchIterator iPrev( maEdges, i ), iNext( maEdges, i );
1804 if ( isConcaNext ) // all null-length segments follow
1806 // look for a VERTEX of the opposite EDGE
1807 // iNext - next after all null-length segments
1808 while (( maE = ++iNext ))
1810 iSeg2 = getBndSegment( maE );
1811 if ( !_boundary->isConcaveSegment( ie1, iSeg2 ))
1814 bool vertexFound = false;
1815 for ( ++iCur; iCur < iNext; ++iCur )
1817 ie2 = getGeomEdge( maEdgesTwin[ iCur.indexMod() ] );
1818 if ( ie2 != edgeIDs2.back() )
1820 // opposite VERTEX found
1821 divisionPnt._iEdge = iCur.indexMod();
1822 divisionPnt._edgeParam = 0;
1823 divPoints.push_back( divisionPnt );
1824 edgeIDs1.push_back( ie1 );
1825 edgeIDs2.push_back( ie2 );
1832 iPrev = iNext; // not to add a BP in the moddle
1833 i = iNext.indexMod();
1837 else if ( isConcaPrev )
1839 // all null-length segments passed, find their beginning
1840 while (( maE = iPrev.edgePrev() ))
1842 iSeg1 = getBndSegment( maE );
1843 if ( _boundary->isConcaveSegment( edgeIDs1.back(), iSeg1 ))
1850 if ( iPrev.index() < i-1 || iNext.index() > i )
1852 // no VERTEX on the opposite EDGE, put the Branch Point in the middle
1853 divisionPnt._iEdge = iPrev.indexMod();
1855 double par1 = _params[ iPrev.indexMod() ], par2 = _params[ iNext.indexMod() ];
1856 double midPar = 0.5 * ( par1 + par2 );
1857 for ( ; _params[ iPrev.indexMod() ] < midPar; ++iPrev )
1858 divisionPnt._iEdge = iPrev.indexMod();
1859 divisionPnt._edgeParam =
1860 ( _params[ iPrev.indexMod() ] - midPar ) /
1861 ( _params[ iPrev.indexMod() ] - _params[ divisionPnt._iEdge ] );
1862 divPoints.push_back( divisionPnt );
1869 //================================================================================
1871 * \brief Return indices of opposite parts of EDGEs equidistant from this branch
1872 * \param [out] edgeIDs1 - EDGE index opposite to the edgeIDs2[i]-th EDGE
1873 * \param [out] edgeIDs2 - EDGE index opposite to the edgeIDs1[i]-th EDGE
1874 * \param [out] divPoints - BranchPoint's located between two successive unique
1875 * pairs of EDGE indices. A \a divPoints[i] can separate e.g. two following pairs
1876 * of EDGE indices < 0, 2 > and < 0, 1 >. Number of \a divPoints is one less
1877 * than number of \a edgeIDs
1879 //================================================================================
1881 void SMESH_MAT2d::Branch::getOppositeGeomEdges( std::vector< std::size_t >& edgeIDs1,
1882 std::vector< std::size_t >& edgeIDs2,
1883 std::vector< BranchPoint >& divPoints) const
1889 std::vector<const TVDEdge*> twins( _maEdges.size() );
1890 for ( size_t i = 0; i < _maEdges.size(); ++i )
1891 twins[i] = _maEdges[i]->twin();
1893 BranchIterator maIter ( _maEdges, 0 );
1894 BranchIterator twIter ( twins, 0 );
1895 // size_t lastConcaE1 = _boundary.nbEdges();
1896 // size_t lastConcaE2 = _boundary.nbEdges();
1898 // if ( maIter._closed ) // closed branch
1900 // edgeIDs1.push_back( getGeomEdge( _maEdges.back() ));
1901 // edgeIDs2.push_back( getGeomEdge( _maEdges.back()->twin() ));
1905 edgeIDs1.push_back( getGeomEdge( maIter.edge() ));
1906 edgeIDs2.push_back( getGeomEdge( twIter.edge() ));
1909 BranchPoint divisionPnt;
1910 divisionPnt._branch = this;
1912 for ( ++maIter, ++twIter; maIter.index() < (int)_maEdges.size(); ++maIter, ++twIter )
1914 size_t ie1 = getGeomEdge( maIter.edge() );
1915 size_t ie2 = getGeomEdge( twIter.edge() );
1917 bool otherE1 = ( edgeIDs1.back() != ie1 );
1918 bool otherE2 = ( edgeIDs2.back() != ie2 );
1920 if ( !otherE1 && !otherE2 && maIter._closed )
1922 int iSegPrev1 = getBndSegment( maIter.edgePrev() );
1923 int iSegCur1 = getBndSegment( maIter.edge() );
1924 otherE1 = Abs( iSegPrev1 - iSegCur1 ) != 1;
1925 int iSegPrev2 = getBndSegment( twIter.edgePrev() );
1926 int iSegCur2 = getBndSegment( twIter.edge() );
1927 otherE2 = Abs( iSegPrev2 - iSegCur2 ) != 1;
1930 if ( otherE1 || otherE2 )
1932 bool isConcaveV = false;
1933 if ( otherE1 && !otherE2 )
1935 isConcaveV = addDivPntForConcaVertex( edgeIDs1, edgeIDs2, divPoints,
1936 _maEdges, twins, maIter._i );
1938 if ( !otherE1 && otherE2 )
1940 isConcaveV = addDivPntForConcaVertex( edgeIDs2, edgeIDs1, divPoints,
1941 twins, _maEdges, maIter._i );
1946 ie1 = getGeomEdge( maIter.edge() );
1947 ie2 = getGeomEdge( twIter.edge() );
1949 if ( !isConcaveV || otherE1 || otherE2 )
1951 edgeIDs1.push_back( ie1 );
1952 edgeIDs2.push_back( ie2 );
1954 if ( divPoints.size() < edgeIDs1.size() - 1 )
1956 divisionPnt._iEdge = maIter.index();
1957 divisionPnt._edgeParam = 0;
1958 divPoints.push_back( divisionPnt );
1961 } // if ( edgeIDs1.back() != ie1 || edgeIDs2.back() != ie2 )
1962 } // loop on _maEdges
1965 //================================================================================
1967 * \brief Store data of boundary segments in TVDEdge
1969 //================================================================================
1971 void SMESH_MAT2d::Branch::setGeomEdge( std::size_t geomIndex, const TVDEdge* maEdge )
1973 if ( maEdge ) maEdge->cell()->color( geomIndex );
1975 std::size_t SMESH_MAT2d::Branch::getGeomEdge( const TVDEdge* maEdge )
1977 return maEdge ? maEdge->cell()->color() : std::string::npos;
1979 void SMESH_MAT2d::Branch::setBndSegment( std::size_t segIndex, const TVDEdge* maEdge )
1981 if ( maEdge ) maEdge->color( segIndex );
1983 std::size_t SMESH_MAT2d::Branch::getBndSegment( const TVDEdge* maEdge )
1985 return maEdge ? maEdge->color() : std::string::npos;
1988 //================================================================================
1990 * \brief Returns a boundary point on a given EDGE
1991 * \param [in] iEdge - index of the EDGE within MedialAxis
1992 * \param [in] iSeg - index of a boundary segment within this Branch
1993 * \param [in] u - [0;1] normalized param within \a iSeg-th segment
1994 * \param [out] bp - the found BoundaryPoint
1995 * \return bool - true if the BoundaryPoint is found
1997 //================================================================================
1999 bool SMESH_MAT2d::Boundary::getPoint( std::size_t iEdge,
2002 BoundaryPoint& bp ) const
2004 if ( iEdge >= _pointsPerEdge.size() )
2006 if ( iSeg+1 >= _pointsPerEdge[ iEdge ]._params.size() )
2009 // This method is called by Branch that can have an opposite orientation,
2010 // hence u is inverted depending on orientation coded as a sign of _maEdge index
2011 bool isReverse = ( _pointsPerEdge[ iEdge ]._maEdges[ iSeg ].second < 0 );
2015 double p0 = _pointsPerEdge[ iEdge ]._params[ iSeg ];
2016 double p1 = _pointsPerEdge[ iEdge ]._params[ iSeg+1 ];
2018 bp._param = p0 * ( 1. - u ) + p1 * u;
2019 bp._edgeIndex = iEdge;