Version: 9.12.0
SMESH_MesherHelper.hxx
Go to the documentation of this file.
1 // Copyright (C) 2007-2023 CEA, EDF, 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: SMESH_MesherHelper.hxx
24 // Created: 15.02.06 14:48:09
25 // Author: Sergey KUUL
26 //
27 #ifndef SMESH_MesherHelper_HeaderFile
28 #define SMESH_MesherHelper_HeaderFile
29 
30 #include "SMESH_SMESH.hxx"
31 
32 #include "SMESH_ComputeError.hxx"
33 #include "SMESH_TypeDefs.hxx"
34 
35 #include <Geom_Surface.hxx>
36 #include <ShapeAnalysis_Surface.hxx>
37 #include <TopoDS_Face.hxx>
38 #include <TopoDS_Shape.hxx>
39 #include <gp_Pnt2d.hxx>
40 
41 #include <map>
42 #include <vector>
43 
44 class GeomAPI_ProjectPointOnCurve;
45 class GeomAPI_ProjectPointOnSurf;
46 class SMDS_MeshNode;
47 class SMESHDS_Hypothesis;
48 class SMESHDS_Mesh;
49 class SMESHDS_SubMesh;
50 class SMESH_Gen;
51 class SMESH_Mesh;
52 class SMESH_ProxyMesh;
53 class SMESH_subMesh;
54 class TopoDS_Edge;
55 class TopoDS_Face;
56 class TopoDS_Vertex;
57 
58 typedef std::map<SMESH_TLink, const SMDS_MeshNode*> TLinkNodeMap;
59 typedef std::map<SMESH_TLink, const SMDS_MeshNode*>::iterator ItTLinkNode;
60 
62 typedef boost::shared_ptr< PShapeIterator > PShapeIteratorPtr;
63 
64 typedef std::vector<const SMDS_MeshNode* > TNodeColumn;
65 typedef std::map< double, TNodeColumn > TParam2ColumnMap;
66 
67 typedef gp_XY (*xyFunPtr)(const gp_XY& uv1, const gp_XY& uv2);
68 
69 //=======================================================================
81 //=======================================================================
82 
84 {
85  public:
86  // ---------- PUBLIC UTILITIES ----------
87 
95  static bool IsSameElemGeometry(const SMESHDS_SubMesh* smDS,
97  const bool nullSubMeshRes = true);
98 
112  static bool LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
113  const TopoDS_Face& theFace,
114  const std::list<TopoDS_Edge>& theBaseSide,
115  SMESHDS_Mesh* theMesh,
116  SMESH_ProxyMesh* theProxyMesh=0);
120  static bool LoadNodeColumns(TParam2ColumnMap & theParam2ColumnMap,
121  const TopoDS_Face& theFace,
122  const TopoDS_Edge& theBaseEdge,
123  SMESHDS_Mesh* theMesh,
124  SMESH_ProxyMesh* theProxyMesh=0);
128  static bool IsStructured( SMESH_subMesh* faceSM );
129 
133  static bool IsCornerOfStructure( const SMDS_MeshNode* n,
134  const SMESHDS_SubMesh* faceSM,
135  SMESH_MesherHelper& faceAnalyser );
136 
140  static bool IsDistorted2D( SMESH_subMesh* faceSM,
141  bool checkUV = false,
142  SMESH_MesherHelper* faceHelper = NULL);
143 
150  static bool IsMedium(const SMDS_MeshNode* node,
151  const SMDSAbs_ElementType typeToCheck = SMDSAbs_All);
159  static TopoDS_Shape GetSubShapeByNode(const SMDS_MeshNode* node,
160  const SMESHDS_Mesh* meshDS);
161 
168  static inline int WrapIndex(int ind, const int nbNodes) {
169  return (( ind %= nbNodes ) < 0 ) ? ind + nbNodes : ind;
170  }
171 
193  inline static gp_XY calcTFI(double x, double y,
194  const gp_XY& a0,const gp_XY& a1,const gp_XY& a2,const gp_XY& a3,
195  const gp_XY& p0,const gp_XY& p1,const gp_XY& p2,const gp_XY& p3);
196 
200  inline static gp_XYZ calcTFI(double x, double y,
201  const gp_XYZ& a0,const gp_XYZ& a1,const gp_XYZ& a2,const gp_XYZ& a3,
202  const gp_XYZ& p0,const gp_XYZ& p1,const gp_XYZ& p2,const gp_XYZ& p3);
210  static int Count(const TopoDS_Shape& shape,
211  const TopAbs_ShapeEnum type,
212  const bool ignoreSame);
213 
217  static int NbAncestors(const TopoDS_Shape& shape,
218  const SMESH_Mesh& mesh,
219  TopAbs_ShapeEnum ancestorType=TopAbs_SHAPE);
223  static PShapeIteratorPtr GetAncestors(const TopoDS_Shape& shape,
224  const SMESH_Mesh& mesh,
225  TopAbs_ShapeEnum ancestorType,
226  const TopoDS_Shape* container = 0);
230  static TopoDS_Shape GetCommonAncestor(const TopoDS_Shape& shape1,
231  const TopoDS_Shape& shape2,
232  const SMESH_Mesh& mesh,
233  TopAbs_ShapeEnum ancestorType);
237  static TopAbs_Orientation GetSubShapeOri(const TopoDS_Shape& shape,
238  const TopoDS_Shape& subShape);
239 
240  static bool IsSubShape( const TopoDS_Shape& shape, const TopoDS_Shape& mainShape );
241 
242  static bool IsSubShape( const TopoDS_Shape& shape, SMESH_Mesh* aMesh );
243 
244  static bool IsBlock( const TopoDS_Shape& shape );
245 
246  static double MaxTolerance( const TopoDS_Shape& shape );
247 
248  static double GetAngle( const TopoDS_Edge & E1, const TopoDS_Edge & E2,
249  const TopoDS_Face & F, const TopoDS_Vertex & V,
250  gp_Vec* faceNormal=0);
251 
252  static bool IsClosedEdge( const TopoDS_Edge& anEdge );
253 
254  static TopoDS_Vertex IthVertex( const bool is2nd, TopoDS_Edge anEdge, const bool CumOri=true );
255 
256  static TopAbs_ShapeEnum GetGroupType(const TopoDS_Shape& group,
257  const bool avoidCompound=false);
258 
259  static TopoDS_Shape GetShapeOfHypothesis( const SMESHDS_Hypothesis * hyp,
260  const TopoDS_Shape& shape,
261  SMESH_Mesh* mesh);
262 
263 
264 public:
265  // ---------- PUBLIC INSTANCE METHODS ----------
266 
267  // constructor
268  SMESH_MesherHelper(SMESH_Mesh& theMesh);
269 
270  SMESH_Gen* GetGen() const;
271 
272  SMESH_Mesh* GetMesh() const { return myMesh; }
273 
274  SMESHDS_Mesh* GetMeshDS() const;
275 
280  bool IsQuadraticSubMesh(const TopoDS_Shape& theShape);
281 
285  void SetIsQuadratic(const bool theBuildQuadratic)
286  { myCreateQuadratic = theBuildQuadratic; }
287 
291  void SetIsBiQuadratic(const bool theBuildBiQuadratic)
292  { myCreateBiQuadratic = theBuildBiQuadratic; }
293 
297  bool GetIsQuadratic() const { return myCreateQuadratic; }
298 
299  /*
300  * \brief Find out elements orientation on a geometrical face
301  */
302  bool IsReversedSubMesh (const TopoDS_Face& theFace);
303 
307  bool GetIsBiQuadratic() const { return myCreateBiQuadratic; }
308 
314  void FixQuadraticElements(SMESH_ComputeErrorPtr& error, bool volumeOnly=true);
315 
321  bool SetElementsOnShape(bool toSet)
322  { bool res = mySetElemOnShape; mySetElemOnShape = toSet; return res; }
323 
327  void SetSubShape(const int subShapeID);
328  void SetSubShape(const TopoDS_Shape& subShape);
333  int GetSubShapeID() const { return myShapeID; }
337  const TopoDS_Shape& GetSubShape() const { return myShape; }
342  void CopySubShapeInfo(const SMESH_MesherHelper& other);
343 
344 
348  int ShapeToIndex( const TopoDS_Shape& S ) const;
349 
353  SMDS_MeshNode* AddNode(double x, double y, double z, smIdType ID = 0, double u=0., double v=0.);
357  SMDS_MeshEdge* AddEdge(const SMDS_MeshNode* n1,
358  const SMDS_MeshNode* n2,
359  const smIdType id = 0,
360  const bool force3d = true);
364  SMDS_MeshFace* AddFace(const SMDS_MeshNode* n1,
365  const SMDS_MeshNode* n2,
366  const SMDS_MeshNode* n3,
367  const smIdType id=0,
368  const bool force3d = false);
372  SMDS_MeshFace* AddFace(const SMDS_MeshNode* n1,
373  const SMDS_MeshNode* n2,
374  const SMDS_MeshNode* n3,
375  const SMDS_MeshNode* n4,
376  const smIdType id = 0,
377  const bool force3d = false);
381  SMDS_MeshFace* AddPolygonalFace (const std::vector<const SMDS_MeshNode*>& nodes,
382  const smIdType id = 0,
383  const bool force3d = false);
387  SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
388  const SMDS_MeshNode* n2,
389  const SMDS_MeshNode* n3,
390  const SMDS_MeshNode* n4,
391  const smIdType id = 0,
392  const bool force3d = true);
396  SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
397  const SMDS_MeshNode* n2,
398  const SMDS_MeshNode* n3,
399  const SMDS_MeshNode* n4,
400  const SMDS_MeshNode* n5,
401  const smIdType id = 0,
402  const bool force3d = true);
406  SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
407  const SMDS_MeshNode* n2,
408  const SMDS_MeshNode* n3,
409  const SMDS_MeshNode* n4,
410  const SMDS_MeshNode* n5,
411  const SMDS_MeshNode* n6,
412  const smIdType id = 0,
413  const bool force3d = true);
417  SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
418  const SMDS_MeshNode* n2,
419  const SMDS_MeshNode* n3,
420  const SMDS_MeshNode* n4,
421  const SMDS_MeshNode* n5,
422  const SMDS_MeshNode* n6,
423  const SMDS_MeshNode* n7,
424  const SMDS_MeshNode* n8,
425  const smIdType id = 0,
426  bool force3d = true);
427 
431  SMDS_MeshVolume* AddVolume(const SMDS_MeshNode* n1,
432  const SMDS_MeshNode* n2,
433  const SMDS_MeshNode* n3,
434  const SMDS_MeshNode* n4,
435  const SMDS_MeshNode* n5,
436  const SMDS_MeshNode* n6,
437  const SMDS_MeshNode* n7,
438  const SMDS_MeshNode* n8,
439  const SMDS_MeshNode* n9,
440  const SMDS_MeshNode* n10,
441  const SMDS_MeshNode* n11,
442  const SMDS_MeshNode* n12,
443  const smIdType id = 0,
444  bool force3d = true);
445 
449  SMDS_MeshVolume* AddPolyhedralVolume (const std::vector<const SMDS_MeshNode*>& nodes,
450  const std::vector<int>& quantities,
451  const smIdType ID=0,
452  const bool force3d = true);
459  void ToFixNodeParameters(bool toFix);
460 
464  double GetNodeU(const TopoDS_Edge& theEdge,
465  const SMDS_MeshNode* theNode,
466  const SMDS_MeshNode* inEdgeNode=0,
467  bool* check=0) const;
473  gp_XY GetNodeUV(const TopoDS_Face& F,
474  const SMDS_MeshNode* n,
475  const SMDS_MeshNode* inFaceNode=0,
476  bool* check=0) const;
483  bool CheckNodeUV(const TopoDS_Face& F,
484  const SMDS_MeshNode* n,
485  gp_XY& uv,
486  const double tol,
487  const bool force=false,
488  double distXYZ[4]=0) const;
495  bool CheckNodeU(const TopoDS_Edge& E,
496  const SMDS_MeshNode* n,
497  double& u,
498  const double tol,
499  const bool force=false,
500  double distXYZ[4]=0) const;
504  static gp_XY GetMiddleUV(const Handle(Geom_Surface)& surface,
505  const gp_XY& uv1,
506  const gp_XY& uv2);
510  static gp_XY GetCenterUV(const gp_XY& uv1,
511  const gp_XY& uv2,
512  const gp_XY& uv3,
513  const gp_XY& uv12,
514  const gp_XY& uv23,
515  const gp_XY& uv31,
516  bool * isBadTria=0);
524 #define gp_XY_FunPtr(meth) \
525  static gp_XY __gpXY_##meth (const gp_XY& uv1, const gp_XY& uv2) { return uv1.meth( uv2 ); } \
526  static xyFunPtr gp_XY_##meth = & __gpXY_##meth
527 
533  static gp_XY ApplyIn2D(Handle(Geom_Surface) surface,
534  const gp_XY& uv1,
535  const gp_XY& uv2,
536  xyFunPtr fun,
537  const bool resultInPeriod=true);
538 
545  void AdjustByPeriod( const TopoDS_Face& face, gp_XY uv[], const int nbUV );
546 
554  bool GetNodeUVneedInFaceNode(const TopoDS_Face& F = TopoDS_Face()) const;
555 
559  GeomAPI_ProjectPointOnSurf& GetProjector(const TopoDS_Face& F,
560  TopLoc_Location& loc,
561  double tol=0 ) const;
565  GeomAPI_ProjectPointOnSurf& GetProjector(const TopoDS_Face& F,
566  double tol=0 ) const;
570  GeomAPI_ProjectPointOnCurve& GetPCProjector(const TopoDS_Edge& E ) const;
574  Handle(ShapeAnalysis_Surface) GetSurface(const TopoDS_Face& F ) const;
575 
583  bool IsDegenShape(const int subShape) const
584  { return myDegenShapeIds.find( subShape ) != myDegenShapeIds.end(); }
590  bool HasDegeneratedEdges() const { return !myDegenShapeIds.empty(); }
596  size_t NbDegeneratedEdges() const { return myDegenShapeIds.size(); }
597 
606  bool IsSeamShape(const int subShape) const
607  { return mySeamShapeIds.find( subShape ) != mySeamShapeIds.end(); }
616  bool IsSeamShape(const TopoDS_Shape& subShape) const
617  { return IsSeamShape( ShapeToIndex( subShape )); }
622  bool IsRealSeam(const int subShape) const
623  { return mySeamShapeIds.find( -subShape ) != mySeamShapeIds.end(); }
628  bool IsRealSeam(const TopoDS_Shape& subShape) const
629  { return IsRealSeam( ShapeToIndex( subShape )); }
636  bool HasSeam() const { return !mySeamShapeIds.empty(); }
642  bool HasRealSeam() const { return HasSeam() && ( *mySeamShapeIds.begin() < 0 ); }
648  size_t NbRealSeam() const;
653  int GetPeriodicIndex() const { return myParIndex; }
657  double GetPeriod(int perioIndex) const { return myPar2[ perioIndex-1 ] - myPar1[ perioIndex-1 ]; }
661  double GetOtherParam(const double param) const;
665  int IsOnSeam(const gp_XY& uv) const;
666 
676  const SMDS_MeshNode* GetMediumNode(const SMDS_MeshNode* n1,
677  const SMDS_MeshNode* n2,
678  const bool force3d,
679  TopAbs_ShapeEnum expectedSupport=TopAbs_SHAPE);
686  const SMDS_MeshNode* GetCentralNode(const SMDS_MeshNode* n1,
687  const SMDS_MeshNode* n2,
688  const SMDS_MeshNode* n3,
689  const SMDS_MeshNode* n4,
690  const SMDS_MeshNode* n12,
691  const SMDS_MeshNode* n23,
692  const SMDS_MeshNode* n34,
693  const SMDS_MeshNode* n41,
694  bool force3d);
701  const SMDS_MeshNode* GetCentralNode(const SMDS_MeshNode* n1,
702  const SMDS_MeshNode* n2,
703  const SMDS_MeshNode* n3,
704  const SMDS_MeshNode* n12,
705  const SMDS_MeshNode* n23,
706  const SMDS_MeshNode* n31,
707  bool force3d);
711  std::pair<int, TopAbs_ShapeEnum> GetMediumPos(const SMDS_MeshNode* n1,
712  const SMDS_MeshNode* n2,
713  const bool useCurSubShape=false,
714  TopAbs_ShapeEnum expectedSupport=TopAbs_SHAPE);
718  void AddTLinkNode(const SMDS_MeshNode* n1,
719  const SMDS_MeshNode* n2,
720  const SMDS_MeshNode* n12);
724  void AddTLinkNodeMap(const TLinkNodeMap& aMap)
725  { myTLinkNodeMap.insert(aMap.begin(), aMap.end()); }
726 
727  bool AddTLinks(const SMDS_MeshEdge* edge);
728  bool AddTLinks(const SMDS_MeshFace* face);
729  bool AddTLinks(const SMDS_MeshVolume* vol);
730 
734  const TLinkNodeMap& GetTLinkNodeMap() const { return myTLinkNodeMap; }
735 
741  enum MType{ LINEAR, QUADRATIC, COMP };
742  MType IsQuadraticMesh();
743 
744  virtual ~SMESH_MesherHelper();
745 
746  static void WriteShape(const TopoDS_Shape& s);
747 
748 
749  protected:
750 
757  gp_Pnt2d getUVOnSeam( const gp_Pnt2d& uv1, const gp_Pnt2d& uv2 ) const;
758 
759  const SMDS_MeshNode* getMediumNodeOnComposedWire(const SMDS_MeshNode* n1,
760  const SMDS_MeshNode* n2,
761  bool force3d);
762 
763  double getFaceMaxTol( const TopoDS_Shape& face ) const;
764 
765 
766  private:
767 
768  // forbidden copy constructor
770 
771  // key of a map of bi-quadratic face to it's central node
772  struct TBiQuad: public std::pair<smIdType, std::pair<smIdType, smIdType> >
773  {
775  const SMDS_MeshNode* n2,
776  const SMDS_MeshNode* n3,
777  const SMDS_MeshNode* n4=0)
778  {
780  s.insert(n1);
781  s.insert(n2);
782  s.insert(n3);
783  if ( n4 ) s.insert(n4);
784  TIDSortedNodeSet::iterator n = s.begin();
785  first = (*n++)->GetID();
786  second.first = (*n++)->GetID();
787  second.second = (*n++)->GetID();
788  }
789  };
790 
791  // maps used during creation of quadratic elements
792  TLinkNodeMap myTLinkNodeMap; // medium nodes on links
793  std::map< TBiQuad, const SMDS_MeshNode* > myMapWithCentralNode; // central nodes of faces
794 
795  std::set< int > myDegenShapeIds;
796  std::set< int > mySeamShapeIds;
797  double myPar1[2], myPar2[2]; // U and V bounds of a closed periodic surface
798  int myParIndex; // bounds' index (1-U, 2-V, 3-both)
799 
800  std::map< int, double > myFaceMaxTol;
801 
802  typedef std::map< int, Handle(ShapeAnalysis_Surface)> TID2Surface;
803  typedef std::map< int, GeomAPI_ProjectPointOnSurf* > TID2ProjectorOnSurf;
804  typedef std::map< int, GeomAPI_ProjectPointOnCurve* > TID2ProjectorOnCurve;
808 
809  TopoDS_Shape myShape;
812 
817 
818  std::map< int,bool > myNodePosShapesValidity;
819  bool toCheckPosOnShape(int shapeID ) const;
820  void setPosOnShapeValidity(int shapeID, bool ok ) const;
821 };
822 
823 //=======================================================================
824 inline gp_XY
825 SMESH_MesherHelper::calcTFI(double x, double y,
826  const gp_XY& a0,const gp_XY& a1,const gp_XY& a2,const gp_XY& a3,
827  const gp_XY& p0,const gp_XY& p1,const gp_XY& p2,const gp_XY& p3)
828 {
829  return
830  ((1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3 ) -
831  ((1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3);
832 }
833 //=======================================================================
834 inline gp_XYZ
835 SMESH_MesherHelper::calcTFI(double x, double y,
836  const gp_XYZ& a0,const gp_XYZ& a1,const gp_XYZ& a2,const gp_XYZ& a3,
837  const gp_XYZ& p0,const gp_XYZ& p1,const gp_XYZ& p2,const gp_XYZ& p3)
838 {
839  return
840  ((1 - y) * p0 + x * p1 + y * p2 + (1 - x) * p3 ) -
841  ((1 - x) * (1 - y) * a0 + x * (1 - y) * a1 + x * y * a2 + (1 - x) * y * a3);
842 }
843 //=======================================================================
844 
845 #endif
Handle(SALOME_InteractiveObject) GeomSelectionTools
Return the first selected Salome Interactive Object (Handle(Salome_InteractiveObject))
Definition: GeomSelectionTools.cxx:97
SMDSAbs_GeometryType
enumeration for element geometry type
Definition: SMDSAbs_ElementType.hxx:47
SMDSAbs_ElementType
Type (node, edge, face or volume) of elements.
Definition: SMDSAbs_ElementType.hxx:34
@ SMDSAbs_All
Definition: SMDSAbs_ElementType.hxx:35
boost::shared_ptr< SMESH_ComputeError > SMESH_ComputeErrorPtr
Definition: SMESH_ComputeError.hxx:36
static bool IsSubShape(const TopTools_IndexedMapOfShape &theMap, const TopoDS_Shape &theShape)
Definition: SMESH_Controls.cxx:5141
std::map< double, TNodeColumn > TParam2ColumnMap
Definition: SMESH_MesherHelper.hxx:65
SMDS_Iterator< const TopoDS_Shape * > PShapeIterator
Definition: SMESH_MesherHelper.hxx:61
boost::shared_ptr< PShapeIterator > PShapeIteratorPtr
Definition: SMESH_MesherHelper.hxx:62
std::map< SMESH_TLink, const SMDS_MeshNode * > TLinkNodeMap
Definition: SMESH_MesherHelper.hxx:56
gp_XY(* xyFunPtr)(const gp_XY &uv1, const gp_XY &uv2)
Definition: SMESH_MesherHelper.hxx:67
std::vector< const SMDS_MeshNode * > TNodeColumn
Definition: SMESH_MesherHelper.hxx:64
std::map< SMESH_TLink, const SMDS_MeshNode * >::iterator ItTLinkNode
Definition: SMESH_MesherHelper.hxx:59
std::set< const SMDS_MeshNode *, TIDCompare > TIDSortedNodeSet
Definition: SMESH_OctreeNode.hxx:51
#define SMESH_EXPORT
Definition: SMESH_SMESH.hxx:37
Abstract class for iterators.
Definition: SMDS_Iterator.hxx:33
Edge mesh element.
Definition: SMDS_MeshEdge.hxx:36
Mesh face.
Definition: SMDS_MeshFace.hxx:41
Definition: SMDS_MeshNode.hxx:36
Mesh volume.
Definition: SMDS_MeshVolume.hxx:40
Definition: SMESHDS_Hypothesis.hxx:37
Definition: SMESHDS_Mesh.hxx:68
Definition: SMESHDS_SubMesh.hxx:48
Definition: SMESH_Gen.hxx:68
Definition: SMESH_Mesh.hxx:80
It helps meshers to add elements and provides other utilities.
Definition: SMESH_MesherHelper.hxx:84
bool IsRealSeam(const int subShape) const
Return true if an edge or a vertex encounters twice in face wire.
Definition: SMESH_MesherHelper.hxx:622
static gp_XY ApplyIn2D(Handle(Geom_Surface) surface, const gp_XY &uv1, const gp_XY &uv2, xyFunPtr fun, const bool resultInPeriod=true)
Perform given operation on two 2d points in parameric space of given surface.
TopoDS_Shape myShape
Definition: SMESH_MesherHelper.hxx:809
bool IsSeamShape(const TopoDS_Shape &subShape) const
Check if shape is a seam edge or it's vertex.
Definition: SMESH_MesherHelper.hxx:616
std::set< int > mySeamShapeIds
Definition: SMESH_MesherHelper.hxx:796
std::set< int > myDegenShapeIds
Definition: SMESH_MesherHelper.hxx:795
bool GetIsQuadratic() const
Return myCreateQuadratic flag.
Definition: SMESH_MesherHelper.hxx:297
bool HasDegeneratedEdges() const
Check if the shape set through IsQuadraticSubMesh() or SetSubShape() has a degenerated edges.
Definition: SMESH_MesherHelper.hxx:590
int myShapeID
Definition: SMESH_MesherHelper.hxx:811
std::map< int, Handle(ShapeAnalysis_Surface)> TID2Surface
Definition: SMESH_MesherHelper.hxx:802
std::map< int, bool > myNodePosShapesValidity
Definition: SMESH_MesherHelper.hxx:818
std::map< TBiQuad, const SMDS_MeshNode * > myMapWithCentralNode
Definition: SMESH_MesherHelper.hxx:793
bool GetIsBiQuadratic() const
Return myCreateBiQuadratic flag.
Definition: SMESH_MesherHelper.hxx:307
bool SetElementsOnShape(bool toSet)
To set created elements on the shape set by IsQuadraticSubMesh() or the next methods.
Definition: SMESH_MesherHelper.hxx:321
double GetPeriod(int perioIndex) const
Return period in given direction [1,2].
Definition: SMESH_MesherHelper.hxx:657
int myParIndex
Definition: SMESH_MesherHelper.hxx:798
static gp_XY calcTFI(double x, double y, const gp_XY &a0, const gp_XY &a1, const gp_XY &a2, const gp_XY &a3, const gp_XY &p0, const gp_XY &p1, const gp_XY &p2, const gp_XY &p3)
Return UV of a point inside a quadrilateral FACE by it's normalized parameters within a unit quadrang...
Definition: SMESH_MesherHelper.hxx:825
std::map< int, double > myFaceMaxTol
Definition: SMESH_MesherHelper.hxx:800
SMESH_Mesh * myMesh
Definition: SMESH_MesherHelper.hxx:810
TID2ProjectorOnSurf myFace2Projector
Definition: SMESH_MesherHelper.hxx:806
SMESH_MesherHelper(const SMESH_MesherHelper &theOther)
bool mySetElemOnShape
Definition: SMESH_MesherHelper.hxx:815
SMESH_Mesh * GetMesh() const
Definition: SMESH_MesherHelper.hxx:272
bool IsRealSeam(const TopoDS_Shape &subShape) const
Return true if an edge or a vertex encounters twice in face wire.
Definition: SMESH_MesherHelper.hxx:628
TID2ProjectorOnCurve myEdge2Projector
Definition: SMESH_MesherHelper.hxx:807
MType
Check mesh without geometry for: if all elements on this shape are quadratic, quadratic elements will...
Definition: SMESH_MesherHelper.hxx:741
int GetSubShapeID() const
Return ID of the shape set by IsQuadraticSubMesh() or SetSubShape()
Definition: SMESH_MesherHelper.hxx:333
bool HasSeam() const
Check if the shape set through IsQuadraticSubMesh() or SetSubShape() has a seam edge,...
Definition: SMESH_MesherHelper.hxx:636
bool IsSeamShape(const int subShape) const
Check if shape is a seam edge or it's vertex.
Definition: SMESH_MesherHelper.hxx:606
bool myCreateQuadratic
Definition: SMESH_MesherHelper.hxx:813
Handle(ShapeAnalysis_Surface) GetSurface(const TopoDS_Face &F) const
Return a cached ShapeAnalysis_Surface of a FACE.
bool myFixNodeParameters
Definition: SMESH_MesherHelper.hxx:816
TLinkNodeMap myTLinkNodeMap
Definition: SMESH_MesherHelper.hxx:792
std::map< int, GeomAPI_ProjectPointOnSurf * > TID2ProjectorOnSurf
Definition: SMESH_MesherHelper.hxx:803
static int WrapIndex(int ind, const int nbNodes)
Return a valid node index, fixing the given one if necessary.
Definition: SMESH_MesherHelper.hxx:168
void SetIsQuadratic(const bool theBuildQuadratic)
Set order of elements to create without calling IsQuadraticSubMesh()
Definition: SMESH_MesherHelper.hxx:285
bool HasRealSeam() const
Check if the shape set through IsQuadraticSubMesh() or SetSubShape() has a seam edge that encounters ...
Definition: SMESH_MesherHelper.hxx:642
int GetPeriodicIndex() const
Return index of periodic parametric direction of a closed face.
Definition: SMESH_MesherHelper.hxx:653
size_t NbDegeneratedEdges() const
Return a number of degenerated edges in the shape set through IsQuadraticSubMesh() or SetSubShape()
Definition: SMESH_MesherHelper.hxx:596
TID2Surface myFace2Surface
Definition: SMESH_MesherHelper.hxx:805
void SetIsBiQuadratic(const bool theBuildBiQuadratic)
Set myCreateBiQuadratic flag.
Definition: SMESH_MesherHelper.hxx:291
std::map< int, GeomAPI_ProjectPointOnCurve * > TID2ProjectorOnCurve
Definition: SMESH_MesherHelper.hxx:804
bool myCreateBiQuadratic
Definition: SMESH_MesherHelper.hxx:814
const TopoDS_Shape & GetSubShape() const
Return the shape set by IsQuadraticSubMesh() or SetSubShape()
Definition: SMESH_MesherHelper.hxx:337
const TLinkNodeMap & GetTLinkNodeMap() const
Returns myTLinkNodeMap.
Definition: SMESH_MesherHelper.hxx:734
bool IsDegenShape(const int subShape) const
Check if shape is a degenerated edge or it's vertex.
Definition: SMESH_MesherHelper.hxx:583
void AddTLinkNodeMap(const TLinkNodeMap &aMap)
Add many links in my data structure.
Definition: SMESH_MesherHelper.hxx:724
Container of xD mesh elements substituting other ones in the input mesh of an (x+1)D algorithm.
Definition: SMESH_ProxyMesh.hxx:51
Definition: SMESH_subMesh.hxx:61
long AddNode(SMESH::SMESH_Mesh_ptr theMesh, float x, float y, float z, const QStringList &theParameters)
Definition: SMESHGUI_NodesDlg.cxx:96
Definition: SMESH_MesherHelper.hxx:773
TBiQuad(const SMDS_MeshNode *n1, const SMDS_MeshNode *n2, const SMDS_MeshNode *n3, const SMDS_MeshNode *n4=0)
Definition: SMESH_MesherHelper.hxx:774