Salome HOME
Merge from V6_5_BR 05/06/2012
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.cxx
1 // Copyright (C) 2007-2012  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                                 std::min( myFirst[i], myLast[i] ),
152                                 std::max( myFirst[i], myLast[i] ));
153       double p2 = myFirst[i]+(myLast[i]-myFirst[i])/2., p4 = myFirst[i]+(myLast[i]-myFirst[i])/4.;
154       double d2 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p2 );
155       double d4 = GCPnts_AbscissaPoint::Length( A2dC, myFirst[i], p4 );
156       //cout<<"len = "<<len<<"  d2 = "<<d2<<"  fabs(2*d2/len-1.0) = "<<fabs(2*d2/len-1.0)<<endl;
157       myIsUniform[i] = !( fabs(2*d2/myEdgeLength[i]-1.0) > 0.01 || fabs(2*d4/d2-1.0) > 0.01 );
158       if ( !myIsUniform[i] )
159       {
160         double fp,lp;
161         TopLoc_Location L;
162         Handle(Geom_Curve) C3d = BRep_Tool::Curve(myEdge[i],L,fp,lp);
163         myC3dAdaptor[i].Load( C3d, fp,lp );
164       }
165     }
166   }
167   vExp.Initialize( theEdges.back() );
168   if ( vExp.Value().Orientation() != TopAbs_REVERSED ) vExp.Next();
169   if ( vExp.More() )
170   {
171     if ( SMESH_Algo::VertexNode( TopoDS::Vertex( vExp.Value()), meshDS ))
172       myNbPonits++; // for the last end
173     else
174       myMissingVertexNodes = true;
175   }
176   if ( nbEdges > 1 && myLength > DBL_MIN ) {
177     const double degenNormLen = 1.e-5;
178     double totLength = myLength;
179     if ( nbDegen )
180       totLength += myLength * degenNormLen * nbDegen;
181     double prevNormPar = 0;
182     for ( int i = 0; i < nbEdges; ++i ) {
183       if ( myEdgeLength[ i ] < DBL_MIN )
184         myEdgeLength[ i ] = myLength * degenNormLen;
185       myNormPar[ i ] = prevNormPar + myEdgeLength[i]/totLength;
186       prevNormPar = myNormPar[ i ];
187     }
188   }
189   myNormPar[nbEdges-1] = 1.;
190   //dump();
191 }
192
193 //================================================================================
194 /*!
195  * \brief Constructor of a side for vertex using data from other FaceSide
196   * \param theVertex - the vertex
197   * \param theSide - the side
198  */
199 //================================================================================
200
201 StdMeshers_FaceSide::StdMeshers_FaceSide(const SMDS_MeshNode* theNode,
202                                          const gp_Pnt2d thePnt2d,
203                                          const StdMeshers_FaceSide* theSide)
204 {
205   myC2d.resize(1);
206   myLength = 0;
207   myMesh = theSide->GetMesh();
208   myDefaultPnt2d = thePnt2d;
209
210   myPoints = theSide->GetUVPtStruct();
211   myNbPonits = myNbSegments = myPoints.size();
212   std::vector<uvPtStruct>::iterator it = myPoints.begin();
213   for(; it!=myPoints.end(); it++) {
214     (*it).u = thePnt2d.X();
215     (*it).v = thePnt2d.Y();
216     (*it).y = 0.0;
217     (*it).node = theNode;
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     SMESH_MesherHelper helper(*myMesh);
237     bool paramOK;
238
239     // sort nodes of all edges putting them into a map
240
241     map< double, const SMDS_MeshNode*> u2node;
242     //int nbOnDegen = 0;
243     for ( int i = 0; i < myEdge.size(); ++i )
244     {
245       // Put 1st vertex node of a current edge
246       TopoDS_Vertex VV[2]; // TopExp::FirstVertex() returns NULL for INTERNAL edge
247       VV[0] = SMESH_MesherHelper::IthVertex( 0, myEdge[i]);
248       VV[1] = SMESH_MesherHelper::IthVertex( 1, myEdge[i]);
249       const SMDS_MeshNode* node = SMESH_Algo::VertexNode( VV[0], meshDS );
250       double prevNormPar = ( i == 0 ? 0 : myNormPar[ i-1 ]); // normalized param
251       if ( node ) { // nodes on internal vertices may be missing
252         u2node.insert( u2node.end(), make_pair( prevNormPar, node ));
253       }
254       else if ( i == 0 ) {
255         MESSAGE(" NO NODE on VERTEX" );
256         return myPoints;
257       }
258
259       // Put internal nodes
260       if ( SMESHDS_SubMesh* sm = meshDS->MeshElements( myEdge[i] ))
261       {
262         vector< pair< double, const SMDS_MeshNode*> > u2nodeVec;
263         u2nodeVec.reserve( sm->NbNodes() );
264         SMDS_NodeIteratorPtr nItr = sm->GetNodes();
265         double paramSize = myLast[i] - myFirst[i];
266         double r = myNormPar[i] - prevNormPar;
267         if ( !myIsUniform[i] )
268           while ( nItr->more() )
269           {
270             const SMDS_MeshNode* node = nItr->next();
271             if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
272               continue;
273             double u = helper.GetNodeU( myEdge[i], node, 0, &paramOK );
274             double aLenU = GCPnts_AbscissaPoint::Length
275               ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), myFirst[i], u );
276             if ( myEdgeLength[i] < aLenU ) // nonregression test "3D_mesh_NETGEN/G6"
277             {
278               u2nodeVec.clear();
279               break;
280             }
281             double normPar = prevNormPar + r*aLenU/myEdgeLength[i];
282             u2nodeVec.push_back( make_pair( normPar, node ));
283           }
284         nItr = sm->GetNodes();
285         if ( u2nodeVec.empty() )
286           while ( nItr->more() )
287           {
288             const SMDS_MeshNode* node = nItr->next();
289             if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
290               continue;
291             double u = helper.GetNodeU( myEdge[i], node, 0, &paramOK );
292
293             // paramSize is signed so orientation is taken into account
294             double normPar = prevNormPar + r * ( u - myFirst[i] ) / paramSize;
295             u2nodeVec.push_back( make_pair( normPar, node ));
296           }
297         for ( size_t j = 0; j < u2nodeVec.size(); ++j )
298         u2node.insert( u2node.end(), u2nodeVec[j] );
299       }
300
301       // Put 2nd vertex node for a last edge
302       if ( i+1 == myEdge.size() ) {
303         node = SMESH_Algo::VertexNode( VV[1], meshDS );
304         if ( !node ) {
305           MESSAGE(" NO NODE on VERTEX" );
306           return myPoints;
307         }
308         u2node.insert( u2node.end(), make_pair( 1., node ));
309       }
310     }
311     if ( u2node.size() != myNbPonits ) {
312       MESSAGE("Wrong node parameters on edges, u2node.size():"
313               <<u2node.size()<<" !=  myNbPonits:"<<myNbPonits);
314       return myPoints;
315     }
316
317     // fill array of UVPtStruct
318
319     vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myPoints );
320     points->resize( myNbPonits );
321
322     int EdgeIndex = 0;
323     double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
324     map< double, const SMDS_MeshNode*>::iterator u_node = u2node.begin();
325     for (int i = 0 ; u_node != u2node.end(); ++u_node, ++i ) {
326       UVPtStruct & uvPt = (*points)[i];
327       uvPt.node = u_node->second;
328       uvPt.x = uvPt.y = uvPt.normParam = u_node->first;
329       if ( isXConst ) uvPt.x = constValue;
330       else            uvPt.y = constValue;
331       const SMDS_EdgePosition* epos =
332         dynamic_cast<const SMDS_EdgePosition*>(uvPt.node->GetPosition());
333       if (( myNormPar[ EdgeIndex ] < uvPt.normParam ) ||
334           ( epos && uvPt.node->getshapeId() != myEdgeID[ EdgeIndex ])) // for myMissingVertexNodes
335       {
336         prevNormPar = myNormPar[ EdgeIndex ];
337         ++EdgeIndex;
338 #ifdef _DEBUG_
339         if ( EdgeIndex >= myEdge.size() ) {
340           dump("DEBUG");
341           MESSAGE ( "WRONg EdgeIndex " << 1+EdgeIndex
342                     << " myNormPar.size()="<<myNormPar.size()
343                     << " myNormPar["<< EdgeIndex<<"]="<< myNormPar[ EdgeIndex ]
344                     << " uvPt.normParam="<<uvPt.normParam );
345         }
346 #endif
347         paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
348       }
349       if ( epos ) {
350         uvPt.param = epos->GetUParameter();
351       }
352       else {
353         double r = ( uvPt.normParam - prevNormPar )/ paramSize;
354 //         uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
355         uvPt.param = ( r > 0.5 ? myLast[EdgeIndex] : myFirst[EdgeIndex] );
356       }
357       if ( !myC2d[ EdgeIndex ].IsNull() ) {
358         gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
359         uvPt.u = p.X();
360         uvPt.v = p.Y();
361       }
362       else {
363         uvPt.u = uvPt.v = 1e+100;
364       }
365     }
366   }
367   return myPoints;
368 }
369
370 //================================================================================
371 /*!
372  * \brief Falsificate info on nodes
373   * \param nbSeg - nb of segments on the side
374   * \retval UVPtStruct* - array of data structures
375  */
376 //================================================================================
377
378 const vector<UVPtStruct>& StdMeshers_FaceSide::SimulateUVPtStruct(int    nbSeg,
379                                                                   bool   isXConst,
380                                                                   double constValue) const
381 {
382   if ( myFalsePoints.empty() ) {
383
384     if ( NbEdges() == 0 ) return myFalsePoints;
385
386     vector<uvPtStruct>* points = const_cast<vector<uvPtStruct>*>( &myFalsePoints );
387     points->resize( nbSeg+1 );
388
389     int EdgeIndex = 0;
390     double prevNormPar = 0, paramSize = myNormPar[ EdgeIndex ];
391     for (int i = 0 ; i < myFalsePoints.size(); ++i ) {
392       double normPar = double(i) / double(nbSeg);
393       UVPtStruct & uvPt = (*points)[i];
394       uvPt.node = 0;
395       uvPt.x = uvPt.y = uvPt.param = uvPt.normParam = normPar;
396       if ( isXConst ) uvPt.x = constValue;
397       else            uvPt.y = constValue;
398       if ( myNormPar[ EdgeIndex ] < normPar ) {
399         prevNormPar = myNormPar[ EdgeIndex ];
400         ++EdgeIndex;
401         paramSize = myNormPar[ EdgeIndex ] - prevNormPar;
402       }
403       double r = ( normPar - prevNormPar )/ paramSize;
404       uvPt.param = myFirst[EdgeIndex] * ( 1 - r ) + myLast[EdgeIndex] * r;
405       if ( !myC2d[ EdgeIndex ].IsNull() ) {
406         gp_Pnt2d p = myC2d[ EdgeIndex ]->Value( uvPt.param );
407         uvPt.u = p.X();
408         uvPt.v = p.Y();
409       }
410       else {
411         uvPt.u = uvPt.v = 1e+100;
412       }
413     }
414   }
415   return myFalsePoints;
416 }
417
418 //=======================================================================
419 //function : GetOrderedNodes
420 //purpose  : Return nodes in the order they encounter while walking along the side
421 //=======================================================================
422
423 std::vector<const SMDS_MeshNode*> StdMeshers_FaceSide::GetOrderedNodes() const
424 {
425   vector<const SMDS_MeshNode*> resultNodes;
426   if ( myPoints.empty() )
427   {
428     if ( NbEdges() == 0 ) return resultNodes;
429
430     SMESHDS_Mesh* meshDS = myMesh->GetMeshDS();
431     SMESH_MesherHelper helper(*myMesh);
432     bool paramOK;
433
434     // Sort nodes of all edges putting them into a map
435
436     map< double, const SMDS_MeshNode*> u2node;
437     for ( int i = 0; i < myEdge.size(); ++i )
438     {
439       // Put 1st vertex node of a current edge
440       TopoDS_Vertex VV[2]; // TopExp::FirstVertex() returns NULL for INTERNAL edge
441       VV[0] = SMESH_MesherHelper::IthVertex( 0, myEdge[i]);
442       VV[1] = SMESH_MesherHelper::IthVertex( 1, myEdge[i]);
443       const SMDS_MeshNode* node = SMESH_Algo::VertexNode( VV[0], meshDS );
444       double prevNormPar = ( i == 0 ? 0 : myNormPar[ i-1 ]); // normalized param
445       if ( node ) { // nodes on internal vertices may be missing
446         u2node.insert( make_pair( prevNormPar, node ));
447       }
448       else if ( i == 0 ) {
449         MESSAGE(" NO NODE on VERTEX" );
450         return resultNodes;
451       }
452
453       // Put internal nodes
454       if ( SMESHDS_SubMesh* sm = meshDS->MeshElements( myEdge[i] ))
455       {
456         SMDS_NodeIteratorPtr nItr = sm->GetNodes();
457         double paramSize = myLast[i] - myFirst[i];
458         double r = myNormPar[i] - prevNormPar;
459         while ( nItr->more() )
460         {
461           const SMDS_MeshNode* node = nItr->next();
462           if ( myIgnoreMediumNodes && SMESH_MeshEditor::IsMedium( node, SMDSAbs_Edge ))
463             continue;
464           double u = helper.GetNodeU( myEdge[i], node, 0, &paramOK );
465
466           // paramSize is signed so orientation is taken into account
467           double normPar = prevNormPar + r * ( u - myFirst[i] ) / paramSize;
468           u2node.insert( u2node.end(), make_pair( normPar, node ));
469         }
470       }
471
472       // Put 2nd vertex node for a last edge
473       if ( i+1 == myEdge.size() ) {
474         node = SMESH_Algo::VertexNode( VV[1], meshDS );
475         if ( !node ) {
476           return resultNodes;
477         }
478         u2node.insert( u2node.end(), make_pair( 1., node ));
479       }
480     }
481
482     // Fill the result vector
483
484     if ( u2node.size() == myNbPonits )
485     {
486       resultNodes.reserve( u2node.size() );
487       map< double, const SMDS_MeshNode*>::iterator u2n = u2node.begin();
488       for ( ; u2n != u2node.end(); ++u2n )
489         resultNodes.push_back( u2n->second );
490     }
491   }
492   else
493   {
494     resultNodes.resize( myPoints.size() );
495     for ( size_t i = 0; i < myPoints.size(); ++i )
496       resultNodes[i] = myPoints[i].node;
497   }
498
499   return resultNodes;
500 }
501
502 //================================================================================
503 /*!
504  * \brief reverse order of vector elements
505   * \param vec - vector to reverse
506  */
507 //================================================================================
508
509 template <typename T > void reverse(vector<T> & vec)
510 {
511   for ( int f=0, r=vec.size()-1; f < r; ++f, --r )
512     std::swap( vec[f], vec[r] );
513 }
514
515 //================================================================================
516 /*!
517  * \brief Change orientation of side geometry
518  */
519 //================================================================================
520
521 void StdMeshers_FaceSide::Reverse()
522 {
523   int nbEdges = myEdge.size();
524   for ( int i = nbEdges-1; i >= 0; --i ) {
525     std::swap( myFirst[i], myLast[i] );
526     myEdge[i].Reverse();
527     if ( i > 0 ) // at the first loop 1. is overwritten
528       myNormPar[i] = 1 - myNormPar[i-1];
529   }
530   if ( nbEdges > 1 ) {
531     reverse( myEdge );
532     reverse( myEdgeID );
533     reverse( myC2d );
534     reverse( myC3dAdaptor );
535     reverse( myFirst );
536     reverse( myLast );
537     reverse( myNormPar );
538     reverse( myEdgeLength );
539     reverse( myIsUniform );
540   }
541   if ( nbEdges > 0 )
542   {
543     myNormPar[nbEdges-1]=1.;
544     myPoints.clear();
545     myFalsePoints.clear();
546   }
547 }
548
549 //================================================================================
550 /*!
551  * \brief Show side features
552  */
553 //================================================================================
554
555 void StdMeshers_FaceSide::dump(const char* msg) const
556 {
557 #ifdef _DEBUG_
558   if (msg) MESSAGE ( std::endl << msg );
559   MESSAGE_BEGIN ("NB EDGES: "<< myEdge.size() );
560   MESSAGE_ADD ( "nbPoints: "<<myNbPonits<<" vecSize: " << myPoints.size()<<" "<<myFalsePoints.size() );
561   for ( int i=0; i<myEdge.size(); ++i)
562   {
563     MESSAGE_ADD ( "\t"<<i+1 );
564     MESSAGE_ADD ( "\tEDGE: " );
565     if (myEdge[i].IsNull()) {
566       MESSAGE_ADD ( "NULL" );
567     }
568     else {
569       TopAbs::Print(myEdge[i].Orientation(),cout)<<" "<<myEdge[i].TShape().operator->()<<endl;
570       MESSAGE_ADD ( "\tV1: " << TopExp::FirstVertex( myEdge[i], 1).TShape().operator->()
571                  << "  V2: " << TopExp::LastVertex( myEdge[i], 1).TShape().operator->() );
572     }
573     MESSAGE_ADD ( "\tC2d: ");
574     
575     if (myC2d[i].IsNull()) {
576       MESSAGE_ADD ( "NULL" );
577     }
578     else {
579       MESSAGE_ADD ( myC2d[i].operator->() );
580     }
581       
582     MESSAGE_ADD ( "\tF: "<<myFirst[i]<< " L: "<< myLast[i] );
583     MESSAGE_END ( "\tnormPar: "<<myNormPar[i]<<endl );
584   }
585 #endif
586 }
587
588 //================================================================================
589 /*!
590  * \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
591   * \retval Adaptor2d_Curve2d* - 
592  */
593 //================================================================================
594
595 struct Adaptor2dCurve2d : public Adaptor2d_Curve2d
596 {
597   const StdMeshers_FaceSide* mySide;
598   Adaptor2dCurve2d(const StdMeshers_FaceSide* faceSide):mySide(faceSide) {}
599   gp_Pnt2d Value(const Standard_Real U) const { return mySide->Value2d( U ); }
600   Standard_Real FirstParameter() const { return 0; }
601   Standard_Real LastParameter() const { return 1; }
602 };
603
604 Adaptor2d_Curve2d* StdMeshers_FaceSide::GetCurve2d() const
605 {
606   return new Adaptor2dCurve2d( this );
607 }
608
609 //================================================================================
610 /*!
611  * \brief Creates a fully functional Adaptor_Curve
612  */
613 //================================================================================
614
615 BRepAdaptor_CompCurve* StdMeshers_FaceSide::GetCurve3d() const
616 {
617   if ( myEdge.empty() )
618     return 0;
619
620   TopoDS_Wire aWire;
621   BRep_Builder aBuilder;
622   aBuilder.MakeWire(aWire);
623   for ( int i=0; i<myEdge.size(); ++i )
624     aBuilder.Add( aWire, myEdge[i] );
625
626   if ( myEdge.size() == 2 && FirstVertex().IsSame( LastVertex() ))
627     aWire.Closed(true); // issue 0021141
628
629   return new BRepAdaptor_CompCurve( aWire );
630 }
631
632
633 //================================================================================
634 /*!
635  * \brief Return 2D point by normalized parameter
636   * \param U - normalized parameter value
637   * \retval gp_Pnt2d - point
638  */
639 //================================================================================
640
641 gp_Pnt2d StdMeshers_FaceSide::Value2d(double U) const
642 {
643   if ( !myC2d[0].IsNull() ) {
644     int i = EdgeIndex( U );
645     double prevU = i ? myNormPar[ i-1 ] : 0;
646     double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
647
648     double par = myFirst[i] * ( 1 - r ) + myLast[i] * r;
649     
650     // check parametrization of curve
651     if( !myIsUniform[i] )
652     {
653       double aLen3dU = r * myEdgeLength[i] * ( myFirst[i]>myLast[i] ? -1. : 1.);
654       GCPnts_AbscissaPoint AbPnt
655         ( const_cast<GeomAdaptor_Curve&>( myC3dAdaptor[i]), aLen3dU, myFirst[i] );
656       if( AbPnt.IsDone() ) {
657         par = AbPnt.Parameter();
658       }
659     }
660     return myC2d[ i ]->Value(par);
661
662   }
663   return myDefaultPnt2d;
664 }
665
666 //================================================================================
667 /*!
668  * \brief Return wires of a face as StdMeshers_FaceSide's
669  */
670 //================================================================================
671
672 TSideVector StdMeshers_FaceSide::GetFaceWires(const TopoDS_Face& theFace,
673                                               SMESH_Mesh &       theMesh,
674                                               const bool         theIgnoreMediumNodes,
675                                               TError &           theError)
676 {
677   TopoDS_Vertex V1;
678   list< TopoDS_Edge > edges, internalEdges;
679   list< int > nbEdgesInWires;
680   int nbWires = SMESH_Block::GetOrderedEdges (theFace, V1, edges, nbEdgesInWires);
681
682   // split list of all edges into separate wires
683   TSideVector wires( nbWires );
684   list< int >::iterator nbE = nbEdgesInWires.begin();
685   list< TopoDS_Edge >::iterator from = edges.begin(), to = from;
686   for ( int iW = 0; iW < nbWires; ++iW, ++nbE )
687   {
688     std::advance( to, *nbE );
689     if ( *nbE == 0 ) // Issue 0020676
690     {
691       --nbWires;
692       --iW;
693       wires.resize( nbWires );
694       continue;
695     }
696     list< TopoDS_Edge > wireEdges( from, to );
697     // assure that there is a node on the first vertex
698     // as StdMeshers_FaceSide::GetUVPtStruct() requires
699     if ( wireEdges.front().Orientation() != TopAbs_INTERNAL ) // Issue 0020676
700     {
701       while ( !SMESH_Algo::VertexNode( TopExp::FirstVertex( wireEdges.front(), true),
702                                        theMesh.GetMeshDS()))
703       {
704         wireEdges.splice(wireEdges.end(), wireEdges,
705                          wireEdges.begin(), ++wireEdges.begin());
706         if ( from->IsSame( wireEdges.front() )) {
707           theError = TError
708             ( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"No nodes on vertices"));
709           return TSideVector(0);
710         }
711       }
712     }
713     else if ( *nbE > 1 ) // Issue 0020676 (Face_pb_netgen.brep) - several internal edges in a wire
714     {
715       internalEdges.splice( internalEdges.end(), wireEdges, ++wireEdges.begin(), wireEdges.end());
716     }
717
718     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, wireEdges, &theMesh,
719                                                          /*isForward=*/true, theIgnoreMediumNodes);
720     wires[ iW ] = StdMeshers_FaceSidePtr( wire );
721     from = to;
722   }
723   while ( !internalEdges.empty() )
724   {
725     StdMeshers_FaceSide* wire = new StdMeshers_FaceSide( theFace, internalEdges.back(), &theMesh,
726                                                          /*isForward=*/true, theIgnoreMediumNodes);
727     wires.push_back( StdMeshers_FaceSidePtr( wire ));
728     internalEdges.pop_back();
729   }
730   return wires;
731 }
732
733 //================================================================================
734 /*!
735  * \brief Return 1st vertex of the i-the edge
736  */
737 //================================================================================
738
739 TopoDS_Vertex StdMeshers_FaceSide::FirstVertex(int i) const
740 {
741   TopoDS_Vertex v;
742   if ( i < NbEdges() )
743   {
744     v = myEdge[i].Orientation() <= TopAbs_REVERSED ? // FORWARD || REVERSED
745         TopExp::FirstVertex( myEdge[i], 1 )        :
746         TopoDS::Vertex( TopoDS_Iterator( myEdge[i] ).Value() );
747   }
748   return v;
749 }
750
751 //================================================================================
752 /*!
753  * \brief Return last vertex of the i-the edge
754  */
755 //================================================================================
756
757 TopoDS_Vertex StdMeshers_FaceSide::LastVertex(int i) const
758 {
759   TopoDS_Vertex v;
760   if ( i < NbEdges() )
761   {
762     const TopoDS_Edge& edge = i<0 ? myEdge[ NbEdges() + i ] : myEdge[i];
763     if ( edge.Orientation() <= TopAbs_REVERSED ) // FORWARD || REVERSED
764       v = TopExp::LastVertex( edge, 1 );
765     else
766       for ( TopoDS_Iterator vIt( edge ); vIt.More(); vIt.Next() )
767         v = TopoDS::Vertex( vIt.Value() );
768   }
769   return v;
770 }
771