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