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