Salome HOME
IPAL54027: Projection algo is very long on a face with many edges
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.hxx
1 // Copyright (C) 2007-2016  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, or (at your option) any later version.
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 #ifndef StdMeshers_FaceSide_HeaderFile
29 #define StdMeshers_FaceSide_HeaderFile
30
31 #include "SMESH_StdMeshers.hxx"
32
33 #include "SMESH_ProxyMesh.hxx"
34
35 #include <Geom2d_Curve.hxx>
36 #include <GeomAdaptor_Curve.hxx>
37 #include <TopoDS_Edge.hxx>
38 #include <TopoDS_Face.hxx>
39 #include <TopoDS_Vertex.hxx>
40 #include <gp_Pnt2d.hxx>
41
42 #include <vector>
43 #include <list>
44 #include <boost/shared_ptr.hpp>
45
46 class Adaptor2d_Curve2d;
47 class Adaptor3d_Curve;
48 class BRepAdaptor_CompCurve;
49 class SMDS_MeshNode;
50 class SMESH_Mesh;
51 class SMESH_MesherHelper;
52 class StdMeshers_FaceSide;
53 struct SMESH_ComputeError;
54
55 typedef boost::shared_ptr< SMESH_ComputeError >  TError;
56 typedef boost::shared_ptr< StdMeshers_FaceSide > StdMeshers_FaceSidePtr;
57 typedef std::vector< StdMeshers_FaceSidePtr >    TSideVector;
58
59 //================================================================================
60 /*!
61  * \brief Represents a side of a quasi quadrilateral face.
62  * It can be composed of several edges. Gives access to geometry and 1D mesh of a side.
63  */
64 //================================================================================
65
66 class STDMESHERS_EXPORT StdMeshers_FaceSide
67 {
68 public:
69
70   enum { ALL_EDGES = -1, LAST_EDGE = -1 }; //!< constants
71
72   /*!
73    * \brief Wrap one edge
74    */
75   StdMeshers_FaceSide(const TopoDS_Face&   theFace,
76                       const TopoDS_Edge&   theEdge,
77                       SMESH_Mesh*          theMesh,
78                       const bool           theIsForward,
79                       const bool           theIgnoreMediumNodes,
80                       SMESH_MesherHelper*  theFaceHelper = NULL,
81                       SMESH_ProxyMesh::Ptr theProxyMesh = SMESH_ProxyMesh::Ptr());
82   /*!
83    * \brief Wrap several edges. Edges must be properly ordered and oriented.
84    */
85   StdMeshers_FaceSide(const TopoDS_Face&      theFace,
86                       std::list<TopoDS_Edge>& theEdges,
87                       SMESH_Mesh*             theMesh,
88                       const bool              theIsForward,
89                       const bool              theIgnoreMediumNodes,
90                       SMESH_MesherHelper*     theFaceHelper = NULL,
91                       SMESH_ProxyMesh::Ptr    theProxyMesh = SMESH_ProxyMesh::Ptr());
92   /*!
93    * \brief Simulate a side from a vertex using data from other FaceSide
94    */
95   StdMeshers_FaceSide(const StdMeshers_FaceSide*  theSide,
96                       const SMDS_MeshNode*        theNode,
97                       const gp_Pnt2d*             thePnt2d1,
98                       const gp_Pnt2d*             thePnt2d2 = NULL,
99                       const Handle(Geom2d_Curve)& theC2d = NULL,
100                       const double                theUFirst = 0.,
101                       const double                theULast = 1.);
102   /*!
103    * \brief Create a side from an UVPtStructVec
104    */
105   StdMeshers_FaceSide(UVPtStructVec&     theSideNodes,
106                       const TopoDS_Face& theFace = TopoDS_Face(),
107                       const TopoDS_Edge& theEdge = TopoDS_Edge(),
108                       SMESH_Mesh*        theMesh = 0);
109
110   ~StdMeshers_FaceSide();
111
112   // static "consrtuctors"
113   static StdMeshers_FaceSidePtr New(const TopoDS_Face&   Face,
114                                     const TopoDS_Edge&   Edge,
115                                     SMESH_Mesh*          Mesh,
116                                     const bool           IsForward,
117                                     const bool           IgnoreMediumNodes,
118                                     SMESH_MesherHelper*  FaceHelper = NULL,
119                                     SMESH_ProxyMesh::Ptr ProxyMesh = SMESH_ProxyMesh::Ptr())
120   { return StdMeshers_FaceSidePtr
121       ( new StdMeshers_FaceSide( Face,Edge,Mesh,IsForward,IgnoreMediumNodes,FaceHelper,ProxyMesh ));
122   }
123   static StdMeshers_FaceSidePtr New (const TopoDS_Face&      Face,
124                                      std::list<TopoDS_Edge>& Edges,
125                                      SMESH_Mesh*             Mesh,
126                                      const bool              IsForward,
127                                      const bool              IgnoreMediumNodes,
128                                      SMESH_MesherHelper*     FaceHelper = NULL,
129                                      SMESH_ProxyMesh::Ptr    ProxyMesh = SMESH_ProxyMesh::Ptr())
130   { return StdMeshers_FaceSidePtr
131       ( new StdMeshers_FaceSide( Face,Edges,Mesh,IsForward,IgnoreMediumNodes,FaceHelper,ProxyMesh ));
132   }
133   static StdMeshers_FaceSidePtr New (const StdMeshers_FaceSide*  Side,
134                                      const SMDS_MeshNode*        Node,
135                                      const gp_Pnt2d*             Pnt2d1,
136                                      const gp_Pnt2d*             Pnt2d2 = NULL,
137                                      const Handle(Geom2d_Curve)& C2d = NULL,
138                                      const double                UFirst = 0.,
139                                      const double                ULast = 1.)
140   { return StdMeshers_FaceSidePtr
141       ( new StdMeshers_FaceSide( Side,Node,Pnt2d1,Pnt2d2,C2d,UFirst,ULast ));
142   }
143   static StdMeshers_FaceSidePtr New (UVPtStructVec&     theSideNodes,
144                                      const TopoDS_Face& theFace = TopoDS_Face())
145   {
146     return StdMeshers_FaceSidePtr( new StdMeshers_FaceSide( theSideNodes, theFace ));
147   }
148
149   /*!
150    * \brief Return wires of a face as StdMeshers_FaceSide's
151    */
152   static TSideVector GetFaceWires(const TopoDS_Face&   theFace,
153                                   SMESH_Mesh &         theMesh,
154                                   const bool           theIgnoreMediumNodes,
155                                   TError &             theError,
156                                   SMESH_MesherHelper*  theFaceHelper = NULL,
157                                   SMESH_ProxyMesh::Ptr theProxyMesh = SMESH_ProxyMesh::Ptr(),
158                                   const bool           theCheckVertexNodes = true);
159   /*!
160    * \brief Change orientation of side geometry
161    */
162   void Reverse();
163   /*!
164    * \brief Make ignore medium nodes
165    */
166   void SetIgnoreMediumNodes(bool toIgnore);
167
168   /*!
169    * \brief Return nb nodes on edges and vertices (+1 to be == GetUVPtStruct().size() ).
170    *        Call it with update == true if mesh of this side can be recomputed
171    *        since creation of this side
172    */
173   int NbPoints(const bool update = false) const;
174   /*!
175    * \brief Return nb edges
176    *        Call it with update == true if mesh of this side can be recomputed
177    *        since creation of this side
178    */
179   int NbSegments(const bool update = false) const;
180   /*!
181    * \brief Return mesh
182    */
183   SMESH_Mesh* GetMesh() const { return myProxyMesh->GetMesh(); }
184   /*!
185    * \brief Return true if there are vertices without nodes
186    */
187   bool MissVertexNode() const { return myMissingVertexNodes; }
188
189   /*!
190    * \brief Return detailed data on nodes
191     * \param isXConst - true if normalized parameter X is constant
192     * \param constValue - constant parameter value
193     *
194     * Missing nodes are allowed only on internal vertices.
195     * For a closed side, the 1st point repeats at end
196    */
197   const UVPtStructVec& GetUVPtStruct( bool isXConst = 0, double constValue = 0 ) const;
198   /*!
199    * \brief Simulates detailed data on nodes
200     * \param isXConst - true if normalized parameter X is constant
201     * \param constValue - constant parameter value
202    */
203   const UVPtStructVec& SimulateUVPtStruct(int    nbSeg,
204                                           bool   isXConst   = 0,
205                                           double constValue = 0) const;
206   /*!
207    * \brief Return nodes in the order they encounter while walking along
208    *  the while side or a specified EDGE. For a closed side, the 1st point repeats at end.
209    *  \param iE - index of the EDGE. Default is "all EDGEs".
210    */
211   std::vector<const SMDS_MeshNode*> GetOrderedNodes( int iE = ALL_EDGES ) const;
212
213   /*!
214    * \brief Return nodes of the i-th EDGE.
215    *        Nodes moved to other geometry by MergeNodes() are also returned.
216    * \retval bool - is OK
217    */
218   bool GetEdgeNodes(const size_t                       i,
219                     std::vector<const SMDS_MeshNode*>& nodes,
220                     bool                               inlude1stVertex=true,
221                     bool                               inludeLastVertex=true) const;
222
223   /*!
224    * \brief Return a node from the i-th VERTEX (count starts from zero)
225    *        Nodes moved to other geometry by MergeNodes() are also returned.
226    */
227   const SMDS_MeshNode* VertexNode(std::size_t i, bool* isMoved = 0) const;
228
229   /*
230    * \brief Return edge and parameter on edge by normalized parameter
231    */
232   inline double Parameter(double U, TopoDS_Edge & edge) const;
233   /*
234    * \brief Return edge ID and parameter on edge by normalized parameter
235    */
236   inline double Parameter(double U, int & edgeID) const;
237   /*!
238    * \brief Return UV by normalized parameter
239    */
240   gp_Pnt2d Value2d(double U) const;
241   /*!
242    * \brief Return XYZ by normalized parameter
243    */
244   gp_Pnt   Value3d(double U) const;
245   /*!
246    * \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
247    */
248   Adaptor2d_Curve2d* GetCurve2d() const;
249   /*!
250    * \brief Creates a fully functional Adaptor_Curve
251    */
252   BRepAdaptor_CompCurve* GetCurve3d() const;
253   /*!
254    * \brief Return nb of wrapped edges
255    */
256   int NbEdges() const { return myEdge.size(); }
257   /*!
258    * \brief Return i-th edge (count starts from zero)
259    */
260   const TopoDS_Edge& Edge(int i) const { return myEdge[i]; }
261   /*!
262    * \brief Return all edges
263    */
264   const std::vector<TopoDS_Edge>& Edges() const { return myEdge; }
265   /*!
266    * \brief Return the FACE
267    */
268   const TopoDS_Face& Face() const { return myFace; }
269   /*!
270    * \brief Return 1st vertex of the i-th edge (count starts from zero)
271    */
272   TopoDS_Vertex FirstVertex(int i=0) const;
273   /*!
274    * \brief Return last vertex of the i-th edge (count starts from zero)
275    */
276   TopoDS_Vertex LastVertex(int i = LAST_EDGE) const;
277   /*!
278    * \brief Return \c true if the chain of EDGEs is closed
279    */
280   bool IsClosed() const;
281   /*!
282    * \brief Return side length
283    */
284   double Length() const { return myLength; }
285   /*!
286    * \brief Return edge index corresponding to normalized parameter
287    */
288   inline int EdgeIndex( double U ) const;
289
290   void dump(const char* msg=0) const;
291   
292   /*!
293    * \brief Return ID of i-th wrapped edge (count starts from zero)
294    */
295   inline int EdgeID(int i) const;
296   /*!
297    * \brief Return p-curve of i-th wrapped edge (count starts from zero)
298    */
299   inline Handle(Geom2d_Curve) Curve2d(int i) const;
300   /*!
301    * \brief Return first normalized parameter of the i-th edge (count starts from zero)
302    */
303   inline double FirstParameter(int i) const;
304   /*!
305    * \brief Return last normalized parameter of the i-th edge (count starts from zero)
306    */
307   inline double LastParameter(int i) const;
308   /*!
309    * \brief Return first parameter of the i-th edge (count starts from zero).
310    *        EDGE orientation is taken into account
311    */
312   inline double FirstU(int i) const;
313   /*!
314    * \brief Return last parameter of the i-th edge (count starts from zero).
315    *        EDGE orientation is taken into account
316    */
317   inline double LastU(int i) const;
318   /*!
319    * \brief Return length of i-th wrapped edge (count starts from zero)
320    */
321   inline double EdgeLength(int i) const;
322   /*!
323    * \brief Return orientation of i-th wrapped edge (count starts from zero)
324    */
325   inline bool IsReversed(int i) const;
326   /*!
327    * \brief Return a helper initialized with the FACE
328    */
329   SMESH_MesherHelper* FaceHelper() const;
330
331 protected:
332
333   void reverseProxySubmesh( const TopoDS_Edge& E );
334
335   // DON't FORGET to update Reverse() when adding one more vector!
336   TopoDS_Face                       myFace;
337   std::vector<uvPtStruct>           myPoints, myFalsePoints;
338   std::vector<TopoDS_Edge>          myEdge;
339   std::vector<int>                  myEdgeID;
340   std::vector<Handle(Geom2d_Curve)> myC2d;
341   std::vector<GeomAdaptor_Curve>    myC3dAdaptor;
342   std::vector<double>               myFirst, myLast;
343   std::vector<double>               myNormPar;
344   std::vector<double>               myEdgeLength;
345   std::vector<int>                  myIsUniform;
346   double                            myLength;
347   int                               myNbPonits, myNbSegments;
348   SMESH_ProxyMesh::Ptr              myProxyMesh;
349   bool                              myMissingVertexNodes, myIgnoreMediumNodes;
350   gp_Pnt2d                          myDefaultPnt2d;
351   SMESH_MesherHelper*               myHelper;
352 };
353
354
355 //================================================================================
356 /*!
357  * \brief Return edge index corresponding to normalized parameter
358   * \param U - the parameter
359   * \retval int - index
360  */
361 //================================================================================
362
363 inline int StdMeshers_FaceSide::EdgeIndex( double U ) const
364 {
365   int i = myNormPar.size() - 1;
366   while ( i > 0 && U < myNormPar[ i-1 ] ) --i;
367   return i;
368 }
369
370 //================================================================================
371 /*!
372  * \brief Return an edge and parameter on the edge by a normalized parameter
373   * \param U - normalized parameter
374   * \retval double - pameter on a curve
375   * \ warning The returned parameter can be inaccurate if the edge is non-uniformly
376   *           parametrized. Use Value2d() to get a precise point on the edge
377  */
378 //================================================================================
379
380 inline double StdMeshers_FaceSide::Parameter(double U, TopoDS_Edge & edge) const
381 {
382   int i = EdgeIndex( U );
383   edge = myEdge[ i ];
384   double prevU = i ? myNormPar[ i-1 ] : 0;
385   double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
386   return myFirst[i] * ( 1 - r ) + myLast[i] * r;
387 }
388
389 //================================================================================
390 /*!
391  * \brief Return an edge ID and parameter on the edge by a normalized parameter
392   * \param U - normalized parameter
393   * \retval double - pameter on a curve
394   * \ warning The returned parameter can be inaccurate if the edge is non-uniformly
395   *           parametrized. Use Value2d() to get a precise point on the edge
396  */
397 //================================================================================
398
399 inline double StdMeshers_FaceSide::Parameter(double U, int & edgeID) const
400 {
401   int i = EdgeIndex( U );
402   edgeID = myEdgeID[ i ];
403   double prevU = i ? myNormPar[ i-1 ] : 0;
404   double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
405   return myFirst[i] * ( 1 - r ) + myLast[i] * r;
406 }
407
408 //================================================================================
409 /*!
410  * \brief Return first normalized parameter of the i-th edge
411  */
412 //================================================================================
413
414 inline double StdMeshers_FaceSide::FirstParameter(int i) const
415 {
416   return i==0 ? 0. : i<(int)myNormPar.size() ? myNormPar[i-1] : 1.;
417 }
418
419 //================================================================================
420 /*!
421  * \brief Return ast normalized parameter of the i-th edge
422  */
423 //================================================================================
424
425 inline double StdMeshers_FaceSide::LastParameter(int i) const
426 {
427   return i < (int)myNormPar.size() ? myNormPar[i] : 1;
428 }
429
430 //================================================================================
431 /*!
432  * \brief Return first parameter of the i-th edge
433  */
434 //================================================================================
435
436 inline double StdMeshers_FaceSide::FirstU(int i) const
437 {
438   return myFirst[ i % myFirst.size() ];
439 }
440
441 //================================================================================
442 /*!
443  * \brief Return last parameter of the i-th edge
444  */
445 //================================================================================
446
447 inline double StdMeshers_FaceSide::LastU(int i) const
448 {
449   return myLast[ i % myLast.size() ];
450 }
451
452 //================================================================================
453   /*!
454    * \brief Return ID of i-th wrapped edge (count starts from zero)
455    */
456 //================================================================================
457
458 inline int StdMeshers_FaceSide::EdgeID(int i) const
459 {
460   return myEdgeID[ i % myEdgeID.size() ];
461 }
462
463 //================================================================================
464 /*!
465    * \brief Return p-curve of i-th wrapped edge (count starts from zero)
466    */
467 //================================================================================
468
469 inline Handle(Geom2d_Curve) StdMeshers_FaceSide::Curve2d(int i) const
470 {
471   return myC2d[ i % myC2d.size() ];
472 }
473
474 //================================================================================
475 /*!
476  * \brief Return length of i-th wrapped edge (count starts from zero)
477  */
478  //================================================================================
479
480 inline double StdMeshers_FaceSide::EdgeLength(int i) const
481 {
482   return myEdgeLength[ i % myEdgeLength.size() ];
483 }
484
485 //================================================================================
486 /*!
487  * \brief Return orientation of i-th wrapped edge (count starts from zero)
488  */
489  //================================================================================
490
491 inline bool StdMeshers_FaceSide::IsReversed(int i) const
492 {
493   return myFirst[i] > myLast[i];
494 }
495
496 #endif