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