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