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