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