Salome HOME
Merge branch 'occ/shaper2smesh'
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.cxx
1 // Copyright (C) 2007-2019  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     } // loop on myEdge's
454
455     // Add 2nd VERTEX node for a last EDGE
456     if ( !proxySubMesh.back() )
457     {
458       if ( u2node.empty() ) return myPoints;
459
460       const SMDS_MeshNode* node;
461       if ( IsClosed() && !proxySubMesh[0] )
462         node = u2node.begin()->second;
463       else
464       {
465         node = VertexNode( iE );
466         if ( myProxyMesh->GetMesh()->HasModificationsToDiscard() )
467           while ( !node && iE > 1 ) // check intermediate VERTEXes
468             node = VertexNode( --iE );
469       }
470       if ( node )
471       {
472         if ( u2node.rbegin()->second == node &&
473              !fHelper.IsRealSeam  ( node->getshapeId() ) &&
474              !fHelper.IsDegenShape( node->getshapeId() ))
475           u2node.erase( --u2node.end() );
476
477         u2node.insert( u2node.end(), make_pair( 1., node ));
478       }
479     }
480
481     if ((int) u2node.size() + nbProxyNodes != myNbPonits &&
482         (int) u2node.size() + nbProxyNodes != NbPoints( /*update=*/true ))
483     {
484       return myPoints;
485     }
486     if (( u2node.size() > 0 ) &&
487         ( u2node.begin()->first < 0 || u2node.rbegin()->first > 1 ))
488     {
489       return myPoints;
490     }
491
492     // fill array of UVPtStruct
493
494     UVPtStructVec& points = me->myPoints;
495     points.resize( myNbPonits );
496
497     int iPt = 0;
498     double prevNormPar = 0, paramSize = myNormPar[ 0 ];
499     map< double, const SMDS_MeshNode*>::iterator u_node = u2node.begin();
500     for ( size_t iE = 0; iE < myEdge.size(); ++iE )
501     {
502       if ( proxySubMesh[ iE ] ) // copy data from a proxy sub-mesh
503       {
504         const UVPtStructVec& edgeUVPtStruct = proxySubMesh[iE]->GetUVPtStructVec();
505         UVPtStruct* pointsPtr = & points[iPt];
506         std::copy( edgeUVPtStruct.begin(), edgeUVPtStruct.end(), pointsPtr );
507         // check orientation
508         double du1 = edgeUVPtStruct.back().param - edgeUVPtStruct[0].param;
509         double du2 = myLast[iE] - myFirst[iE];
510         if ( du1 * du2 < 0 )
511         {
512           std::reverse( pointsPtr, pointsPtr + edgeUVPtStruct.size());
513           for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
514             pointsPtr[i].normParam = 1. - pointsPtr[i].normParam;
515         }
516         // update normalized params
517         if ( myEdge.size() > 1 ) {
518           for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
519           {
520             UVPtStruct & uvPt = pointsPtr[i];
521             uvPt.normParam    = prevNormPar + uvPt.normParam * paramSize;
522             uvPt.x = uvPt.y   = uvPt.normParam;
523           }
524           iPt += edgeUVPtStruct.size() - 1; // to point to the 1st VERTEX of the next EDGE
525         }
526         // update UV on a seam EDGE
527         if ( fHelper.IsRealSeam( myEdgeID[ iE ]))
528         {
529           // check if points lye on the EDGE
530           const UVPtStruct& pm = edgeUVPtStruct[ edgeUVPtStruct.size()/2 ];
531           gp_Pnt         pNode = SMESH_TNodeXYZ( pm.node );
532           gp_Pnt         pCurv = myC3dAdaptor[ iE ].Value( pm.param );
533           double           tol = BRep_Tool::Tolerance( myEdge[ iE ]) * 10;
534           bool   isPointOnEdge = ( pNode.SquareDistance( pCurv ) < tol * tol );
535           if ( isPointOnEdge )
536             for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
537               pointsPtr[i].SetUV( myC2d[ iE ]->Value( pointsPtr[i].param ).XY() );
538         }
539       }
540       else
541       {
542         for ( ; u_node != u2node.end(); ++u_node, ++iPt )
543         {
544           if ( myNormPar[ iE ]-eps < u_node->first )
545             break; // u_node is at VERTEX of the next EDGE
546
547           UVPtStruct & uvPt = points[iPt];
548           uvPt.node       = u_node->second;
549           // -- normParam, x, y --------------------------------
550           uvPt.normParam  = u_node->first;
551           uvPt.x = uvPt.y = uvPt.normParam;
552           // -- U ----------------------------------------------
553           SMDS_EdgePositionPtr epos = uvPt.node->GetPosition();
554           if ( epos && uvPt.node->getshapeId() == myEdgeID[iE] ) {
555             uvPt.param = epos->GetUParameter();
556           }
557           else {
558             double r = ( uvPt.normParam - prevNormPar )/ paramSize;
559             uvPt.param = ( r > 0.5 ? myLast[iE] : myFirst[iE] );
560           }
561           // -- UV ---------------------------------------------
562           if ( !myC2d[ iE ].IsNull() ) {
563             gp_Pnt2d p = myC2d[ iE ]->Value( uvPt.param );
564             uvPt.u = p.X();
565             uvPt.v = p.Y();
566           }
567           else {
568             uvPt.u = uvPt.v = 1e+100;
569           }
570         }
571       }
572       // prepare for the next EDGE
573       if ( iE+1 < myEdge.size() )
574       {
575         prevNormPar = myNormPar[ iE ];
576         paramSize   = myNormPar[ iE+1 ] - prevNormPar;
577       }
578     } // loop on myEdge's
579
580     // set <constValue>
581     if ( isXConst )
582       for ( iPt = 0; iPt < (int)points.size(); ++iPt ) points[ iPt ].x = constValue;
583     else
584       for ( iPt = 0; iPt < (int)points.size(); ++iPt ) points[ iPt ].y = constValue;
585
586   } // if ( myPoints.empty())
587
588   return myPoints;
589 }
590
591 //================================================================================
592 /*!
593  * \brief Falsificate info on nodes
594  * \param nbSeg - nb of segments on the side
595  * \retval UVPtStruct* - array of data structures
596  */
597 //================================================================================
598
599 const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int    nbSeg,
600                                                                   bool   isXConst,
601                                                                   double constValue) const
602 {
603   if ( myFalsePoints.empty() )
604   {
605     if ( NbEdges() == 0 ) return myFalsePoints;
606
607     vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myFalsePoints );
608     points->resize( nbSeg+1 );
609
610     int EdgeIndex = 0;
611     double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
612     gp_Pnt2d p;
613     for ( size_t i = 0 ; i < myFalsePoints.size(); ++i )
614     {
615       double normPar = double(i) / double(nbSeg);
616       UVPtStruct & uvPt = (*points)[i];
617       uvPt.node = 0;
618       uvPt.x = uvPt.y = uvPt.param = uvPt.normParam = normPar;
619       if ( isXConst ) uvPt.x = constValue;
620       else            uvPt.y = constValue;
621       if ( myNormPar[ EdgeIndex ] < normPar )
622       {
623         prevNormPar = myNormPar[ EdgeIndex ];
624         ++EdgeIndex;
625         paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
626       }
627       double r = ( normPar - prevNormPar )/ paramSize;
628       uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
629       if ( !myC2d[ EdgeIndex ].IsNull() )
630         p = myC2d[ EdgeIndex ]->Value( uvPt.param );
631       else
632         p = Value2d( normPar );
633       uvPt.u = p.X();
634       uvPt.v = p.Y();
635     }
636   }
637   return myFalsePoints;
638 }
639
640 //=======================================================================
641 //function : GetOrderedNodes
642 //purpose  : Return nodes in the order they encounter while walking along the side
643 //=======================================================================
644
645 std::vector<const SMDS_MeshNode*> StdMeshers_FaceSide::GetOrderedNodes(int theEdgeInd) const
646 {
647   vector<const SMDS_MeshNode*> resultNodes;
648   if ( myPoints.empty() || ( theEdgeInd >= 0 && NbEdges() > 0 ))
649   {
650     if ( NbEdges() == 0 ) return resultNodes;
651
652     //SMESHDS_Mesh* meshDS = myProxyMesh->GetMeshDS();
653     SMESH_MesherHelper  eHelper( *myProxyMesh->GetMesh() );
654     SMESH_MesherHelper& fHelper = * FaceHelper();
655     bool paramOK = true;
656
657     // Sort nodes of all edges putting them into a map
658
659     map< double, const SMDS_MeshNode*> u2node;
660     vector<const SMDS_MeshNode*>       nodes;
661     set<const SMDS_MeshNode*>          vertexNodes;
662     int iE = 0, iEnd = myEdge.size();
663     if ( theEdgeInd >= 0 )
664     {
665       iE   = theEdgeInd % NbEdges();
666       iEnd = iE + 1;
667     }
668     for ( ; iE < iEnd; ++iE )
669     {
670       double prevNormPar = ( iE == 0 ? 0 : myNormPar[ iE-1 ]); // normalized param
671
672       const SMESH_ProxyMesh::SubMesh* proxySM = myProxyMesh->GetProxySubMesh( myEdge[iE] );
673       if ( proxySM )
674       {
675         const UVPtStructVec& points = proxySM->GetUVPtStructVec();
676         for ( size_t i = 0; i < points.size(); ++i )
677           u2node.insert( make_pair( prevNormPar + points[i].normParam, points[i].node ));
678         continue;
679       }
680
681       // Add 1st vertex node of a current EDGE
682       const SMDS_MeshNode* node = VertexNode( iE );
683       if ( node ) { // nodes on internal vertices may be missing
684         if ( vertexNodes.insert( node ).second ||
685              fHelper.IsRealSeam  ( node->getshapeId() ) ||
686              fHelper.IsDegenShape( node->getshapeId() ))
687           u2node.insert( make_pair( prevNormPar, node ));
688       }
689       else if ( iE == 0 )
690       {
691         if ( nodes.empty() ) {
692           for ( ++iE; iE < iEnd; ++iE )
693             if (( node = VertexNode( iE ))) {
694               u2node.insert( make_pair( prevNormPar, node ));
695               break;
696             }
697           --iE;
698         }
699         if ( !node )
700           return resultNodes;
701         vertexNodes.insert( node );
702       }
703
704       // Add internal nodes
705       nodes.clear();
706       if ( !GetEdgeNodes( iE, nodes, /*v0=*/false, /*v1=*/false ))
707         return resultNodes;
708       if ( !nodes.empty() )
709       {
710         double paramSize = myLast[iE] - myFirst[iE];
711         double r         = myNormPar[iE] - prevNormPar;
712         eHelper.SetSubShape( myEdge[iE] );
713         eHelper.ToFixNodeParameters( true );
714         for ( size_t i = 0; i < nodes.size(); ++i )
715         {
716           double u = eHelper.GetNodeU( myEdge[iE], nodes[i], 0, &paramOK );
717           // paramSize is signed so orientation is taken into account
718           double normPar = prevNormPar + r * ( u - myFirst[iE] ) / paramSize;
719           u2node.insert( u2node.end(), make_pair( normPar, nodes[i] ));
720         }
721       }
722
723     } // loop on myEdges
724
725     if ( u2node.empty() ) return resultNodes;
726
727     // Add 2nd vertex node for a last EDGE
728     {
729       const SMDS_MeshNode* node;
730       if ( IsClosed() && theEdgeInd < 0 )
731         node = u2node.begin()->second;
732       else
733       {
734         node = VertexNode( iE );
735         while ( !node && iE > 0 )
736           node = VertexNode( --iE );
737         if ( !node )
738           return resultNodes;
739       }
740       if ( u2node.rbegin()->second == node &&
741            !fHelper.IsRealSeam  ( node->getshapeId() ) &&
742            !fHelper.IsDegenShape( node->getshapeId() ))
743         u2node.erase( --u2node.end() );
744
745       u2node.insert( u2node.end(), make_pair( 1., node ));
746     }
747
748     // Fill the result vector
749
750     if ( theEdgeInd < 0 &&
751          (int) u2node.size() != myNbPonits &&
752          (int) u2node.size() != NbPoints( /*update=*/true ))
753     {
754       u2node.clear();
755     }
756     resultNodes.reserve( u2node.size() );
757     map< double, const SMDS_MeshNode*>::iterator u2n = u2node.begin();
758     for ( ; u2n != u2node.end(); ++u2n )
759       resultNodes.push_back( u2n->second );
760   }
761   else
762   {
763     resultNodes.resize( myPoints.size() );
764     for ( size_t i = 0; i < myPoints.size(); ++i )
765       resultNodes[i] = myPoints[i].node;
766   }
767
768   return resultNodes;
769 }
770
771 //================================================================================
772 /*!
773  * \brief Return (unsorted) nodes of the i-th EDGE.
774  *        Nodes moved to other geometry by MergeNodes() are also returned.
775  *  \retval bool - is OK
776  */
777 //================================================================================
778
779 bool StdMeshers_FaceSide::GetEdgeNodes(size_t                        i,
780                                        vector<const SMDS_MeshNode*>& nodes,
781                                        bool                          inlude1stVertex,
782                                        bool                          inludeLastVertex) const
783 {
784   if ( i >= myEdge.size() )
785     return false;
786
787   SMESH_Mesh*     mesh = myProxyMesh->GetMesh();
788   SMESHDS_Mesh* meshDS = mesh->GetMeshDS();
789   SMESHDS_SubMesh*  sm = meshDS->MeshElements( myEdge[i] );
790
791   if ( inlude1stVertex )
792   {
793     if ( const SMDS_MeshNode* n0 = VertexNode( i ))
794       nodes.push_back( n0 );
795   }
796
797   if ( sm && ( sm->NbElements() > 0 || sm->NbNodes() > 0 ))
798   {
799     if ( mesh->HasModificationsToDiscard() ) // check nb of nodes on the EDGE sub-mesh
800     {
801       int iQuad    = sm->NbElements() ? sm->GetElements()->next()->IsQuadratic() : 0;
802       int nbExpect = sm->NbElements() - 1 + iQuad * sm->NbElements();
803       if ( nbExpect != sm->NbNodes() ) // some nodes are moved from the EDGE by MergeNodes()
804       {
805         // add nodes of all segments
806         typedef set< const SMDS_MeshNode* > TNodeSet;
807         TNodeSet sharedNodes;
808         SMDS_ElemIteratorPtr segIt = sm->GetElements();
809         while ( segIt->more() )
810         {
811           const SMDS_MeshElement* seg = segIt->next();
812           if ( seg->GetType() != SMDSAbs_Edge )
813             continue;
814           for ( int i = 0; i < 3-myIgnoreMediumNodes; ++i )
815           {
816             const SMDS_MeshNode* n = seg->GetNode( i );
817             if ( i == 2 ) // medium node
818             {
819               nodes.push_back( n );
820             }
821             else
822             {
823               pair<TNodeSet::iterator, bool> it2new = sharedNodes.insert( n );
824               if ( !it2new.second ) // n encounters twice == it's on EDGE, not on VERTEX
825               {
826                 nodes.push_back( n );
827                 sharedNodes.erase( it2new.first );
828               }
829             }
830           }
831         }
832       }
833     }
834     if ( nodes.size() < 2 ) // add nodes assigned to the EDGE
835     {
836       SMDS_NodeIteratorPtr nItr = sm->GetNodes();
837       while ( nItr->more() )
838       {
839         const SMDS_MeshNode* n = nItr->next();
840         if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( n, SMDSAbs_Edge ))
841           continue;
842         nodes.push_back( n );
843       }
844     }
845   } // if ( sm && sm->NbElements() > 0 )
846
847   if ( inludeLastVertex )
848   {
849     if ( const SMDS_MeshNode* n1 = VertexNode( i+1 ))
850       nodes.push_back( n1 );
851   }
852   return true;
853 }
854
855 //================================================================================
856 /*!
857  * \brief Return a node from the i-th VERTEX (count starts from zero)
858  *        Nodes moved to other geometry by MergeNodes() are also returned.
859  *  \param [in] i - the VERTEX index
860  *  \param [out] isMoved - returns \c true if the found node is moved by MergeNodes()
861  *  \return const SMDS_MeshNode* - the found node
862  */
863 //================================================================================
864
865 const SMDS_MeshNode* StdMeshers_FaceSide::VertexNode(std::size_t i, bool* isMoved) const
866 {
867   TopoDS_Vertex V = ( i >= myEdge.size() ) ? LastVertex() : FirstVertex(i);
868
869   const SMDS_MeshNode* n = SMESH_Algo::VertexNode( V, myProxyMesh->GetMeshDS() );
870
871   if ( !n && !myEdge.empty() && myProxyMesh->GetMesh()->HasModificationsToDiscard() )
872   {
873     size_t iE = ( i < myEdge.size() ) ? i : myEdge.size()-1;
874     SMESHDS_SubMesh* sm = myProxyMesh->GetMeshDS()->MeshElements( myEdgeID[ iE ]);
875
876     n = SMESH_Algo::VertexNode( V, sm, myProxyMesh->GetMesh(), /*checkV=*/false );
877
878     if (( !n ) &&
879         (( i > 0 && (int) i < NbEdges() ) || IsClosed() ))
880     {
881       iE = SMESH_MesherHelper::WrapIndex( int(i)-1, NbEdges() );
882       sm = myProxyMesh->GetMeshDS()->MeshElements( myEdgeID[ iE ]);
883       n  = SMESH_Algo::VertexNode( V, sm, myProxyMesh->GetMesh(), /*checkV=*/false );
884     }
885
886     if ( n && n->GetPosition()->GetDim() == 1 ) // check that n does not lie on an EDGE of myFace
887     {
888       TopoDS_Shape S = SMESH_MesherHelper::GetSubShapeByNode( n, myProxyMesh->GetMeshDS() );
889       if ( SMESH_MesherHelper::IsSubShape( S, myFace ))
890         n = 0; // VERTEX ignored by e.g. Composite Wire Discretization algo
891     }
892     if ( isMoved )
893       *isMoved = n;
894   }
895   return n;
896 }
897
898 //================================================================================
899 /*!
900  * \brief reverse order of vector elements
901   * \param vec - vector to reverse
902  */
903 //================================================================================
904
905 template <typename T > void reverse(vector<T> & vec)
906 {
907   std::reverse( vec.begin(), vec.end() );
908 }
909
910 //================================================================================
911 /*!
912  * \brief Change orientation of side geometry
913  */
914 //================================================================================
915
916 void StdMeshers_FaceSide::Reverse()
917 {
918   int nbEdges = myEdge.size();
919   for ( int i = nbEdges-1; i >= 0; --i ) {
920     std::swap( myFirst[i], myLast[i] );
921     if ( !myEdge[i].IsNull() )
922       myEdge[i].Reverse();
923     if ( i > 0 ) // at the first loop 1. is overwritten
924       myNormPar[i] = 1 - myNormPar[i-1];
925   }
926   if ( nbEdges > 1 ) {
927     reverse( myEdge );
928     reverse( myEdgeID );
929     reverse( myC2d );
930     //reverse( myC3dAdaptor );
931     reverse( myFirst );
932     reverse( myLast );
933     reverse( myNormPar );
934     reverse( myEdgeLength );
935     reverse( myIsUniform );
936   }
937   if ( nbEdges > 0 )
938   {
939     myNormPar[nbEdges-1]=1.;
940     if ( !myEdge[0].IsNull() )
941     {
942       for ( size_t i = 0; i < myEdge.size(); ++i )
943         reverseProxySubmesh( myEdge[i] );
944       myPoints.clear();
945       myFalsePoints.clear();
946     }
947     else
948     {
949       for ( size_t i = 0; i < myPoints.size(); ++i )
950       {
951         UVPtStruct & uvPt = myPoints[i];
952         uvPt.normParam = 1 - uvPt.normParam;
953         uvPt.x         = 1 - uvPt.x;
954         uvPt.y         = 1 - uvPt.y;
955       }
956       reverse( myPoints );
957
958       for ( size_t i = 0; i < myFalsePoints.size(); ++i )
959       {
960         UVPtStruct & uvPt = myFalsePoints[i];
961         uvPt.normParam = 1 - uvPt.normParam;
962         uvPt.x         = 1 - uvPt.x;
963         uvPt.y         = 1 - uvPt.y;
964       }
965       reverse( myFalsePoints );
966     }
967   }
968   for ( size_t i = 0; i < myEdge.size(); ++i )
969   {
970     if ( myEdge[i].IsNull() ) continue; // for a side on points only
971     double fp,lp;
972     Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],fp,lp);
973     if ( !C3d.IsNull() )
974       myC3dAdaptor[i].Load( C3d, fp,lp );
975   }
976 }
977
978 //=======================================================================
979 //function : SetIgnoreMediumNodes
980 //purpose  : Make ignore medium nodes
981 //=======================================================================
982
983 void StdMeshers_FaceSide::SetIgnoreMediumNodes(bool toIgnore)
984 {
985   if ( myIgnoreMediumNodes != toIgnore )
986   {
987     myIgnoreMediumNodes = toIgnore;
988
989     if ( !myPoints.empty() )
990     {
991       UVPtStructVec newPoints;
992       newPoints.reserve( myPoints.size()/2 + 1 );
993       for ( size_t i = 0; i < myPoints.size(); i += 2 )
994         newPoints.push_back( myPoints[i] );
995
996       myPoints.swap( newPoints );
997     }
998     else
999     {
1000       NbPoints( /*update=*/true );
1001     }
1002   }
1003 }
1004
1005 //=======================================================================
1006 //function : NbPoints
1007 //purpose  : Return nb nodes on edges and vertices (+1 to be == GetUVPtStruct().size() )
1008 //           Call it with update == true if mesh of this side can be recomputed
1009 //           since creation of this side
1010 //=======================================================================
1011
1012 int StdMeshers_FaceSide::NbPoints(const bool update) const
1013 {
1014   if ( !myPoints.empty() )
1015     return myPoints.size();
1016
1017   // if ( !myFalsePoints.empty() )
1018   //   return myFalsePoints.size();
1019
1020   if ( update && myEdge.size() > 0 )
1021   {
1022     StdMeshers_FaceSide* me = (StdMeshers_FaceSide*) this;
1023     me->myNbPonits = 0;
1024     me->myNbSegments = 0;
1025     me->myMissingVertexNodes = false;
1026
1027     vector<const SMDS_MeshNode*> nodes;
1028     for ( int i = 0; i < NbEdges(); ++i )
1029     {
1030       if ( const SMESHDS_SubMesh* sm = myProxyMesh->GetSubMesh( Edge(i) ))
1031       {
1032         if ( sm->NbNodes() == sm->NbElements()-1 || sm->NbElements() == 0 )
1033         {
1034           me->myNbPonits += sm->NbNodes();
1035           if ( myIgnoreMediumNodes && sm->IsQuadratic() )
1036             me->myNbPonits -= sm->NbElements();
1037         }
1038         else // nodes can be moved to other shapes by MergeNodes()
1039         {
1040           nodes.clear();
1041           GetEdgeNodes( i, nodes, /*v1=*/false, /*v2=*/false );
1042           me->myNbPonits += nodes.size();
1043         }
1044         me->myNbSegments += sm->NbElements();
1045       }
1046     }
1047
1048     SMESH_MesherHelper* helper = FaceHelper();
1049
1050     std::set< const SMDS_MeshNode* > vNodes;
1051     const int nbV = NbEdges() + !IsClosed();
1052     for ( int i = 0; i < nbV; ++i )
1053       if ( const SMDS_MeshNode* n = VertexNode( i ))
1054       {
1055         if ( !vNodes.insert( n ).second &&
1056              ( helper->IsRealSeam  ( n->getshapeId() ) ||
1057                helper->IsDegenShape( n->getshapeId() )))
1058           me->myNbPonits++;
1059       }
1060       else
1061       {
1062         me->myMissingVertexNodes = true;
1063       }
1064     me->myNbPonits += vNodes.size();
1065
1066     if ( IsClosed() )
1067       me->myNbPonits++; // closing node is repeated
1068   }
1069   return myNbPonits;
1070 }
1071
1072 //=======================================================================
1073 //function : NbSegments
1074 //purpose  : Return nb edges
1075 //           Call it with update == true if mesh of this side can be recomputed
1076 //           since creation of this side
1077 //=======================================================================
1078
1079 int StdMeshers_FaceSide::NbSegments(const bool update) const
1080 {
1081   return NbPoints( update ), myNbSegments;
1082 }
1083
1084 //================================================================================
1085 /*!
1086  * \brief Reverse UVPtStructVec if a proxy sub-mesh of E
1087  */
1088 //================================================================================
1089
1090 void StdMeshers_FaceSide::reverseProxySubmesh( const TopoDS_Edge& E )
1091 {
1092   if ( !myProxyMesh ) return;
1093   if ( const SMESH_ProxyMesh::SubMesh* sm = myProxyMesh->GetProxySubMesh( E ))
1094   {
1095     UVPtStructVec& edgeUVPtStruct = (UVPtStructVec& ) sm->GetUVPtStructVec();
1096     for ( size_t i = 0; i < edgeUVPtStruct.size(); ++i )
1097     {
1098       UVPtStruct & uvPt = edgeUVPtStruct[i];
1099       uvPt.normParam = 1 - uvPt.normParam;
1100       uvPt.x         = 1 - uvPt.x;
1101       uvPt.y         = 1 - uvPt.y;
1102     }
1103     reverse( edgeUVPtStruct );
1104   }
1105 }
1106
1107 //================================================================================
1108 /*!
1109  * \brief Show side features
1110  */
1111 //================================================================================
1112
1113 void StdMeshers_FaceSide::dump(const char* msg) const
1114 {
1115 #ifdef _DEBUG_
1116   if (msg) MESSAGE ( std::endl << msg );
1117   MESSAGE_BEGIN ("NB EDGES: "<< myEdge.size() );
1118   MESSAGE_ADD ( "nbPoints: "<<myNbPonits<<" vecSize: " << myPoints.size()<<" "<<myFalsePoints.size() );
1119   for ( size_t i = 0; i < myEdge.size(); ++i )
1120   {
1121     MESSAGE_ADD ( "\t"<<i+1 );
1122     MESSAGE_ADD ( "\tEDGE: " );
1123     if (myEdge[i].IsNull()) {
1124       MESSAGE_ADD ( "NULL" );
1125     }
1126     else {
1127       TopAbs::Print(myEdge[i].Orientation(),cout)<<" "<<myEdge[i].TShape().operator->()<<endl;
1128       MESSAGE_ADD ( "\tV1: " << TopExp::FirstVertex( myEdge[i], 1).TShape().operator->()
1129                     << "  V2: " << TopExp::LastVertex( myEdge[i], 1).TShape().operator->() );
1130     }
1131     MESSAGE_ADD ( "\tC2d: ");
1132
1133     if (myC2d[i].IsNull()) {
1134       MESSAGE_ADD ( "NULL" );
1135     }
1136     else {
1137       MESSAGE_ADD ( myC2d[i].operator->() );
1138     }
1139
1140     MESSAGE_ADD ( "\tF: "<<myFirst[i]<< " L: "<< myLast[i] );
1141     MESSAGE_END ( "\tnormPar: "<<myNormPar[i]<<endl );
1142   }
1143 #endif
1144 }
1145
1146 //================================================================================
1147 /*!
1148  * \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
1149   * \retval Adaptor2d_Curve2d* - 
1150  */
1151 //================================================================================
1152
1153 struct Adaptor2dCurve2d : public Adaptor2d_Curve2d
1154 {
1155   const StdMeshers_FaceSide* mySide;
1156   Adaptor2dCurve2d(const StdMeshers_FaceSide* faceSide):mySide(faceSide) {}
1157   gp_Pnt2d Value(const Standard_Real U) const { return mySide->Value2d( U ); }
1158   Standard_Real FirstParameter() const { return 0; }
1159   Standard_Real LastParameter() const { return 1; }
1160 };
1161
1162 Adaptor2d_Curve2d* StdMeshers_FaceSide::GetCurve2d() const
1163 {
1164   return new Adaptor2dCurve2d( this );
1165 }
1166
1167 //================================================================================
1168 /*!
1169  * \brief Creates a fully functional Adaptor_Curve
1170  */
1171 //================================================================================
1172
1173 BRepAdaptor_CompCurve* StdMeshers_FaceSide::GetCurve3d() const
1174 {
1175   if ( myEdge.empty() )
1176     return 0;
1177
1178   TopoDS_Wire aWire;
1179   BRep_Builder aBuilder;
1180   aBuilder.MakeWire(aWire);
1181   for ( size_t i = 0; i < myEdge.size(); ++i )
1182     aBuilder.Add( aWire, myEdge[i] );
1183
1184   if ( myEdge.size() == 2 && IsClosed() )
1185     aWire.Closed(true); // issue 0021141
1186
1187   return new BRepAdaptor_CompCurve( aWire );
1188 }
1189
1190 //================================================================================
1191 /*!
1192  * \brief Return 2D point by normalized parameter
1193   * \param U - normalized parameter value
1194   * \retval gp_Pnt2d - point
1195  */
1196 //================================================================================
1197
1198 gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
1199 {
1200   if ( !myC2d[0].IsNull() ) {
1201     int i = EdgeIndex( U );
1202     double prevU = i ? myNormPar[ i-1 ] : 0;
1203     double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
1204
1205     double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
1206     
1207     // check parametrization of curve
1208     if( !myIsUniform[i] )
1209     {
1210       double aLen3dU = r * myEdgeLength[i] * ( myFirst[i]>myLast[i] ? -1. : 1.);
1211       GCPnts_AbscissaPoint AbPnt
1212         ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
1213       if( AbPnt.IsDone() ) {
1214         par = AbPnt.Parameter();
1215       }
1216     }
1217     return myC2d[ i ]->Value(par);
1218
1219   }
1220   else if ( !myPoints.empty() )
1221   {
1222     int i = U * double( myPoints.size()-1 );
1223     while ( i > 0 && myPoints[ i ].normParam > U )
1224       --i;
1225     while ( i+1 < (int)myPoints.size() && myPoints[ i+1 ].normParam < U )
1226       ++i;
1227     double r = (( U - myPoints[ i ].normParam ) /
1228                 ( myPoints[ i+1 ].normParam - myPoints[ i ].normParam ));
1229     return ( myPoints[ i   ].UV() * ( 1 - r ) +
1230              myPoints[ i+1 ].UV() * r );
1231   }
1232   return myDefaultPnt2d;
1233 }
1234
1235 //================================================================================
1236 /*!
1237  * \brief Return XYZ by normalized parameter
1238   * \param U - normalized parameter value
1239   * \retval gp_Pnt - point
1240  */
1241 //================================================================================
1242
1243 gp_Pnt StdMeshers_FaceSide::Value3d(double U) const
1244 {
1245   int        i = EdgeIndex( U );
1246   double prevU = i ? myNormPar[ i-1 ] : 0;
1247   double     r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
1248
1249   double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
1250
1251   // check parametrization of curve
1252   if( !myIsUniform[i] )
1253   {
1254     double aLen3dU = r * myEdgeLength[i] * ( myFirst[i] > myLast[i] ? -1. : 1. );
1255     GCPnts_AbscissaPoint AbPnt
1256       ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
1257     if( AbPnt.IsDone() ) {
1258       par = AbPnt.Parameter();
1259     }
1260   }
1261   return myC3dAdaptor[ i ].Value(par);
1262 }
1263
1264 //================================================================================
1265 /*!
1266  * \brief Return wires of a face as StdMeshers_FaceSide's
1267  */
1268 //================================================================================
1269
1270 TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face&   theFace,
1271                                               SMESH_Mesh &         theMesh,
1272                                               const bool           theIgnoreMediumNodes,
1273                                               TError &             theError,
1274                                               SMESH_MesherHelper*  theFaceHelper,
1275                                               SMESH_ProxyMesh::Ptr theProxyMesh,
1276                                               const bool           theCheckVertexNodes)
1277 {
1278   SMESH_MesherHelper helper( theMesh );
1279   if ( theFaceHelper && theFaceHelper->GetSubShape() == theFace )
1280     helper.CopySubShapeInfo( *theFaceHelper );
1281
1282   list< TopoDS_Edge > edges, internalEdges;
1283   list< int > nbEdgesInWires;
1284   int nbWires = SMESH_Block::GetOrderedEdges (theFace, edges, nbEdgesInWires);
1285
1286   // split list of all edges into separate wires
1287   TSideVector wires( nbWires );
1288   list< int >::iterator nbE = nbEdgesInWires.begin();
1289   list< TopoDS_Edge >::iterator from = edges.begin(), to = from;
1290   for ( int iW = 0; iW < nbWires; ++iW, ++nbE )
1291   {
1292     std::advance( to, *nbE );
1293     if ( *nbE == 0 ) // Issue 0020676
1294     {
1295       --nbWires;
1296       --iW;
1297       wires.resize( nbWires );
1298       continue;
1299     }
1300     list< TopoDS_Edge > wireEdges( from, to );
1301     // assure that there is a node on the first vertex
1302     // as StdMeshers_FaceSide::GetUVPtStruct() requires
1303     if ( wireEdges.front().Orientation() != TopAbs_INTERNAL ) // Issue 0020676
1304     {
1305       if ( theCheckVertexNodes )
1306         while ( !SMESH_Algo::VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
1307                                          theMesh.GetMeshDS()))
1308         {
1309           wireEdges.splice(wireEdges.end(), wireEdges,
1310                            wireEdges.begin(), ++wireEdges.begin());
1311           if ( from->IsSame( wireEdges.front() )) {
1312             theError = TError
1313               ( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"No nodes on vertices"));
1314             return TSideVector(0);
1315           }
1316         }
1317     }
1318     else if ( *nbE > 1 ) // Issue 0020676 (Face_pb_netgen.brep) - several internal edges in a wire
1319     {
1320       internalEdges.splice( internalEdges.end(), wireEdges, ++wireEdges.begin(), wireEdges.end());
1321     }
1322
1323     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, wireEdges, &theMesh,
1324                                                          /*isForward=*/true, theIgnoreMediumNodes,
1325                                                          &helper, theProxyMesh );
1326     wires[ iW ] = StdMeshers_FaceSidePtr( wire );
1327     from = to;
1328   }
1329   while ( !internalEdges.empty() )
1330   {
1331     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, internalEdges.back(), &theMesh,
1332                                                          /*isForward=*/true, theIgnoreMediumNodes,
1333                                                          &helper, theProxyMesh );
1334     wires.push_back( StdMeshers_FaceSidePtr( wire ));
1335     internalEdges.pop_back();
1336   }
1337   return wires;
1338 }
1339
1340 //================================================================================
1341 /*!
1342  * \brief Return 1st vertex of the i-the edge
1343  */
1344 //================================================================================
1345
1346 TopoDS_Vertex StdMeshers_FaceSide::FirstVertex(int i) const
1347 {
1348   TopoDS_Vertex v;
1349   if ( i < NbEdges() )
1350   {
1351     v = myEdge[i].Orientation() <= TopAbs_REVERSED ? // FORWARD || REVERSED
1352         TopExp::FirstVertex( myEdge[i], 1 )        :
1353         TopoDS::Vertex( TopoDS_Iterator( myEdge[i] ).Value() );
1354   }
1355   return v;
1356 }
1357
1358 //================================================================================
1359 /*!
1360  * \brief Return last vertex of the i-the edge
1361  */
1362 //================================================================================
1363
1364 TopoDS_Vertex StdMeshers_FaceSide::LastVertex(int i) const
1365 {
1366   TopoDS_Vertex v;
1367   if ( i < NbEdges() )
1368   {
1369     const TopoDS_Edge& edge = i<0 ? myEdge[ NbEdges() + i ] : myEdge[i];
1370     if ( edge.Orientation() <= TopAbs_REVERSED ) // FORWARD || REVERSED
1371       v = TopExp::LastVertex( edge, 1 );
1372     else
1373       for ( TopoDS_Iterator vIt( edge ); vIt.More(); vIt.Next() )
1374         v = TopoDS::Vertex( vIt.Value() );
1375   }
1376   return v;
1377 }
1378
1379 //================================================================================
1380 /*!
1381  * \brief Return \c true if the chain of EDGEs is closed
1382  */
1383 //================================================================================
1384
1385 bool StdMeshers_FaceSide::IsClosed() const
1386 {
1387   return myEdge.empty() ? false : FirstVertex().IsSame( LastVertex() );
1388 }
1389
1390 //================================================================================
1391 /*!
1392  * \brief Return a helper initialized with the FACE
1393  */
1394 //================================================================================
1395
1396 SMESH_MesherHelper* StdMeshers_FaceSide::FaceHelper() const
1397 {
1398   StdMeshers_FaceSide* me = const_cast< StdMeshers_FaceSide* >( this );
1399   if ( !myHelper && myProxyMesh )
1400   {
1401     me->myHelper = new SMESH_MesherHelper( *myProxyMesh->GetMesh() );
1402     me->myHelper->SetSubShape( myFace );
1403   }
1404   return me->myHelper;
1405 }