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