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