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