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