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