]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_FaceSide.hxx
Salome HOME
MeshFormatWriter and MeshFormaReader are in MEDCOUPLING namespace
[modules/smesh.git] / src / StdMeshers / StdMeshers_FaceSide.hxx
1 // Copyright (C) 2007-2020  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                       const 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 "constructors"
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                                      const 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                                      const TopoDS_Edge& theEdge = TopoDS_Edge(),
146                                      SMESH_Mesh*        theMesh = 0)
147   {
148     return StdMeshers_FaceSidePtr( new StdMeshers_FaceSide( theSideNodes, theFace, theEdge, theMesh ));
149   }
150
151   /*!
152    * \brief Return wires of a face as StdMeshers_FaceSide's
153    */
154   static TSideVector GetFaceWires(const TopoDS_Face&   theFace,
155                                   SMESH_Mesh &         theMesh,
156                                   const bool           theIgnoreMediumNodes,
157                                   TError &             theError,
158                                   SMESH_MesherHelper*  theFaceHelper = NULL,
159                                   SMESH_ProxyMesh::Ptr theProxyMesh = SMESH_ProxyMesh::Ptr(),
160                                   const bool           theCheckVertexNodes = true);
161   /*!
162    * \brief Change orientation of side geometry
163    */
164   void Reverse();
165   /*!
166    * \brief Make ignore medium nodes
167    */
168   void SetIgnoreMediumNodes(bool toIgnore);
169
170   /*!
171    * \brief Return nb nodes on edges and vertices (+1 to be == GetUVPtStruct().size() ).
172    *        Call it with update == true if mesh of this side can be recomputed
173    *        since creation of this side
174    */
175   int NbPoints(const bool update = false) const;
176   /*!
177    * \brief Return nb edges
178    *        Call it with update == true if mesh of this side can be recomputed
179    *        since creation of this side
180    */
181   int NbSegments(const bool update = false) const;
182   /*!
183    * \brief Return mesh
184    */
185   SMESH_Mesh* GetMesh() const { return myProxyMesh->GetMesh(); }
186   /*!
187    * \brief Return true if there are vertices without nodes
188    */
189   bool MissVertexNode() const { return myMissingVertexNodes; }
190
191   /*!
192    * \brief Return detailed data on nodes
193     * \param isXConst - true if normalized parameter X is constant
194     * \param constValue - constant parameter value
195     *
196     * Missing nodes are allowed only on internal vertices.
197     * For a closed side, the 1st point repeats at end
198    */
199   const UVPtStructVec& GetUVPtStruct( bool isXConst = 0, double constValue = 0 ) const;
200   /*!
201    * \brief Simulates detailed data on nodes
202     * \param isXConst - true if normalized parameter X is constant
203     * \param constValue - constant parameter value
204    */
205   const UVPtStructVec& SimulateUVPtStruct(int    nbSeg,
206                                           bool   isXConst   = 0,
207                                           double constValue = 0) const;
208   /*!
209    * \brief Return nodes in the order they encounter while walking along
210    *  the whole side or a specified EDGE. For a closed side, the 1st point repeats at end.
211    *  \param iE - index of the EDGE. Default is "all EDGEs".
212    */
213   std::vector<const SMDS_MeshNode*> GetOrderedNodes( int iE = ALL_EDGES ) const;
214
215   /*!
216    * \brief Return nodes of the i-th EDGE.
217    *        Nodes moved to other geometry by MergeNodes() are also returned.
218    * \retval bool - is OK
219    */
220   bool GetEdgeNodes(const size_t                       i,
221                     std::vector<const SMDS_MeshNode*>& nodes,
222                     bool                               inlude1stVertex=true,
223                     bool                               inludeLastVertex=true) const;
224
225   /*!
226    * \brief Return a node from the i-th VERTEX (count starts from zero)
227    *        Nodes moved to other geometry by MergeNodes() are also returned.
228    */
229   const SMDS_MeshNode* VertexNode(std::size_t i, bool* isMoved = 0) const;
230
231   /*
232    * \brief Return edge and parameter on edge by normalized parameter
233    */
234   inline double Parameter(double U, TopoDS_Edge & edge) const;
235   /*
236    * \brief Return edge ID and parameter on edge by normalized parameter
237    */
238   inline double Parameter(double U, int & edgeID) const;
239   /*!
240    * \brief Return UV by normalized parameter
241    */
242   gp_Pnt2d Value2d(double U) const;
243   /*!
244    * \brief Return XYZ by normalized parameter
245    */
246   gp_Pnt   Value3d(double U) const;
247   /*!
248    * \brief Creates a Adaptor2d_Curve2d to be used in SMESH_Block
249    */
250   Adaptor2d_Curve2d* GetCurve2d() const;
251   /*!
252    * \brief Creates a fully functional Adaptor_Curve
253    */
254   BRepAdaptor_CompCurve* GetCurve3d() const;
255   /*!
256    * \brief Return nb of wrapped edges
257    */
258   int NbEdges() const { return myEdge.size(); }
259   /*!
260    * \brief Return i-th edge (count starts from zero)
261    */
262   const TopoDS_Edge& Edge(int i) const { return myEdge[i]; }
263   /*!
264    * \brief Return all edges
265    */
266   const std::vector<TopoDS_Edge>& Edges() const { return myEdge; }
267   /*!
268    * \brief Return the FACE
269    */
270   const TopoDS_Face& Face() const { return myFace; }
271   /*!
272    * \brief Return 1st vertex of the i-th edge (count starts from zero)
273    */
274   TopoDS_Vertex FirstVertex(int i=0) const;
275   /*!
276    * \brief Return last vertex of the i-th edge (count starts from zero)
277    */
278   TopoDS_Vertex LastVertex(int i = LAST_EDGE) const;
279   /*!
280    * \brief Return \c true if the chain of EDGEs is closed
281    */
282   bool IsClosed() const;
283   /*!
284    * \brief Return side length
285    */
286   double Length() const { return myLength; }
287   /*!
288    * \brief Return edge index corresponding to normalized parameter
289    */
290   inline int EdgeIndex( double U ) const;
291
292   void dump(const char* msg=0) const;
293   
294   /*!
295    * \brief Return ID of i-th wrapped edge (count starts from zero)
296    */
297   inline int EdgeID(int i) const;
298   /*!
299    * \brief Return p-curve of i-th wrapped edge (count starts from zero)
300    */
301   inline Handle(Geom2d_Curve) Curve2d(int i) const;
302   /*!
303    * \brief Return first normalized parameter of the i-th edge (count starts from zero)
304    */
305   inline double FirstParameter(int i) const;
306   /*!
307    * \brief Return last normalized parameter of the i-th edge (count starts from zero)
308    */
309   inline double LastParameter(int i) const;
310   /*!
311    * \brief Return first parameter of the i-th edge (count starts from zero).
312    *        EDGE orientation is taken into account
313    */
314   inline double FirstU(int i) const;
315   /*!
316    * \brief Return last parameter of the i-th edge (count starts from zero).
317    *        EDGE orientation is taken into account
318    */
319   inline double LastU(int i) const;
320   /*!
321    * \brief Return length of i-th wrapped edge (count starts from zero)
322    */
323   inline double EdgeLength(int i) const;
324   /*!
325    * \brief Return orientation of i-th wrapped edge (count starts from zero)
326    */
327   inline bool IsReversed(int i) const;
328   /*!
329    * \brief Return a helper initialized with the FACE
330    */
331   SMESH_MesherHelper* FaceHelper() const;
332
333 protected:
334
335   void reverseProxySubmesh( const TopoDS_Edge& E );
336
337   // DON't FORGET to update Reverse() when adding one more vector!
338   TopoDS_Face                       myFace;
339   std::vector<uvPtStruct>           myPoints, myFalsePoints;
340   std::vector<TopoDS_Edge>          myEdge;
341   std::vector<int>                  myEdgeID;
342   std::vector<Handle(Geom2d_Curve)> myC2d;
343   std::vector<GeomAdaptor_Curve>    myC3dAdaptor;
344   std::vector<double>               myFirst, myLast;
345   std::vector<double>               myNormPar;
346   std::vector<double>               myEdgeLength;
347   std::vector<int>                  myIsUniform;
348   double                            myLength;
349   int                               myNbPonits, myNbSegments;
350   SMESH_ProxyMesh::Ptr              myProxyMesh;
351   bool                              myMissingVertexNodes, myIgnoreMediumNodes;
352   gp_Pnt2d                          myDefaultPnt2d;
353   SMESH_MesherHelper*               myHelper;
354 };
355
356
357 //================================================================================
358 /*!
359  * \brief Return edge index corresponding to normalized parameter
360   * \param U - the parameter
361   * \retval int - index
362  */
363 //================================================================================
364
365 inline int StdMeshers_FaceSide::EdgeIndex( double U ) const
366 {
367   int i = myNormPar.size() - 1;
368   while ( i > 0 && U < myNormPar[ i-1 ] ) --i;
369   return i;
370 }
371
372 //================================================================================
373 /*!
374  * \brief Return an edge and parameter on the edge by a normalized parameter
375   * \param U - normalized parameter
376   * \retval double - pameter on a curve
377   * \ warning The returned parameter can be inaccurate if the edge is non-uniformly
378   *           parametrized. Use Value2d() to get a precise point on the edge
379  */
380 //================================================================================
381
382 inline double StdMeshers_FaceSide::Parameter(double U, TopoDS_Edge & edge) const
383 {
384   int i = EdgeIndex( U );
385   edge = myEdge[ i ];
386   double prevU = i ? myNormPar[ i-1 ] : 0;
387   double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
388   return myFirst[i] * ( 1 - r ) + myLast[i] * r;
389 }
390
391 //================================================================================
392 /*!
393  * \brief Return an edge ID and parameter on the edge by a normalized parameter
394   * \param U - normalized parameter
395   * \retval double - pameter on a curve
396   * \ warning The returned parameter can be inaccurate if the edge is non-uniformly
397   *           parametrized. Use Value2d() to get a precise point on the edge
398  */
399 //================================================================================
400
401 inline double StdMeshers_FaceSide::Parameter(double U, int & edgeID) const
402 {
403   int i = EdgeIndex( U );
404   edgeID = myEdgeID[ i ];
405   double prevU = i ? myNormPar[ i-1 ] : 0;
406   double r = ( U - prevU )/ ( myNormPar[ i ] - prevU );
407   return myFirst[i] * ( 1 - r ) + myLast[i] * r;
408 }
409
410 //================================================================================
411 /*!
412  * \brief Return first normalized parameter of the i-th edge
413  */
414 //================================================================================
415
416 inline double StdMeshers_FaceSide::FirstParameter(int i) const
417 {
418   return i==0 ? 0. : i<(int)myNormPar.size() ? myNormPar[i-1] : 1.;
419 }
420
421 //================================================================================
422 /*!
423  * \brief Return ast normalized parameter of the i-th edge
424  */
425 //================================================================================
426
427 inline double StdMeshers_FaceSide::LastParameter(int i) const
428 {
429   return i < (int)myNormPar.size() ? myNormPar[i] : 1;
430 }
431
432 //================================================================================
433 /*!
434  * \brief Return first parameter of the i-th edge
435  */
436 //================================================================================
437
438 inline double StdMeshers_FaceSide::FirstU(int i) const
439 {
440   return myFirst[ i % myFirst.size() ];
441 }
442
443 //================================================================================
444 /*!
445  * \brief Return last parameter of the i-th edge
446  */
447 //================================================================================
448
449 inline double StdMeshers_FaceSide::LastU(int i) const
450 {
451   return myLast[ i % myLast.size() ];
452 }
453
454 //================================================================================
455   /*!
456    * \brief Return ID of i-th wrapped edge (count starts from zero)
457    */
458 //================================================================================
459
460 inline int StdMeshers_FaceSide::EdgeID(int i) const
461 {
462   return myEdgeID[ i % myEdgeID.size() ];
463 }
464
465 //================================================================================
466 /*!
467    * \brief Return p-curve of i-th wrapped edge (count starts from zero)
468    */
469 //================================================================================
470
471 inline Handle(Geom2d_Curve) StdMeshers_FaceSide::Curve2d(int i) const
472 {
473   return myC2d[ i % myC2d.size() ];
474 }
475
476 //================================================================================
477 /*!
478  * \brief Return length of i-th wrapped edge (count starts from zero)
479  */
480  //================================================================================
481
482 inline double StdMeshers_FaceSide::EdgeLength(int i) const
483 {
484   return myEdgeLength[ i % myEdgeLength.size() ];
485 }
486
487 //================================================================================
488 /*!
489  * \brief Return orientation of i-th wrapped edge (count starts from zero)
490  */
491  //================================================================================
492
493 inline bool StdMeshers_FaceSide::IsReversed(int i) const
494 {
495   return myFirst[i] > myLast[i];
496 }
497
498 #endif