]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_FaceSide.cxx
Salome HOME
32f798dee5dfd9677dbaade884f861961b691119
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.cxx
1 //  Copyright (C) 2007-2010  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 //  SMESH SMESH : implementaion of SMESH idl descriptions
24 // File      : StdMeshers_FaceSide.hxx
25 // Created   : Wed Jan 31 18:41:25 2007
26 // Author    : Edward AGAPOV (eap)
27 // Module    : SMESH
28 //
29 #include "StdMeshers_FaceSide.hxx"
30
31 #include "SMDS_EdgePosition.hxx"
32 #include "SMDS_MeshNode.hxx"
33 #include "SMESHDS_Mesh.hxx"
34 #include "SMESHDS_SubMesh.hxx"
35 #include "SMESH_Algo.hxx"
36 #include "SMESH_Mesh.hxx"
37 #include "SMESH_MesherHelper.hxx"
38 #include "SMESH_ComputeError.hxx"
39 #include "SMESH_Block.hxx"
40
41 #include <Adaptor2d_Curve2d.hxx>
42 #include <BRepAdaptor_CompCurve.hxx>
43 #include <BRepAdaptor_Curve.hxx>
44 #include <BRep_Builder.hxx>
45 #include <BRep_Tool.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 <GCPnts_AbscissaPoint.hxx>
54 #include <Geom2dAdaptor_Curve.hxx>
55
56 #include <map>
57
58 #include "utilities.h"
59
60 //================================================================================
61 /*!
62  * \brief Constructor of a side of one edge
63   * \param theFace - the face
64   * \param theEdge - the edge
65  */
66 //================================================================================
67
68 StdMeshers_FaceSide::StdMeshers_FaceSide(const TopoDS_Face& theFace,
69                                          const TopoDS_Edge& theEdge,
70                                          SMESH_Mesh*        theMesh,
71                                          const bool         theIsForward,
72                                          const bool         theIgnoreMediumNodes)
73 {
74   list<TopoDS_Edge> edges(1,theEdge);
75   *this = StdMeshers_FaceSide( theFace, edges, theMesh, theIsForward, theIgnoreMediumNodes );
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 {
92   int nbEdges = theEdges.size();
93   myEdge.resize( nbEdges );
94   myC2d.resize( nbEdges );
95   myC3dAdaptor.resize( nbEdges );
96   myFirst.resize( nbEdges );
97   myLast.resize( nbEdges );
98   myNormPar.resize( nbEdges );
99   myEdgeLength.resize( nbEdges );
100   myIsUniform.resize( nbEdges );
101   myLength = 0;
102   myNbPonits = myNbSegments = 0;
103   myMesh = theMesh;
104   myMissingVertexNodes = false;
105   myIgnoreMediumNodes = theIgnoreMediumNodes;
106   myDefaultPnt2d = gp_Pnt2d( 1e+100, 1e+100 );
107   if ( nbEdges == 0 ) return;
108
109   SMESHDS_Mesh* meshDS = theMesh->GetMeshDS();
110
111   int nbDegen = 0;
112   list<TopoDS_Edge>::iterator edge = theEdges.begin();
113   TopoDS_Iterator vExp;
114   for ( int index = 0; edge != theEdges.end(); ++index, ++edge )
115   {
116     int i = theIsForward ? index : nbEdges - index - 1;
117     myEdgeLength[i] = SMESH_Algo::EdgeLength( *edge );
118     if ( myEdgeLength[i] < DBL_MIN ) nbDegen++;
119     myLength += myEdgeLength[i];
120     myEdge[i] = *edge;
121     if ( !theIsForward ) myEdge[i].Reverse();
122
123     if ( theFace.IsNull() )
124       BRep_Tool::Range( *edge, myFirst[i], myLast[i] );
125     else
126       myC2d[i] = BRep_Tool::CurveOnSurface( *edge, theFace, myFirst[i], myLast[i] );
127     if ( myEdge[i].Orientation() == TopAbs_REVERSED )
128       std::swap( myFirst[i], myLast[i] );
129
130     if ( SMESHDS_SubMesh* sm = meshDS->MeshElements( *edge )) {
131       int nbN = sm->NbNodes();
132       if ( theIgnoreMediumNodes ) {
133         SMDS_ElemIteratorPtr elemIt = sm->GetElements();
134         if ( elemIt->more() && elemIt->next()->IsQuadratic() )
135           nbN -= sm->NbElements();
136       }
137       myNbPonits += nbN;
138       myNbSegments += sm->NbElements();
139     }
140     // TopExp::FirstVertex() and TopExp::LastVertex() return NULL from INTERNAL edge
141     vExp.Initialize( *edge );
142     if ( vExp.Value().Orientation() == TopAbs_REVERSED ) vExp.Next();
143     if ( SMESH_Algo::VertexNode( TopoDS::Vertex( vExp.Value()), meshDS ))
144       myNbPonits += 1; // for the first end
145     else
146       myMissingVertexNodes = true;
147
148     // check if edge has non-uniform parametrization (issue 0020705)
149     if ( !myC2d[i].IsNull() )
150     {
151       Geom2dAdaptor_Curve A2dC( myC2d[i] );
152       double p2 = myFirst[i]+(myLast[i]-myFirst[i])/2., p4 = myFirst[i]+(myLast[i]-myFirst[i])/4.;
153       double d2 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p2 );
154       double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
155       //cout<<"len = "<<len<<"  d2 = "<<d2<<"  fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
156       myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
157       if ( !myIsUniform[i] )
158       {
159         double fp,lp;
160         TopLoc_Location L;
161         Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],L,fp,lp);
162         myC3dAdaptor[i].Load( C3d, fp,lp );
163       }
164     }
165   }
166   vExp.Initialize( theEdges.back() );
167   if ( vExp.Value().Orientation() != TopAbs_REVERSED ) vExp.Next();
168   if ( vExp.More() )
169   {
170     if ( SMESH_Algo::VertexNode( TopoDS::Vertex( vExp.Value()), meshDS ))
171       myNbPonits++; // for the last end
172     else
173       myMissingVertexNodes = true;
174   }
175   if ( nbEdges > 1 && myLength > DBL_MIN ) {
176     const double degenNormLen = 1.e-5;
177     double totLength = myLength;
178     if ( nbDegen )
179       totLength += myLength * degenNormLen * nbDegen;
180     double prevNormPar = 0;
181     for ( int i = 0; i < nbEdges; ++i ) {
182       if ( myEdgeLength[ i ] < DBL_MIN )
183         myEdgeLength[ i ] = myLength * degenNormLen;
184       myNormPar[ i ] = prevNormPar + myEdgeLength[i]/totLength;
185       prevNormPar = myNormPar[ i ];
186     }
187   }
188   myNormPar[nbEdges-1] = 1.;
189   //dump();
190 }
191
192 //================================================================================
193 /*!
194  * \brief Constructor of a side for vertex using data from other FaceSide
195   * \param theVertex - the vertex
196   * \param theSide - the side
197  */
198 //================================================================================
199
200 StdMeshers_FaceSide::StdMeshers_FaceSide(const SMDS_MeshNode* theNode,
201                                          const gp_Pnt2d thePnt2d,
202                                          const StdMeshers_FaceSide* theSide)
203 {
204   myC2d.resize(1);
205   myLength = 0;
206   myMesh = theSide->GetMesh();
207   myDefaultPnt2d = thePnt2d;
208
209   myPoints = theSide->GetUVPtStruct();
210   myNbPonits = myNbSegments = myPoints.size();
211   std::vector<uvPtStruct>::iterator it = myPoints.begin();
212   for(; it!=myPoints.end(); it++) {
213     (*it).u = thePnt2d.X();
214     (*it).v = thePnt2d.Y();
215     (*it).y = 0.0;
216     (*it).node = theNode;
217   }
218 }
219
220 //================================================================================
221 /*!
222  * \brief Return info on nodes on the side
223   * \retval UVPtStruct* - array of data structures
224  */
225 //================================================================================
226
227 const vector<UVPtStruct>& StdMeshers_FaceSide::GetUVPtStruct(bool   isXConst,
228                                                              double constValue) const
229 {
230   if ( myPoints.empty() ) {
231
232     if ( NbEdges() == 0 ) return myPoints;
233
234     SMESHDS_Mesh* meshDS = myMesh->GetMeshDS();
235     SMESH_MesherHelper helper(*myMesh);
236     bool paramOK;
237
238     // sort nodes of all edges putting them into a map
239
240     map< double, const SMDS_MeshNode*> u2node;
241     //int nbOnDegen = 0;
242     for ( int i = 0; i < myEdge.size(); ++i ) {
243       // put 1st vertex node
244       TopoDS_Vertex VV[2]; // TopExp::FirstVertex() returns NULL for INTERNAL edge
245       for ( TopoDS_Iterator vIt(myEdge[i]); vIt.More(); vIt.Next() )
246         VV[ VV[0].IsNull() ? 0 : 1 ] = TopoDS::Vertex(vIt.Value());
247       if ( VV[0].Orientation() == TopAbs_REVERSED ) std::swap ( VV[0], VV[1] );
248       const SMDS_MeshNode* node = SMESH_Algo::VertexNode( VV[0], meshDS );
249       double prevNormPar = ( i == 0 ? 0 : myNormPar[ i-1 ]); // normalized param
250       if ( node ) { // internal nodes may be missing
251         u2node.insert( make_pair( prevNormPar, node ));
252       }
253       else if ( i == 0 ) {
254         MESSAGE(" NO NODE on VERTEX" );
255         return myPoints;
256       }
257
258       // put 2nd vertex node for a last edge
259       if ( i+1 == myEdge.size() ) {
260         node = SMESH_Algo::VertexNode( VV[1], meshDS );
261         if ( !node ) {
262           MESSAGE(" NO NODE on VERTEX" );
263           return myPoints;
264         }
265         u2node.insert( make_pair( 1., node ));
266       }
267
268       // put internal nodes
269       SMESHDS_SubMesh* sm = meshDS->MeshElements( myEdge[i] );
270       if ( !sm ) continue;
271       SMDS_NodeIteratorPtr nItr = sm->GetNodes();
272       double paramSize = myLast[i] - myFirst[i];
273       double r = myNormPar[i] - prevNormPar;
274       while ( nItr->more() )
275       {
276         const SMDS_MeshNode* node = nItr->next();
277         if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
278           continue;
279         double u = helper.GetNodeU( myEdge[i], node, &paramOK );
280
281         // paramSize is signed so orientation is taken into account
282         double normPar = prevNormPar + r * ( u - myFirst[i] ) / paramSize;
283         if(!myIsUniform[i])
284         {
285           double aLenU = GCPnts_AbscissaPoint::Length
286             ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), myFirst[i], u );
287           if ( myEdgeLength[i] > aLenU ) // nonregression test "3D_mesh_NETGEN/G6"
288             normPar = prevNormPar + r*aLenU/myEdgeLength[i];
289         }
290 #ifdef _DEBUG_
291         if ( normPar > 1 || normPar < 0) {
292           dump("DEBUG");
293           MESSAGE ( "WRONG normPar: "<<normPar<< " prevNormPar="<<prevNormPar
294                     << " u="<<u << " myFirst[i]="<<myFirst[i]<< " myLast[i]="<<myLast[i]
295                     << " paramSize="<<paramSize );
296         }
297 #endif
298         u2node.insert( make_pair( normPar, node ));
299       }
300     }
301     if ( u2node.size() != myNbPonits ) {
302       MESSAGE("Wrong node parameters on edges, u2node.size():"
303               <<u2node.size()<<" !=  myNbPonits:"<<myNbPonits);
304       return myPoints;
305     }
306
307     // fill array of UVPtStruct
308
309     vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myPoints );
310     points->resize( myNbPonits );
311
312     int EdgeIndex = 0;
313     double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
314     map< double, const SMDS_MeshNode*>::iterator u_node = u2node.begin();
315     for (int i = 0 ; u_node != u2node.end(); ++u_node, ++i ) {
316       UVPtStruct & uvPt = (*points)[i];
317       uvPt.node = u_node->second;
318       uvPt.x = uvPt.y = uvPt.normParam = u_node->first;
319       if ( isXConst ) uvPt.x = constValue;
320       else            uvPt.y = constValue;
321       if ( myNormPar[ EdgeIndex ] < uvPt.normParam ) {
322         prevNormPar = myNormPar[ EdgeIndex ];
323         ++EdgeIndex;
324 #ifdef _DEBUG_
325         if ( EdgeIndex >= myEdge.size() ) {
326           dump("DEBUG");
327           MESSAGE ( "WRONg EdgeIndex " << 1+EdgeIndex
328                     << " myNormPar.size()="<<myNormPar.size()
329                     << " myNormPar["<< EdgeIndex<<"]="<< myNormPar[ EdgeIndex ]
330                     << " uvPt.normParam="<<uvPt.normParam );
331         }
332 #endif
333         paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
334       }
335       const SMDS_EdgePosition* epos =
336         dynamic_cast<const SMDS_EdgePosition*>(uvPt.node->GetPosition().get());
337       if ( epos ) {
338         uvPt.param = epos->GetUParameter();
339       }
340       else {
341         double r = ( uvPt.normParam - prevNormPar )/ paramSize;
342 //         uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
343         uvPt.param = ( r > 0.5 ? myLast[EdgeIndex] : myFirst[EdgeIndex] );
344       }
345       if ( !myC2d[ EdgeIndex ].IsNull() ) {
346         gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
347         uvPt.u = p.X();
348         uvPt.v = p.Y();
349       }
350       else {
351         uvPt.u = uvPt.v = 1e+100;
352       }
353     }
354   }
355   return myPoints;
356 }
357
358 //================================================================================
359 /*!
360  * \brief Falsificate info on nodes
361   * \param nbSeg - nb of segments on the side
362   * \retval UVPtStruct* - array of data structures
363  */
364 //================================================================================
365
366 const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int    nbSeg,
367                                                                   bool   isXConst,
368                                                                   double constValue) const
369 {
370   if ( myFalsePoints.empty() ) {
371
372     if ( NbEdges() == 0 ) return myFalsePoints;
373
374     vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myFalsePoints );
375     points->resize( nbSeg+1 );
376
377     int EdgeIndex = 0;
378     double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
379     for (int i = 0 ; i < myFalsePoints.size(); ++i ) {
380       double normPar = double(i) / double(nbSeg);
381       UVPtStruct & uvPt = (*points)[i];
382       uvPt.node = 0;
383       uvPt.x = uvPt.y = uvPt.param = uvPt.normParam = normPar;
384       if ( isXConst ) uvPt.x = constValue;
385       else            uvPt.y = constValue;
386       if ( myNormPar[ EdgeIndex ] < normPar ) {
387         prevNormPar = myNormPar[ EdgeIndex ];
388         ++EdgeIndex;
389         paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
390       }
391       double r = ( normPar - prevNormPar )/ paramSize;
392       uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
393       if ( !myC2d[ EdgeIndex ].IsNull() ) {
394         gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
395         uvPt.u = p.X();
396         uvPt.v = p.Y();
397       }
398       else {
399         uvPt.u = uvPt.v = 1e+100;
400       }
401     }
402   }
403   return myFalsePoints;
404 }
405 // gp_Pnt StdMeshers_FaceSide::Value(double U) const
406 // {
407 // }
408
409 //================================================================================
410 /*!
411  * \brief reverse order of vector elements
412   * \param vec - vector to reverse
413  */
414 //================================================================================
415
416 template <typename T > void reverse(vector<T> & vec)
417 {
418   for ( int f=0, r=vec.size()-1; f < r; ++f, --r )
419     std::swap( vec[f], vec[r] );
420 }
421
422 //================================================================================
423 /*!
424  * \brief Change orientation of side geometry
425  */
426 //================================================================================
427
428 void StdMeshers_FaceSide::Reverse()
429 {
430   int nbEdges = myEdge.size();
431   for ( int i = nbEdges-1; i >= 0; --i ) {
432     std::swap( myFirst[i], myLast[i] );
433     myEdge[i].Reverse();
434     if ( i > 0 ) // at the first loop 1. is overwritten
435       myNormPar[i] = 1 - myNormPar[i-1];
436   }
437   if ( nbEdges > 1 ) {
438     reverse( myEdge );
439     reverse( myC2d );
440     reverse( myFirst );
441     reverse( myLast );
442     reverse( myNormPar );
443   }
444   myNormPar[nbEdges-1]=1.;
445   myPoints.clear();
446   myFalsePoints.clear();
447 }
448
449 //================================================================================
450 /*!
451  * \brief Show side features
452  */
453 //================================================================================
454
455 void StdMeshers_FaceSide::dump(const char* msg) const
456 {
457 #ifdef _DEBUG_
458   if (msg) MESSAGE ( std::endl << msg );
459   MESSAGE_BEGIN ("NB EDGES: "<< myEdge.size() );
460   MESSAGE_ADD ( "nbPoints: "<<myNbPonits<<" vecSize: " << myPoints.size()<<" "<<myFalsePoints.size() );
461   for ( int i=0; i<myEdge.size(); ++i)
462   {
463     MESSAGE_ADD ( "\t"<<i+1 );
464     MESSAGE_ADD ( "\tEDGE: " );
465     if (myEdge[i].IsNull()) {
466       MESSAGE_ADD ( "NULL" );
467     }
468     else {
469       TopAbs::Print(myEdge[i].Orientation(),cout)<<" "<<myEdge[i].TShape().operator->()<<endl;
470       MESSAGE_ADD ( "\tV1: " << TopExp::FirstVertex( myEdge[i], 1).TShape().operator->()
471                  << "  V2: " << TopExp::LastVertex( myEdge[i], 1).TShape().operator->() );
472     }
473     MESSAGE_ADD ( "\tC2d: ");
474     
475     if (myC2d[i].IsNull()) {
476       MESSAGE_ADD ( "NULL" );
477     }
478     else {
479       MESSAGE_ADD ( myC2d[i].operator->() );
480     }
481       
482     MESSAGE_ADD ( "\tF: "<<myFirst[i]<< " L: "<< myLast[i] );
483     MESSAGE_END ( "\tnormPar: "<<myNormPar[i]<<endl );
484   }
485 #endif
486 }
487
488 //================================================================================
489 /*!
490  * \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
491   * \retval Adaptor2d_Curve2d* - 
492  */
493 //================================================================================
494
495 struct Adaptor2dCurve2d : public Adaptor2d_Curve2d
496 {
497   const StdMeshers_FaceSide* mySide;
498   Adaptor2dCurve2d(const StdMeshers_FaceSide* faceSide):mySide(faceSide) {}
499   gp_Pnt2d Value(const Standard_Real U) const { return mySide->Value2d( U ); }
500   Standard_Real FirstParameter() const { return 0; }
501   Standard_Real LastParameter() const { return 1; }
502 };
503
504 Adaptor2d_Curve2d* StdMeshers_FaceSide::GetCurve2d() const
505 {
506   return new Adaptor2dCurve2d( this );
507 }
508
509 //================================================================================
510 /*!
511  * \brief Creates a fully functional Adaptor_Curve
512  */
513 //================================================================================
514
515 BRepAdaptor_CompCurve* StdMeshers_FaceSide::GetCurve3d() const
516 {
517   if ( myEdge.empty() )
518     return 0;
519
520 //   if ( myEdge.size() == 1 )
521 //     return new BRepAdaptor_Curve( myEdge[0] );
522
523   TopoDS_Wire aWire;
524   BRep_Builder aBuilder;
525   aBuilder.MakeWire(aWire);
526   for ( int i=0; i<myEdge.size(); ++i )
527     aBuilder.Add( aWire, myEdge[i] );
528   return new BRepAdaptor_CompCurve( aWire );
529 }
530
531
532 //================================================================================
533 /*!
534  * \brief Return 2D point by normalized parameter
535   * \param U - normalized parameter value
536   * \retval gp_Pnt2d - point
537  */
538 //================================================================================
539
540 gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
541 {
542   if ( !myC2d[0].IsNull() ) {
543     int i = EdgeIndex( U );
544     double prevU = i ? myNormPar[ i-1 ] : 0;
545     double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
546
547     double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
548     
549     // check parametrization of curve
550     if( !myIsUniform[i] )
551     {
552       double aLen3dU = r * myEdgeLength[i] * ( myFirst[i]>myLast[i] ? -1. : 1.);
553       GCPnts_AbscissaPoint AbPnt
554         ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
555       if( AbPnt.IsDone() ) {
556         par = AbPnt.Parameter();
557       }
558     }
559     return myC2d[ i ]->Value(par);
560
561   }
562   return myDefaultPnt2d;
563 }
564
565 //================================================================================
566 /*!
567  * \brief Return wires of a face as StdMeshers_FaceSide's
568  */
569 //================================================================================
570
571 TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face& theFace,
572                                               SMESH_Mesh &       theMesh,
573                                               const bool         theIgnoreMediumNodes,
574                                               TError &           theError)
575 {
576   TopoDS_Vertex V1;
577   list< TopoDS_Edge > edges;
578   list< int > nbEdgesInWires;
579   int nbWires = SMESH_Block::GetOrderedEdges (theFace, V1, edges, nbEdgesInWires);
580
581   // split list of all edges into separate wires
582   TSideVector wires( nbWires );
583   list< int >::iterator nbE = nbEdgesInWires.begin();
584   list< TopoDS_Edge >::iterator from, to;
585   from = to = edges.begin();
586   for ( int iW = 0; iW < nbWires; ++iW )
587   {
588     std::advance( to, *nbE++ );
589     if ( *nbE == 0 ) // Issue 0020676
590     {
591       --nbWires;
592       --iW;
593       wires.resize( nbWires );
594       continue;
595     }
596     list< TopoDS_Edge > wireEdges( from, to );
597     // assure that there is a node on the first vertex
598     // as StdMeshers_FaceSide::GetUVPtStruct() requires
599     if ( wireEdges.front().Orientation() != TopAbs_INTERNAL ) // Issue 0020676
600       while ( !SMESH_Algo::VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
601                                        theMesh.GetMeshDS()))
602       {
603         wireEdges.splice(wireEdges.end(), wireEdges,
604                          wireEdges.begin(), ++wireEdges.begin());
605         if ( from->IsSame( wireEdges.front() )) {
606           theError = TError
607             ( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"No nodes on vertices"));
608           return TSideVector(0);
609         }
610       }
611     const bool isForward = true;
612     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, wireEdges, &theMesh,
613                                                          isForward, theIgnoreMediumNodes);
614     wires[ iW ] = StdMeshers_FaceSidePtr( wire );
615     from = to;
616   }
617   return wires;
618 }
619
620 //================================================================================
621 /*!
622  * \brief Return 1st vertex of the i-the edge
623  */
624 //================================================================================
625
626 TopoDS_Vertex StdMeshers_FaceSide::FirstVertex(int i) const
627 {
628   TopoDS_Vertex v;
629   if ( i < NbEdges() )
630   {
631     v = myEdge[i].Orientation() <= TopAbs_REVERSED ? // FORWARD || REVERSED
632         TopExp::FirstVertex( myEdge[i], 1 )        :
633         TopoDS::Vertex( TopoDS_Iterator( myEdge[i] ).Value() );
634   }
635   return v;
636 }
637
638 //================================================================================
639 /*!
640  * \brief Return last vertex of the i-the edge
641  */
642 //================================================================================
643
644 TopoDS_Vertex StdMeshers_FaceSide::LastVertex(int i) const
645 {
646   TopoDS_Vertex v;
647   if ( i < NbEdges() )
648   {
649     const TopoDS_Edge& edge = i<0 ? myEdge[ NbEdges() + i ] : myEdge[i];
650     if ( edge.Orientation() <= TopAbs_REVERSED ) // FORWARD || REVERSED
651       v = TopExp::LastVertex( edge, 1 );
652     else
653       for ( TopoDS_Iterator vIt( edge ); vIt.More(); vIt.Next() )
654         v = TopoDS::Vertex( vIt.Value() );
655   }
656   return v;
657 }
658