Salome HOME
Fix regressions caused by improvements
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.cxx
1 // Copyright (C) 2007-2015  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
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.
10 //
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.
15 //
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
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 // File      : StdMeshers_FaceSide.hxx
24 // Created   : Wed Jan 31 18:41:25 2007
25 // Author    : Edward AGAPOV (eap)
26 // Module    : SMESH
27 //
28 #include "StdMeshers_FaceSide.hxx"
29
30 #include "SMDS_EdgePosition.hxx"
31 #include "SMDS_MeshNode.hxx"
32 #include "SMESHDS_Mesh.hxx"
33 #include "SMESHDS_SubMesh.hxx"
34 #include "SMESH_Algo.hxx"
35 #include "SMESH_Mesh.hxx"
36 #include "SMESH_MesherHelper.hxx"
37 #include "SMESH_ComputeError.hxx"
38 #include "SMESH_Block.hxx"
39
40 #include <Adaptor2d_Curve2d.hxx>
41 #include <BRepAdaptor_CompCurve.hxx>
42 #include <BRep_Builder.hxx>
43 #include <BRep_Tool.hxx>
44 #include <GCPnts_AbscissaPoint.hxx>
45 #include <Geom2dAdaptor_Curve.hxx>
46 #include <Geom_Line.hxx>
47 #include <TopExp.hxx>
48 #include <TopExp_Explorer.hxx>
49 #include <TopoDS.hxx>
50 #include <TopoDS_Face.hxx>
51 #include <TopoDS_Vertex.hxx>
52 #include <TopoDS_Wire.hxx>
53
54 #include <map>
55 #include <limits>
56
57 #include "utilities.h"
58
59 //================================================================================
60 /*!
61  * \brief Constructor of a side of one edge
62   * \param theFace - the face
63   * \param theEdge - the edge
64  */
65 //================================================================================
66
67 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
68                                          const TopoDS_Edge&   theEdge,
69                                          SMESH_Mesh*          theMesh,
70                                          const bool           theIsForward,
71                                          const bool           theIgnoreMediumNodes,
72                                          SMESH_ProxyMesh::Ptr theProxyMesh)
73 {
74   list<TopoDS_Edge> edges(1,theEdge);
75   *this = StdMeshers_FaceSide( theFace, edges, theMesh, theIsForward,
76                                theIgnoreMediumNodes, theProxyMesh );
77 }
78
79 //================================================================================
80 /*!
81  * \brief Constructor of a side of several edges
82  */
83 //================================================================================
84
85 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face&   theFace,
86                                          list<TopoDS_Edge>&   theEdges,
87                                          SMESH_Mesh*          theMesh,
88                                          const bool           theIsForward,
89                                          const bool           theIgnoreMediumNodes,
90                                          SMESH_ProxyMesh::Ptr theProxyMesh)
91 {
92   int nbEdges = theEdges.size();
93   myEdge.resize      ( nbEdges );
94   myEdgeID.resize    ( nbEdges );
95   myC2d.resize       ( nbEdges );
96   myC3dAdaptor.resize( nbEdges );
97   myFirst.resize     ( nbEdges );
98   myLast.resize      ( nbEdges );
99   myNormPar.resize   ( nbEdges );
100   myEdgeLength.resize( nbEdges );
101   myIsUniform.resize ( nbEdges, true );
102   myFace               = theFace;
103   myLength             = 0;
104   myNbPonits           = myNbSegments = 0;
105   myProxyMesh          = theProxyMesh;
106   myMissingVertexNodes = false;
107   myIgnoreMediumNodes  = theIgnoreMediumNodes;
108   myDefaultPnt2d       = gp_Pnt2d( 1e+100, 1e+100 );
109   if ( !myProxyMesh ) myProxyMesh.reset( new SMESH_ProxyMesh( *theMesh ));
110   if ( nbEdges == 0 ) return;
111
112   SMESHDS_Mesh* meshDS = myProxyMesh->GetMeshDS();
113
114   int nbDegen = 0;
115   list<TopoDS_Edge>::iterator edge = theEdges.begin();
116   for ( int index = 0; edge != theEdges.end(); ++index, ++edge )
117   {
118     int i = theIsForward ? index : nbEdges-index-1;
119     myEdgeLength[i] = SMESH_Algo::EdgeLength( *edge );
120     if ( myEdgeLength[i] < DBL_MIN ) nbDegen++;
121     myLength += myEdgeLength[i];
122     myEdge  [i] = *edge;
123     myEdgeID[i] = meshDS->ShapeToIndex( *edge );
124     if ( !theIsForward ) myEdge[i].Reverse();
125
126     if ( theFace.IsNull() )
127       BRep_Tool::Range( *edge, myFirst[i], myLast[i] );
128     else
129       myC2d[i] = BRep_Tool::CurveOnSurface( *edge, theFace, myFirst[i], myLast[i] );
130     if ( myEdge[i].Orientation() == TopAbs_REVERSED )
131       std::swap( myFirst[i], myLast[i] );
132
133     // check if the edge has a non-uniform parametrization (issue 0020705)
134     if ( !myC2d[i].IsNull() )
135     {
136       if ( myEdgeLength[i] > DBL_MIN )
137       {
138         Geom2dAdaptor_Curve A2dC( myC2d[i],
139                                   std::min( myFirst[i], myLast[i] ),
140                                   std::max( myFirst[i], myLast[i] ));
141         double p2 = myFirst[i]+(myLast[i]-myFirst[i])/2.;
142         double p4 = myFirst[i]+(myLast[i]-myFirst[i])/4.;
143         double d2 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p2 );
144         double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
145         //cout<<"len = "<<len<<"  d2 = "<<d2<<"  fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
146         myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
147         Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],d2,d4);
148         myC3dAdaptor[i].Load( C3d, d2,d4 );
149       }
150       else
151       {
152         const TopoDS_Vertex& V = SMESH_MesherHelper::IthVertex( 0, *edge );
153         Handle(Geom_Curve) C3d = new Geom_Line( BRep_Tool::Pnt( V ), gp::DX() );
154         myC3dAdaptor[i].Load( C3d, 0, 0.5 * BRep_Tool::Tolerance( V ));
155       }
156     }
157     // reverse a proxy sub-mesh
158     if ( !theIsForward )
159       reverseProxySubmesh( myEdge[i] );
160
161   } // loop on edges
162
163   // count nodes and segments
164   NbPoints( /*update=*/true );
165
166   if ( nbEdges > 1 && myLength > DBL_MIN ) {
167     const double degenNormLen = 1.e-5;
168     double totLength = myLength;
169     if ( nbDegen )
170       totLength += myLength * degenNormLen * nbDegen;
171     double prevNormPar = 0;
172     for ( int i = 0; i < nbEdges; ++i ) {
173       if ( myEdgeLength[ i ] < DBL_MIN )
174         myEdgeLength[ i ] = myLength * degenNormLen;
175       myNormPar[ i ] = prevNormPar + myEdgeLength[i]/totLength;
176       prevNormPar = myNormPar[ i ];
177     }
178   }
179   myNormPar[nbEdges-1] = 1.;
180   //dump();
181 }
182
183 //================================================================================
184 /*!
185  * \brief Constructor of a side for vertex using data from other FaceSide
186  */
187 //================================================================================
188
189 StdMeshers_FaceSide::StdMeshers_FaceSide(const StdMeshers_FaceSide*  theSide,
190                                          const SMDS_MeshNode*        theNode,
191                                          const gp_Pnt2d*             thePnt2d1,
192                                          const gp_Pnt2d*             thePnt2d2,
193                                          const Handle(Geom2d_Curve)& theC2d,
194                                          const double                theUFirst,
195                                          const double                theULast)
196 {
197   myC2d.push_back      ( theC2d );
198   myFirst.push_back    ( theUFirst );
199   myLast.push_back     ( theULast );
200   myNormPar.push_back  ( 1. );
201   myIsUniform.push_back( true );
202   myEdgeID.push_back   ( 0 );
203   myLength       = 0;
204   myProxyMesh    = theSide->myProxyMesh;
205   myDefaultPnt2d = *thePnt2d1;
206   myPoints       = theSide->GetUVPtStruct();
207   myNbPonits     = myPoints.size();
208   myNbSegments   = theSide->myNbSegments;
209   if ( thePnt2d2 )
210     for ( size_t i = 0; i < myPoints.size(); ++i )
211     {
212       double r = i / ( myPoints.size() - 1. );
213       myPoints[i].u = (1-r) * thePnt2d1->X() + r * thePnt2d2->X();
214       myPoints[i].v = (1-r) * thePnt2d1->Y() + r * thePnt2d2->Y();
215       myPoints[i].node = theNode;
216     }
217   else
218     for ( size_t i = 0; i < myPoints.size(); ++i )
219     {
220       myPoints[i].u = thePnt2d1->X();
221       myPoints[i].v = thePnt2d1->Y();
222       myPoints[i].node = theNode;
223     }
224 }
225
226 //================================================================================
227 /*
228  * Create a side from an UVPtStructVec
229  */
230 //================================================================================
231
232 StdMeshers_FaceSide::StdMeshers_FaceSide(UVPtStructVec&     theSideNodes,
233                                          const TopoDS_Face& theFace)
234 {
235   myEdge.resize( 1 );
236   myEdgeID.resize( 1, -1 );
237   myC2d.resize( 1 );
238   myC3dAdaptor.resize( 1 );
239   myFirst.resize( 1, 0. );
240   myLast.resize( 1, 1. );
241   myNormPar.resize( 1, 1. );
242   myIsUniform.resize( 1, 1 );
243   myMissingVertexNodes = myIgnoreMediumNodes = false;
244   myDefaultPnt2d.SetCoord( 1e100, 1e100 );
245
246   myFace       = theFace;
247   myPoints     = theSideNodes;
248   myNbPonits   = myPoints.size();
249   myNbSegments = myNbPonits + 1;
250
251   myLength = 0;
252   if ( !myPoints.empty() )
253   {
254     myPoints[0].normParam = 0;
255     if ( myPoints[0].node &&
256          myPoints.back().node &&
257          myPoints[ myNbPonits/2 ].node )
258     {
259       gp_Pnt pPrev = SMESH_TNodeXYZ( myPoints[0].node );
260       for ( size_t i = 1; i < myPoints.size(); ++i )
261       {
262         gp_Pnt p = SMESH_TNodeXYZ( myPoints[i].node );
263         myLength += p.Distance( pPrev );
264         myPoints[i].normParam = myLength;
265         pPrev = p;
266       }
267     }
268     else if ( !theFace.IsNull() )
269     {
270       TopLoc_Location loc;
271       Handle(Geom_Surface) surf = BRep_Tool::Surface( theFace, loc );
272       gp_Pnt pPrev = surf->Value( myPoints[0].u, myPoints[0].v );
273       for ( size_t i = 1; i < myPoints.size(); ++i )
274       {
275         gp_Pnt p = surf->Value( myPoints[i].u, myPoints[i].v );
276         myLength += p.Distance( pPrev );
277         myPoints[i].normParam = myLength;
278         pPrev = p;
279       }
280     }
281     else
282     {
283       gp_Pnt2d pPrev = myPoints[0].UV();
284       for ( size_t i = 1; i < myPoints.size(); ++i )
285       {
286         gp_Pnt2d p = myPoints[i].UV();
287         myLength += p.Distance( pPrev );
288         myPoints[i].normParam = myLength;
289         pPrev = p;
290       }
291     }
292     if ( myLength > std::numeric_limits<double>::min() )
293       for ( size_t i = 1; i < myPoints.size(); ++i )
294         myPoints[i].normParam /= myLength;
295   }
296   myEdgeLength.resize( 1, myLength );
297 }
298
299 //================================================================================
300 /*
301  * Return info on nodes on the side
302  */
303 //================================================================================
304
305 const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
306                                                              double constValue) const
307 {
308   if ( myPoints.empty() )
309   {
310     if ( NbEdges() == 0 ) return myPoints;
311
312     StdMeshers_FaceSide* me = const_cast< StdMeshers_FaceSide* >( this );
313     SMESHDS_Mesh*    meshDS = myProxyMesh->GetMeshDS();
314     SMESH_MesherHelper helper( *myProxyMesh->GetMesh() );
315     bool paramOK;
316     double eps = 1e-100;
317
318     // sort nodes of all edges putting them into a map
319
320     map< double, const SMDS_MeshNode*>            u2node;
321     vector< pair< double, const SMDS_MeshNode*> > u2nodeVec;
322     vector<const SMDS_MeshNode*>                  nodes;
323     set<const SMDS_MeshNode*>                     vertexNodes;
324     vector< const SMESH_ProxyMesh::SubMesh* >     proxySubMesh( myEdge.size() );
325     int nbProxyNodes = 0;
326     size_t iE;
327
328     for ( iE = 0; iE < myEdge.size(); ++iE )
329     {
330       proxySubMesh[iE] = myProxyMesh->GetProxySubMesh( myEdge[iE] );
331       if ( proxySubMesh[iE] )
332       {
333         if ( proxySubMesh[iE]->GetUVPtStructVec().empty() ) {
334           proxySubMesh[iE] = 0;
335         }
336         else {
337           nbProxyNodes += proxySubMesh[iE]->GetUVPtStructVec().size() - 1;
338           if ( iE+1 == myEdge.size() )
339             ++nbProxyNodes;
340           continue;
341         }
342       }
343
344       // Add 1st vertex node of a current edge
345       const SMDS_MeshNode* node = VertexNode( iE );
346       const double  prevNormPar = ( iE == 0 ? 0 : myNormPar[ iE-1 ]); // normalized param
347       if ( node ) // nodes on internal vertices may be missing
348       {
349         if ( vertexNodes.insert( node ).second )
350           u2node.insert( u2node.end(), make_pair( prevNormPar, node ));
351       }
352       else if ( iE == 0 )
353       {
354         if ( nodes.empty() ) {
355           for ( ++iE; iE < myEdge.size(); ++iE )
356             if (( node = VertexNode( iE ))) {
357               u2node.insert( make_pair( prevNormPar, node ));
358               break;
359             }
360           --iE;
361         }
362         if ( !node )
363           return myPoints;
364         vertexNodes.insert( node );
365       }
366
367       // Add internal nodes
368       nodes.clear();
369       if ( !GetEdgeNodes( iE, nodes, /*v0=*/false, /*v1=*/false ))
370         return myPoints;
371       if ( !nodes.empty() )
372       {
373         u2nodeVec.clear();
374         double paramSize = myLast[iE] - myFirst[iE];
375         double r         = myNormPar[iE] - prevNormPar;
376         helper.SetSubShape( myEdge[iE] );
377         helper.ToFixNodeParameters( true );
378         if ( !myIsUniform[iE] )
379           for ( size_t i = 0; i < nodes.size(); ++i )
380           {
381             double     u = helper.GetNodeU( myEdge[iE], nodes[i], 0, &paramOK );
382             double aLenU = GCPnts_AbscissaPoint::Length( me->myC3dAdaptor[iE], myFirst[iE], u );
383             if ( myEdgeLength[iE] < aLenU ) // nonregression test "3D_mesh_NETGEN/G6"
384             {
385               u2nodeVec.clear();
386               break;
387             }
388             double normPar = prevNormPar + r * aLenU / myEdgeLength[iE];
389             u2nodeVec.push_back( make_pair( normPar, nodes[i] ));
390           }
391         if ( u2nodeVec.empty() )
392           for ( size_t i = 0; i < nodes.size(); ++i )
393           {
394             double u = helper.GetNodeU( myEdge[iE], nodes[i], 0, &paramOK );
395             // paramSize is signed so orientation is taken into account
396             double normPar = prevNormPar + r * ( u - myFirst[iE] ) / paramSize;
397             u2nodeVec.push_back( make_pair( normPar, nodes[i] ));
398           }
399         for ( size_t j = 0; j < u2nodeVec.size(); ++j )
400           u2node.insert( u2node.end(), u2nodeVec[j] );
401       }
402     } // loop on myEdge's
403
404     if ( u2node.empty() ) return myPoints;
405
406     // Add 2nd VERTEX node for a last EDGE
407     {
408       const SMDS_MeshNode* node;
409       if ( IsClosed() )
410         node = u2node.begin()->second;
411       else
412       {
413         node = VertexNode( iE );
414         while ( !node && iE > 0 )
415           node = VertexNode( --iE );
416         if ( !node )
417           return myPoints;
418       }
419       if ( u2node.rbegin()->second == node )
420         u2node.erase( --u2node.end() );
421       u2node.insert( u2node.end(), make_pair( 1., node ));
422     }
423
424     if ( u2node.size() + nbProxyNodes != myNbPonits &&
425          u2node.size() + nbProxyNodes != NbPoints( /*update=*/true ))
426     {
427       MESSAGE("Wrong node parameters on edges, u2node.size():"
428               <<u2node.size()<<" !=  myNbPonits:"<<myNbPonits);
429       return myPoints;
430     }
431
432     // fill array of UVPtStruct
433
434     UVPtStructVec& points = me->myPoints;
435     points.resize( myNbPonits );
436
437     int iPt = 0;
438     double prevNormPar = 0, paramSize = myNormPar[ 0 ];
439     map< double, const SMDS_MeshNode*>::iterator u_node = u2node.begin();
440     for ( size_t iE = 0; iE < myEdge.size(); ++iE )
441     {
442       if ( proxySubMesh[ iE ] ) // copy data from a proxy sub-mesh
443       {
444         const UVPtStructVec& edgeUVPtStruct = proxySubMesh[iE]->GetUVPtStructVec();
445         std::copy( edgeUVPtStruct.begin(), edgeUVPtStruct.end(), & points[iPt] );
446         // check orientation
447         double du1 = edgeUVPtStruct.back().param - edgeUVPtStruct[0].param;
448         double du2 = myLast[iE] - myFirst[iE];
449         if ( du1 * du2 < 0 )
450         {
451           std::reverse( & points[iPt], & points[iPt + edgeUVPtStruct.size()]);
452           for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
453             points[iPt+i].normParam = 1. - points[iPt+i].normParam;
454         }
455         // update normalized params
456         if ( myEdge.size() > 1 ) {
457           for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i, ++iPt )
458           {
459             UVPtStruct & uvPt = points[iPt];
460             uvPt.normParam    = prevNormPar + uvPt.normParam * paramSize;
461             uvPt.x = uvPt.y   = uvPt.normParam;
462           }
463           --iPt; // to point to the 1st VERTEX of the next EDGE
464         }
465       }
466       else
467       {
468         for ( ; u_node != u2node.end(); ++u_node, ++iPt )
469         {
470           if ( myNormPar[ iE ]-eps < u_node->first )
471             break; // u_node is at VERTEX of the next EDGE 
472
473           UVPtStruct & uvPt = points[iPt];
474           uvPt.node       = u_node->second;
475           // -- normParam, x, y --------------------------------
476           uvPt.normParam  = u_node->first;
477           uvPt.x = uvPt.y = uvPt.normParam;
478           // -- U ----------------------------------------------
479           const SMDS_EdgePosition* epos =
480             dynamic_cast<const SMDS_EdgePosition*>(uvPt.node->GetPosition());
481           if ( epos && uvPt.node->getshapeId() == myEdgeID[iE] ) {
482             uvPt.param = epos->GetUParameter();
483           }
484           else {
485             double r = ( uvPt.normParam - prevNormPar )/ paramSize;
486             uvPt.param = ( r > 0.5 ? myLast[iE] : myFirst[iE] );
487           }
488           // -- UV ---------------------------------------------
489           if ( !myC2d[ iE ].IsNull() ) {
490             gp_Pnt2d p = myC2d[ iE ]->Value( uvPt.param );
491             uvPt.u = p.X();
492             uvPt.v = p.Y();
493           }
494           else {
495             uvPt.u = uvPt.v = 1e+100;
496           }
497         }
498       }
499       // prepare for the next EDGE
500       if ( iE+1 < myEdge.size() )
501       {
502         prevNormPar = myNormPar[ iE ];
503         paramSize   = myNormPar[ iE+1 ] - prevNormPar;
504       }
505     } // loop on myEdge's
506
507     // set <constValue>
508     if ( isXConst )
509       for ( iPt = 0; iPt < points.size(); ++iPt ) points[ iPt ].x = constValue;
510     else
511       for ( iPt = 0; iPt < points.size(); ++iPt ) points[ iPt ].y = constValue;
512
513   } // if ( myPoints.empty())
514
515   return myPoints;
516 }
517
518 //================================================================================
519 /*!
520  * \brief Falsificate info on nodes
521  * \param nbSeg - nb of segments on the side
522  * \retval UVPtStruct* - array of data structures
523  */
524 //================================================================================
525
526 const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int    nbSeg,
527                                                                   bool   isXConst,
528                                                                   double constValue) const
529 {
530   if ( myFalsePoints.empty() ) {
531
532     if ( NbEdges() == 0 ) return myFalsePoints;
533
534     vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myFalsePoints );
535     points->resize( nbSeg+1 );
536
537     int EdgeIndex = 0;
538     double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
539     for (int i = 0 ; i < myFalsePoints.size(); ++i ) {
540       double normPar = double(i) / double(nbSeg);
541       UVPtStruct & uvPt = (*points)[i];
542       uvPt.node = 0;
543       uvPt.x = uvPt.y = uvPt.param = uvPt.normParam = normPar;
544       if ( isXConst ) uvPt.x = constValue;
545       else            uvPt.y = constValue;
546       if ( myNormPar[ EdgeIndex ] < normPar ) {
547         prevNormPar = myNormPar[ EdgeIndex ];
548         ++EdgeIndex;
549         paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
550       }
551       double r = ( normPar - prevNormPar )/ paramSize;
552       uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
553       if ( !myC2d[ EdgeIndex ].IsNull() ) {
554         gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
555         uvPt.u = p.X();
556         uvPt.v = p.Y();
557       }
558       else {
559         uvPt.u = uvPt.v = 1e+100;
560       }
561     }
562   }
563   return myFalsePoints;
564 }
565
566 //=======================================================================
567 //function : GetOrderedNodes
568 //purpose  : Return nodes in the order they encounter while walking along the side
569 //=======================================================================
570
571 std::vector<const SMDS_MeshNode*> StdMeshers_FaceSide::GetOrderedNodes(int theEdgeInd) const
572 {
573   vector<const SMDS_MeshNode*> resultNodes;
574   if ( myPoints.empty() || ( theEdgeInd >= 0 && NbEdges() > 0 ))
575   {
576     if ( NbEdges() == 0 ) return resultNodes;
577
578     SMESHDS_Mesh* meshDS = myProxyMesh->GetMeshDS();
579     SMESH_MesherHelper helper(*myProxyMesh->GetMesh());
580     bool paramOK = true;
581
582     // Sort nodes of all edges putting them into a map
583
584     map< double, const SMDS_MeshNode*>        u2node;
585     vector<const SMDS_MeshNode*>              nodes;
586     set<const SMDS_MeshNode*>                 vertexNodes;
587     int iE = 0, iEnd = myEdge.size();
588     if ( theEdgeInd >= 0 )
589     {
590       iE   = theEdgeInd % NbEdges();
591       iEnd = iE + 1;
592     }
593     for ( iE = 0; iE < iEnd; ++iE )
594     {
595       double prevNormPar = ( iE == 0 ? 0 : myNormPar[ iE-1 ]); // normalized param
596
597       const SMESH_ProxyMesh::SubMesh* proxySM = myProxyMesh->GetProxySubMesh( myEdge[iE] );
598       if ( proxySM )
599       {
600         const UVPtStructVec& points = proxySM->GetUVPtStructVec();
601         for ( size_t i = 0; i < points.size(); ++i )
602           u2node.insert( make_pair( prevNormPar + points[i].normParam, points[i].node ));
603         continue;
604       }
605
606       // Add 1st vertex node of a current EDGE
607       const SMDS_MeshNode* node = VertexNode( iE );
608       if ( node ) { // nodes on internal vertices may be missing
609         if ( vertexNodes.insert( node ).second )
610           u2node.insert( make_pair( prevNormPar, node ));
611       }
612       else if ( iE == 0 )
613       {
614         if ( nodes.empty() ) {
615           for ( ++iE; iE < iEnd; ++iE )
616             if (( node = VertexNode( iE ))) {
617               u2node.insert( make_pair( prevNormPar, node ));
618               break;
619             }
620           --iE;
621         }
622         if ( !node )
623           return resultNodes;
624         vertexNodes.insert( node );
625       }
626
627       // Add internal nodes
628       nodes.clear();
629       if ( !GetEdgeNodes( iE, nodes, /*v0=*/false, /*v1=*/false ))
630         return resultNodes;
631       if ( !nodes.empty() )
632       {
633         double paramSize = myLast[iE] - myFirst[iE];
634         double r         = myNormPar[iE] - prevNormPar;
635         helper.SetSubShape( myEdge[iE] );
636         helper.ToFixNodeParameters( true );
637         for ( size_t i = 0; i < nodes.size(); ++i )
638         {
639           double u = helper.GetNodeU( myEdge[iE], nodes[i], 0, &paramOK );
640           // paramSize is signed so orientation is taken into account
641           double normPar = prevNormPar + r * ( u - myFirst[iE] ) / paramSize;
642           u2node.insert( u2node.end(), make_pair( normPar, nodes[i] ));
643         }
644       }
645
646     } // loop on myEdges
647
648     if ( u2node.empty() ) return resultNodes;
649
650     // Add 2nd vertex node for a last EDGE
651     {
652       const SMDS_MeshNode* node;
653       if ( IsClosed() && theEdgeInd < 0 )
654         node = u2node.begin()->second;
655       else
656       {
657         node = VertexNode( iE );
658         while ( !node && iE > 0 )
659           node = VertexNode( --iE );
660         if ( !node )
661           return resultNodes;
662       }
663       if ( u2node.rbegin()->second == node )
664         u2node.erase( --u2node.end() );
665       u2node.insert( u2node.end(), make_pair( 1., node ));
666     }
667
668     // Fill the result vector
669
670     if ( u2node.size() == myNbPonits )
671     {
672       resultNodes.reserve( u2node.size() );
673       map< double, const SMDS_MeshNode*>::iterator u2n = u2node.begin();
674       for ( ; u2n != u2node.end(); ++u2n )
675         resultNodes.push_back( u2n->second );
676     }
677   }
678   else
679   {
680     resultNodes.resize( myPoints.size() );
681     for ( size_t i = 0; i < myPoints.size(); ++i )
682       resultNodes[i] = myPoints[i].node;
683   }
684
685   return resultNodes;
686 }
687
688 //================================================================================
689 /*!
690  * \brief Return (unsorted) nodes of the i-th EDGE.
691  *        Nodes moved to other geometry by MergeNodes() are also returned.
692  *  \retval bool - is OK
693  */
694 //================================================================================
695
696 bool StdMeshers_FaceSide::GetEdgeNodes(size_t                        i,
697                                        vector<const SMDS_MeshNode*>& nodes,
698                                        bool                          inlude1stVertex,
699                                        bool                          inludeLastVertex) const
700 {
701   if ( i >= myEdge.size() )
702     return false;
703
704   SMESH_Mesh*     mesh = myProxyMesh->GetMesh();
705   SMESHDS_Mesh* meshDS = mesh->GetMeshDS();
706   SMESHDS_SubMesh*  sm = meshDS->MeshElements( myEdge[i] );
707
708   if ( inlude1stVertex )
709   {
710     const SMDS_MeshNode* n0 = VertexNode( i );
711     if ( !n0 ) return false;
712     nodes.push_back( n0 );
713   }
714
715   if ( sm && ( sm->NbElements() > 0 || sm->NbNodes() > 0 ))
716   {
717     if ( mesh->HasModificationsToDiscard() ) // check nb of nodes on the EDGE sub-mesh
718     {
719       int iQuad    = sm->NbElements() ? sm->GetElements()->next()->IsQuadratic() : 0;
720       int nbExpect = sm->NbElements() - 1 + iQuad * sm->NbElements();
721       if ( nbExpect != sm->NbNodes() ) // some nodes are moved from the EDGE by MergeNodes()
722       {
723         // add nodes of all segments
724         typedef set< const SMDS_MeshNode* > TNodeSet;
725         TNodeSet sharedNodes;
726         SMDS_ElemIteratorPtr segIt = sm->GetElements();
727         while ( segIt->more() )
728         {
729           const SMDS_MeshElement* seg = segIt->next();
730           if ( seg->GetType() != SMDSAbs_Edge )
731             continue;
732           for ( int i = 0; i < 3-myIgnoreMediumNodes; ++i )
733           {
734             const SMDS_MeshNode* n = seg->GetNode( i );
735             if ( i == 2 ) // medium node
736             {
737               nodes.push_back( n );
738             }
739             else
740             {
741               pair<TNodeSet::iterator, bool> it2new = sharedNodes.insert( n );
742               if ( !it2new.second ) // n encounters twice == it's on EDGE, not on VERTEX
743               {
744                 nodes.push_back( n );
745                 sharedNodes.erase( it2new.first );
746               }
747             }
748           }
749         }
750       }
751     }
752     if ( nodes.size() < 2 ) // add nodes assigned to the EDGE
753     {
754       SMDS_NodeIteratorPtr nItr = sm->GetNodes();
755       while ( nItr->more() )
756       {
757         const SMDS_MeshNode* n = nItr->next();
758         if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( n, SMDSAbs_Edge ))
759           continue;
760         nodes.push_back( n );
761       }
762     }
763   } // if ( sm && sm->NbElements() > 0 )
764
765   if ( inludeLastVertex )
766   {
767     const SMDS_MeshNode* n1 = VertexNode( i+1 );
768     if ( !n1 ) return false;
769     nodes.push_back( n1 );
770   }
771   return true;
772 }
773
774 //================================================================================
775 /*!
776  * \brief Return a node from the i-th VERTEX (count starts from zero)
777  *        Nodes moved to other geometry by MergeNodes() are also returned.
778  *  \param [in] i - the VERTEX index
779  *  \param [out] isMoved - returns \c true if the found node is moved by MergeNodes()
780  *  \return const SMDS_MeshNode* - the found node
781  */
782 //================================================================================
783
784 const SMDS_MeshNode* StdMeshers_FaceSide::VertexNode(std::size_t i, bool* isMoved) const
785 {
786   TopoDS_Vertex V = ( i >= myEdge.size() ) ? LastVertex() : FirstVertex(i);
787
788   const SMDS_MeshNode* n = SMESH_Algo::VertexNode( V, myProxyMesh->GetMeshDS() );
789
790   if ( !n && !myEdge.empty() && myProxyMesh->GetMesh()->HasModificationsToDiscard() )
791   {
792     size_t iE = ( i < myEdge.size() ) ? i : myEdge.size()-1;
793     SMESHDS_SubMesh* sm = myProxyMesh->GetMeshDS()->MeshElements( myEdgeID[ iE ]);
794
795     n = SMESH_Algo::VertexNode( V, sm, myProxyMesh->GetMesh(), /*checkV=*/false );
796
797     if (( !n ) &&
798         (( i > 0 && i < NbEdges() ) || IsClosed() ))
799     {
800       iE = SMESH_MesherHelper::WrapIndex( int(i)-1, NbEdges() );
801       sm = myProxyMesh->GetMeshDS()->MeshElements( myEdgeID[ iE ]);
802       n  = SMESH_Algo::VertexNode( V, sm, myProxyMesh->GetMesh(), /*checkV=*/false );
803     }
804
805     if ( n && n->GetPosition()->GetDim() == 1 ) // check that n does not lie on an EDGE of myFace
806     {
807       TopoDS_Shape S = SMESH_MesherHelper::GetSubShapeByNode( n, myProxyMesh->GetMeshDS() );
808       if ( SMESH_MesherHelper::IsSubShape( S, myFace ))
809         n = 0; // VERTEX ignored by e.g. Composite Wire Discretization algo
810     }
811     if ( isMoved )
812       *isMoved = n;
813   }
814   return n;
815 }
816
817 //================================================================================
818 /*!
819  * \brief reverse order of vector elements
820   * \param vec - vector to reverse
821  */
822 //================================================================================
823
824 template <typename T > void reverse(vector<T> & vec)
825 {
826   std::reverse( vec.begin(), vec.end() );
827 }
828
829 //================================================================================
830 /*!
831  * \brief Change orientation of side geometry
832  */
833 //================================================================================
834
835 void StdMeshers_FaceSide::Reverse()
836 {
837   int nbEdges = myEdge.size();
838   for ( int i = nbEdges-1; i >= 0; --i ) {
839     std::swap( myFirst[i], myLast[i] );
840     if ( !myEdge[i].IsNull() )
841       myEdge[i].Reverse();
842     if ( i > 0 ) // at the first loop 1. is overwritten
843       myNormPar[i] = 1 - myNormPar[i-1];
844   }
845   if ( nbEdges > 1 ) {
846     reverse( myEdge );
847     reverse( myEdgeID );
848     reverse( myC2d );
849     //reverse( myC3dAdaptor );
850     reverse( myFirst );
851     reverse( myLast );
852     reverse( myNormPar );
853     reverse( myEdgeLength );
854     reverse( myIsUniform );
855   }
856   if ( nbEdges > 0 )
857   {
858     myNormPar[nbEdges-1]=1.;
859     if ( !myEdge[0].IsNull() )
860     {
861       for ( size_t i = 0; i < myEdge.size(); ++i )
862         reverseProxySubmesh( myEdge[i] );
863       myPoints.clear();
864       myFalsePoints.clear();
865     }
866     else
867     {
868       for ( size_t i = 0; i < myPoints.size(); ++i )
869       {
870         UVPtStruct & uvPt = myPoints[i];
871         uvPt.normParam = 1 - uvPt.normParam;
872         uvPt.x         = 1 - uvPt.x;
873         uvPt.y         = 1 - uvPt.y;
874       }
875       reverse( myPoints );
876
877       for ( size_t i = 0; i < myFalsePoints.size(); ++i )
878       {
879         UVPtStruct & uvPt = myFalsePoints[i];
880         uvPt.normParam = 1 - uvPt.normParam;
881         uvPt.x         = 1 - uvPt.x;
882         uvPt.y         = 1 - uvPt.y;
883       }
884       reverse( myFalsePoints );
885     }
886   }
887   for ( size_t i = 0; i < myEdge.size(); ++i )
888   {
889     if ( myEdge[i].IsNull() ) continue; // for a side on points only
890     double fp,lp;
891     Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],fp,lp);
892     if ( !C3d.IsNull() )
893       myC3dAdaptor[i].Load( C3d, fp,lp );
894   }
895 }
896
897 //=======================================================================
898 //function : SetIgnoreMediumNodes
899 //purpose  : Make ignore medium nodes
900 //=======================================================================
901
902 void StdMeshers_FaceSide::SetIgnoreMediumNodes(bool toIgnore)
903 {
904   if ( myIgnoreMediumNodes != toIgnore )
905   {
906     myIgnoreMediumNodes = toIgnore;
907
908     if ( !myPoints.empty() )
909     {
910       UVPtStructVec newPoints;
911       newPoints.reserve( myPoints.size()/2 + 1 );
912       for ( size_t i = 0; i < myPoints.size(); i += 2 )
913         newPoints.push_back( myPoints[i] );
914
915       myPoints.swap( newPoints );
916     }
917     else
918     {
919       NbPoints( /*update=*/true );
920     }
921   }
922 }
923
924 //=======================================================================
925 //function : NbPoints
926 //purpose  : Return nb nodes on edges and vertices (+1 to be == GetUVPtStruct().size() )
927 //           Call it with update == true if mesh of this side can be recomputed
928 //           since creation of this side
929 //=======================================================================
930
931 int StdMeshers_FaceSide::NbPoints(const bool update) const
932 {
933   if ( !myPoints.empty() )
934     return myPoints.size();
935
936   // if ( !myFalsePoints.empty() )
937   //   return myFalsePoints.size();
938
939   if ( update && myEdge.size() > 0 )
940   {
941     StdMeshers_FaceSide* me = (StdMeshers_FaceSide*) this;
942     me->myNbPonits = 0;
943     me->myNbSegments = 0;
944     me->myMissingVertexNodes = false;
945
946     for ( int i = 0; i < NbEdges(); ++i )
947     {
948       if ( const SMESHDS_SubMesh* sm = myProxyMesh->GetSubMesh( Edge(i) )) {
949         int nbN = sm->NbNodes();
950         if ( sm->NbElements() > 0 ) {
951           nbN = sm->NbElements() - 1; // nodes can be moved to other shapes by MergeNodes()
952           if ( !myIgnoreMediumNodes &&
953                sm->GetElements()->next()->IsQuadratic() )
954             nbN += sm->NbElements();
955         }
956         me->myNbPonits   += nbN;
957         me->myNbSegments += sm->NbElements();
958       }
959     }
960
961     std::set< const SMDS_MeshNode* > vNodes;
962     for ( int i = 0; i <= NbEdges(); ++i )
963       if ( const SMDS_MeshNode* n = VertexNode( i ))
964         vNodes.insert( n );
965       else
966         me->myMissingVertexNodes = true;
967     me->myNbPonits += vNodes.size();
968
969     if ( IsClosed() )
970       me->myNbPonits++; // closing node is repeated
971   }
972   return myNbPonits;
973 }
974
975 //=======================================================================
976 //function : NbSegments
977 //purpose  : Return nb edges
978 //           Call it with update == true if mesh of this side can be recomputed
979 //           since creation of this side
980 //=======================================================================
981
982 int StdMeshers_FaceSide::NbSegments(const bool update) const
983 {
984   return NbPoints( update ), myNbSegments;
985 }
986
987 //================================================================================
988 /*!
989  * \brief Reverse UVPtStructVec if a proxy sub-mesh of E
990  */
991 //================================================================================
992
993 void StdMeshers_FaceSide::reverseProxySubmesh( const TopoDS_Edge& E )
994 {
995   if ( !myProxyMesh ) return;
996   if ( const SMESH_ProxyMesh::SubMesh* sm = myProxyMesh->GetProxySubMesh( E ))
997   {
998     UVPtStructVec& edgeUVPtStruct = (UVPtStructVec& ) sm->GetUVPtStructVec();
999     for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
1000     {
1001       UVPtStruct & uvPt = edgeUVPtStruct[i];
1002       uvPt.normParam = 1 - uvPt.normParam;
1003       uvPt.x         = 1 - uvPt.x;
1004       uvPt.y         = 1 - uvPt.y;
1005     }
1006     reverse( edgeUVPtStruct );
1007   }
1008 }
1009
1010 //================================================================================
1011 /*!
1012  * \brief Show side features
1013  */
1014 //================================================================================
1015
1016 void StdMeshers_FaceSide::dump(const char* msg) const
1017 {
1018 #ifdef _DEBUG_
1019   if (msg) MESSAGE ( std::endl << msg );
1020   MESSAGE_BEGIN ("NB EDGES: "<< myEdge.size() );
1021   MESSAGE_ADD ( "nbPoints: "<<myNbPonits<<" vecSize: " << myPoints.size()<<" "<<myFalsePoints.size() );
1022   for ( int i=0; i<myEdge.size(); ++i)
1023   {
1024     MESSAGE_ADD ( "\t"<<i+1 );
1025     MESSAGE_ADD ( "\tEDGE: " );
1026     if (myEdge[i].IsNull()) {
1027       MESSAGE_ADD ( "NULL" );
1028     }
1029     else {
1030       TopAbs::Print(myEdge[i].Orientation(),cout)<<" "<<myEdge[i].TShape().operator->()<<endl;
1031       MESSAGE_ADD ( "\tV1: " << TopExp::FirstVertex( myEdge[i], 1).TShape().operator->()
1032                  << "  V2: " << TopExp::LastVertex( myEdge[i], 1).TShape().operator->() );
1033     }
1034     MESSAGE_ADD ( "\tC2d: ");
1035     
1036     if (myC2d[i].IsNull()) {
1037       MESSAGE_ADD ( "NULL" );
1038     }
1039     else {
1040       MESSAGE_ADD ( myC2d[i].operator->() );
1041     }
1042       
1043     MESSAGE_ADD ( "\tF: "<<myFirst[i]<< " L: "<< myLast[i] );
1044     MESSAGE_END ( "\tnormPar: "<<myNormPar[i]<<endl );
1045   }
1046 #endif
1047 }
1048
1049 //================================================================================
1050 /*!
1051  * \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
1052   * \retval Adaptor2d_Curve2d* - 
1053  */
1054 //================================================================================
1055
1056 struct Adaptor2dCurve2d : public Adaptor2d_Curve2d
1057 {
1058   const StdMeshers_FaceSide* mySide;
1059   Adaptor2dCurve2d(const StdMeshers_FaceSide* faceSide):mySide(faceSide) {}
1060   gp_Pnt2d Value(const Standard_Real U) const { return mySide->Value2d( U ); }
1061   Standard_Real FirstParameter() const { return 0; }
1062   Standard_Real LastParameter() const { return 1; }
1063 };
1064
1065 Adaptor2d_Curve2d* StdMeshers_FaceSide::GetCurve2d() const
1066 {
1067   return new Adaptor2dCurve2d( this );
1068 }
1069
1070 //================================================================================
1071 /*!
1072  * \brief Creates a fully functional Adaptor_Curve
1073  */
1074 //================================================================================
1075
1076 BRepAdaptor_CompCurve* StdMeshers_FaceSide::GetCurve3d() const
1077 {
1078   if ( myEdge.empty() )
1079     return 0;
1080
1081   TopoDS_Wire aWire;
1082   BRep_Builder aBuilder;
1083   aBuilder.MakeWire(aWire);
1084   for ( int i=0; i<myEdge.size(); ++i )
1085     aBuilder.Add( aWire, myEdge[i] );
1086
1087   if ( myEdge.size() == 2 && IsClosed() )
1088     aWire.Closed(true); // issue 0021141
1089
1090   return new BRepAdaptor_CompCurve( aWire );
1091 }
1092
1093 //================================================================================
1094 /*!
1095  * \brief Return 2D point by normalized parameter
1096   * \param U - normalized parameter value
1097   * \retval gp_Pnt2d - point
1098  */
1099 //================================================================================
1100
1101 gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
1102 {
1103   if ( !myC2d[0].IsNull() ) {
1104     int i = EdgeIndex( U );
1105     double prevU = i ? myNormPar[ i-1 ] : 0;
1106     double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
1107
1108     double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
1109     
1110     // check parametrization of curve
1111     if( !myIsUniform[i] )
1112     {
1113       double aLen3dU = r * myEdgeLength[i] * ( myFirst[i]>myLast[i] ? -1. : 1.);
1114       GCPnts_AbscissaPoint AbPnt
1115         ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
1116       if( AbPnt.IsDone() ) {
1117         par = AbPnt.Parameter();
1118       }
1119     }
1120     return myC2d[ i ]->Value(par);
1121
1122   }
1123   else if ( !myPoints.empty() )
1124   {
1125     int i = U * double( myPoints.size()-1 );
1126     while ( i > 0 && myPoints[ i ].normParam > U )
1127       --i;
1128     while ( i+1 < myPoints.size() && myPoints[ i+1 ].normParam < U )
1129       ++i;
1130     double r = (( U - myPoints[ i ].normParam ) /
1131                 ( myPoints[ i+1 ].normParam - myPoints[ i ].normParam ));
1132     return ( myPoints[ i   ].UV() * ( 1 - r ) +
1133              myPoints[ i+1 ].UV() * r );
1134   }
1135   return myDefaultPnt2d;
1136 }
1137
1138 //================================================================================
1139 /*!
1140  * \brief Return XYZ by normalized parameter
1141   * \param U - normalized parameter value
1142   * \retval gp_Pnt - point
1143  */
1144 //================================================================================
1145
1146 gp_Pnt StdMeshers_FaceSide::Value3d(double U) const
1147 {
1148   int        i = EdgeIndex( U );
1149   double prevU = i ? myNormPar[ i-1 ] : 0;
1150   double     r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
1151
1152   double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
1153
1154   // check parametrization of curve
1155   if( !myIsUniform[i] )
1156   {
1157     double aLen3dU = r * myEdgeLength[i] * ( myFirst[i]>myLast[i] ? -1. : 1.);
1158     GCPnts_AbscissaPoint AbPnt
1159       ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
1160     if( AbPnt.IsDone() ) {
1161       par = AbPnt.Parameter();
1162     }
1163   }
1164   return myC3dAdaptor[ i ].Value(par);
1165 }
1166
1167 //================================================================================
1168 /*!
1169  * \brief Return wires of a face as StdMeshers_FaceSide's
1170  */
1171 //================================================================================
1172
1173 TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face&   theFace,
1174                                               SMESH_Mesh &         theMesh,
1175                                               const bool           theIgnoreMediumNodes,
1176                                               TError &             theError,
1177                                               SMESH_ProxyMesh::Ptr theProxyMesh,
1178                                               const bool           theCheckVertexNodes)
1179 {
1180   list< TopoDS_Edge > edges, internalEdges;
1181   list< int > nbEdgesInWires;
1182   int nbWires = SMESH_Block::GetOrderedEdges (theFace, edges, nbEdgesInWires);
1183
1184   // split list of all edges into separate wires
1185   TSideVector wires( nbWires );
1186   list< int >::iterator nbE = nbEdgesInWires.begin();
1187   list< TopoDS_Edge >::iterator from = edges.begin(), to = from;
1188   for ( int iW = 0; iW < nbWires; ++iW, ++nbE )
1189   {
1190     std::advance( to, *nbE );
1191     if ( *nbE == 0 ) // Issue 0020676
1192     {
1193       --nbWires;
1194       --iW;
1195       wires.resize( nbWires );
1196       continue;
1197     }
1198     list< TopoDS_Edge > wireEdges( from, to );
1199     // assure that there is a node on the first vertex
1200     // as StdMeshers_FaceSide::GetUVPtStruct() requires
1201     if ( wireEdges.front().Orientation() != TopAbs_INTERNAL ) // Issue 0020676
1202     {
1203       if ( theCheckVertexNodes )
1204         while ( !SMESH_Algo::VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
1205                                          theMesh.GetMeshDS()))
1206         {
1207           wireEdges.splice(wireEdges.end(), wireEdges,
1208                            wireEdges.begin(), ++wireEdges.begin());
1209           if ( from->IsSame( wireEdges.front() )) {
1210             theError = TError
1211               ( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"No nodes on vertices"));
1212             return TSideVector(0);
1213           }
1214         }
1215     }
1216     else if ( *nbE > 1 ) // Issue 0020676 (Face_pb_netgen.brep) - several internal edges in a wire
1217     {
1218       internalEdges.splice( internalEdges.end(), wireEdges, ++wireEdges.begin(), wireEdges.end());
1219     }
1220
1221     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, wireEdges, &theMesh,
1222                                                          /*isForward=*/true, theIgnoreMediumNodes,
1223                                                          theProxyMesh );
1224     wires[ iW ] = StdMeshers_FaceSidePtr( wire );
1225     from = to;
1226   }
1227   while ( !internalEdges.empty() )
1228   {
1229     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, internalEdges.back(), &theMesh,
1230                                                          /*isForward=*/true, theIgnoreMediumNodes,
1231                                                          theProxyMesh );
1232     wires.push_back( StdMeshers_FaceSidePtr( wire ));
1233     internalEdges.pop_back();
1234   }
1235   return wires;
1236 }
1237
1238 //================================================================================
1239 /*!
1240  * \brief Return 1st vertex of the i-the edge
1241  */
1242 //================================================================================
1243
1244 TopoDS_Vertex StdMeshers_FaceSide::FirstVertex(int i) const
1245 {
1246   TopoDS_Vertex v;
1247   if ( i < NbEdges() )
1248   {
1249     v = myEdge[i].Orientation() <= TopAbs_REVERSED ? // FORWARD || REVERSED
1250         TopExp::FirstVertex( myEdge[i], 1 )        :
1251         TopoDS::Vertex( TopoDS_Iterator( myEdge[i] ).Value() );
1252   }
1253   return v;
1254 }
1255
1256 //================================================================================
1257 /*!
1258  * \brief Return last vertex of the i-the edge
1259  */
1260 //================================================================================
1261
1262 TopoDS_Vertex StdMeshers_FaceSide::LastVertex(int i) const
1263 {
1264   TopoDS_Vertex v;
1265   if ( i < NbEdges() )
1266   {
1267     const TopoDS_Edge& edge = i<0 ? myEdge[ NbEdges() + i ] : myEdge[i];
1268     if ( edge.Orientation() <= TopAbs_REVERSED ) // FORWARD || REVERSED
1269       v = TopExp::LastVertex( edge, 1 );
1270     else
1271       for ( TopoDS_Iterator vIt( edge ); vIt.More(); vIt.Next() )
1272         v = TopoDS::Vertex( vIt.Value() );
1273   }
1274   return v;
1275 }
1276
1277 //================================================================================
1278 /*!
1279  * \brief Return \c true if the chain of EDGEs is closed
1280  */
1281 //================================================================================
1282
1283 bool StdMeshers_FaceSide::IsClosed() const
1284 {
1285   return myEdge.empty() ? false : FirstVertex().IsSame( LastVertex() );
1286 }