SALOME - SMESH
SMESH_ControlsDef.hxx
Go to the documentation of this file.
1 // Copyright (C) 2007-2008 CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 #ifndef _SMESH_CONTROLSDEF_HXX_
23 #define _SMESH_CONTROLSDEF_HXX_
24 
25 #include <set>
26 #include <map>
27 #include <vector>
28 
29 #include <boost/shared_ptr.hpp>
30 
31 #include <gp_XYZ.hxx>
32 #include <GeomAPI_ProjectPointOnSurf.hxx>
33 #include <GeomAPI_ProjectPointOnCurve.hxx>
34 #include <TColStd_SequenceOfInteger.hxx>
35 #include <TColStd_MapOfInteger.hxx>
36 #include <TCollection_AsciiString.hxx>
37 #include <TopAbs.hxx>
38 #include <TopoDS_Face.hxx>
39 #include <TopTools_MapOfShape.hxx>
40 #include <BRepClass3d_SolidClassifier.hxx>
41 #include <Quantity_Color.hxx>
42 
43 #include "SMDSAbs_ElementType.hxx"
44 #include "SMDS_MeshNode.hxx"
45 
46 #include "SMESH_Controls.hxx"
47 
48 #ifdef WNT
49  #if defined SMESHCONTROLS_EXPORTS || defined SMESHControls_EXPORTS
50  #define SMESHCONTROLS_EXPORT __declspec( dllexport )
51  #else
52  #define SMESHCONTROLS_EXPORT __declspec( dllimport )
53  #endif
54 #else
55  #define SMESHCONTROLS_EXPORT
56 #endif
57 
58 class SMDS_MeshElement;
59 class SMDS_MeshFace;
60 class SMDS_MeshNode;
61 class SMDS_Mesh;
62 
63 class SMESHDS_Mesh;
64 class SMESHDS_SubMesh;
65 
66 class gp_Pnt;
67 
68 namespace SMESH{
69  namespace Controls{
70 
72  {
73  typedef std::vector<gp_XYZ>::size_type size_type;
74 
75  public:
77 
79 
80  TSequenceOfXYZ(size_type n, const gp_XYZ& t);
81 
82  TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ);
83 
84  template <class InputIterator>
85  TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd);
86 
87  ~TSequenceOfXYZ();
88 
89  TSequenceOfXYZ& operator=(const TSequenceOfXYZ& theSequenceOfXYZ);
90 
91  gp_XYZ& operator()(size_type n);
92 
93  const gp_XYZ& operator()(size_type n) const;
94 
95  void clear();
96 
97  void reserve(size_type n);
98 
99  void push_back(const gp_XYZ& v);
100 
101  size_type size() const;
102 
103  private:
104  std::vector<gp_XYZ> myArray;
105  };
106 
107  /*
108  Class : Functor
109  Description : Root of all Functors
110  */
112  {
113  public:
115  virtual void SetMesh( const SMDS_Mesh* theMesh ) = 0;
116  virtual SMDSAbs_ElementType GetType() const = 0;
117  };
118 
119  /*
120  Class : NumericalFunctor
121  Description : Root of all Functors returning numeric value
122  */
124  public:
126  virtual void SetMesh( const SMDS_Mesh* theMesh );
127  virtual double GetValue( long theElementId );
128  virtual double GetValue(const TSequenceOfXYZ& thePoints) { return -1.0;};
129  virtual SMDSAbs_ElementType GetType() const = 0;
130  virtual double GetBadRate( double Value, int nbNodes ) const = 0;
131  long GetPrecision() const;
132  void SetPrecision( const long thePrecision );
133 
134  bool GetPoints(const int theId,
135  TSequenceOfXYZ& theRes) const;
136  static bool GetPoints(const SMDS_MeshElement* theElem,
137  TSequenceOfXYZ& theRes);
138  protected:
142  };
143 
144 
145  /*
146  Class : Volume
147  Description : Functor calculating volume of 3D mesh element
148  */
150  public:
151  virtual double GetValue( long theElementId );
152  //virtual double GetValue( const TSequenceOfXYZ& thePoints );
153  virtual double GetBadRate( double Value, int nbNodes ) const;
154  virtual SMDSAbs_ElementType GetType() const;
155  };
156 
157 
158  /*
159  Class : SMESH_MinimumAngle
160  Description : Functor for calculation of minimum angle
161  */
163  public:
164  virtual double GetValue( const TSequenceOfXYZ& thePoints );
165  virtual double GetBadRate( double Value, int nbNodes ) const;
166  virtual SMDSAbs_ElementType GetType() const;
167  };
168 
169 
170  /*
171  Class : AspectRatio
172  Description : Functor for calculating aspect ratio
173  */
175  public:
176  virtual double GetValue( const TSequenceOfXYZ& thePoints );
177  virtual double GetBadRate( double Value, int nbNodes ) const;
178  virtual SMDSAbs_ElementType GetType() const;
179  };
180 
181 
182  /*
183  Class : AspectRatio3D
184  Description : Functor for calculating aspect ratio of 3D elems.
185  */
187  public:
188  virtual double GetValue( const TSequenceOfXYZ& thePoints );
189  virtual double GetBadRate( double Value, int nbNodes ) const;
190  virtual SMDSAbs_ElementType GetType() const;
191  };
192 
193 
194  /*
195  Class : Warping
196  Description : Functor for calculating warping
197  */
199  public:
200  virtual double GetValue( const TSequenceOfXYZ& thePoints );
201  virtual double GetBadRate( double Value, int nbNodes ) const;
202  virtual SMDSAbs_ElementType GetType() const;
203 
204  private:
205  double ComputeA( const gp_XYZ&, const gp_XYZ&, const gp_XYZ&, const gp_XYZ& ) const;
206  };
207 
208 
209  /*
210  Class : Taper
211  Description : Functor for calculating taper
212  */
214  public:
215  virtual double GetValue( const TSequenceOfXYZ& thePoints );
216  virtual double GetBadRate( double Value, int nbNodes ) const;
217  virtual SMDSAbs_ElementType GetType() const;
218  };
219 
220 
221  /*
222  Class : Skew
223  Description : Functor for calculating skew in degrees
224  */
226  public:
227  virtual double GetValue( const TSequenceOfXYZ& thePoints );
228  virtual double GetBadRate( double Value, int nbNodes ) const;
229  virtual SMDSAbs_ElementType GetType() const;
230  };
231 
232 
233  /*
234  Class : Area
235  Description : Functor for calculating area
236  */
238  public:
239  virtual double GetValue( const TSequenceOfXYZ& thePoints );
240  virtual double GetBadRate( double Value, int nbNodes ) const;
241  virtual SMDSAbs_ElementType GetType() const;
242  };
243 
244 
245  /*
246  Class : Length
247  Description : Functor for calculating length of edge
248  */
250  public:
251  virtual double GetValue( const TSequenceOfXYZ& thePoints );
252  virtual double GetBadRate( double Value, int nbNodes ) const;
253  virtual SMDSAbs_ElementType GetType() const;
254  };
255 
256  /*
257  Class : Length2D
258  Description : Functor for calculating length of edge
259  */
261  public:
262  virtual double GetValue( long theElementId );
263  virtual double GetBadRate( double Value, int nbNodes ) const;
264  virtual SMDSAbs_ElementType GetType() const;
265  struct Value{
266  double myLength;
267  long myPntId[2];
268  Value(double theLength, long thePntId1, long thePntId2);
269  bool operator<(const Value& x) const;
270  };
271  typedef std::set<Value> TValues;
272  void GetValues(TValues& theValues);
273  };
274  typedef boost::shared_ptr<Length2D> Length2DPtr;
275 
276  /*
277  Class : MultiConnection
278  Description : Functor for calculating number of faces conneted to the edge
279  */
281  public:
282  virtual double GetValue( long theElementId );
283  virtual double GetValue( const TSequenceOfXYZ& thePoints );
284  virtual double GetBadRate( double Value, int nbNodes ) const;
285  virtual SMDSAbs_ElementType GetType() const;
286  };
287 
288  /*
289  Class : MultiConnection2D
290  Description : Functor for calculating number of faces conneted to the edge
291  */
293  public:
294  virtual double GetValue( long theElementId );
295  virtual double GetValue( const TSequenceOfXYZ& thePoints );
296  virtual double GetBadRate( double Value, int nbNodes ) const;
297  virtual SMDSAbs_ElementType GetType() const;
298  struct Value{
299  long myPntId[2];
300  Value(long thePntId1, long thePntId2);
301  bool operator<(const Value& x) const;
302  };
303  typedef std::map<Value,int> MValues;
304 
305  void GetValues(MValues& theValues);
306  };
307  typedef boost::shared_ptr<MultiConnection2D> MultiConnection2DPtr;
308  /*
309  PREDICATES
310  */
311  /*
312  Class : Predicate
313  Description : Base class for all predicates
314  */
315  class SMESHCONTROLS_EXPORT Predicate: public virtual Functor{
316  public:
317  virtual bool IsSatisfy( long theElementId ) = 0;
318  virtual SMDSAbs_ElementType GetType() const = 0;
319  };
320 
321 
322  /*
323  Class : FreeBorders
324  Description : Predicate for free borders
325  */
327  public:
328  FreeBorders();
329  virtual void SetMesh( const SMDS_Mesh* theMesh );
330  virtual bool IsSatisfy( long theElementId );
331  virtual SMDSAbs_ElementType GetType() const;
332 
333  protected:
335  };
336 
337 
338  /*
339  Class : BadOrientedVolume
340  Description : Predicate bad oriented volumes
341  */
343  public:
345  virtual void SetMesh( const SMDS_Mesh* theMesh );
346  virtual bool IsSatisfy( long theElementId );
347  virtual SMDSAbs_ElementType GetType() const;
348 
349  protected:
351  };
352 
353 
354  /*
355  Class : FreeEdges
356  Description : Predicate for free Edges
357  */
358  class SMESHCONTROLS_EXPORT FreeEdges: public virtual Predicate{
359  public:
360  FreeEdges();
361  virtual void SetMesh( const SMDS_Mesh* theMesh );
362  virtual bool IsSatisfy( long theElementId );
363  virtual SMDSAbs_ElementType GetType() const;
364  static bool IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId );
365  typedef long TElemId;
366  struct Border{
368  TElemId myPntId[2];
369  Border(long theElemId, long thePntId1, long thePntId2);
370  bool operator<(const Border& x) const;
371  };
372  typedef std::set<Border> TBorders;
373  void GetBoreders(TBorders& theBorders);
374 
375  protected:
377  };
378  typedef boost::shared_ptr<FreeEdges> FreeEdgesPtr;
379 
380 
381  /*
382  Class : FreeNodes
383  Description : Predicate for free nodes
384  */
385  class SMESHCONTROLS_EXPORT FreeNodes: public virtual Predicate{
386  public:
387  FreeNodes();
388  virtual void SetMesh( const SMDS_Mesh* theMesh );
389  virtual bool IsSatisfy( long theNodeId );
390  virtual SMDSAbs_ElementType GetType() const;
391 
392  protected:
394  };
395 
396 
397  /*
398  Class : RangeOfIds
399  Description : Predicate for Range of Ids.
400  Range may be specified with two ways.
401  1. Using AddToRange method
402  2. With SetRangeStr method. Parameter of this method is a string
403  like as "1,2,3,50-60,63,67,70-"
404  */
406  {
407  public:
408  RangeOfIds();
409  virtual void SetMesh( const SMDS_Mesh* theMesh );
410  virtual bool IsSatisfy( long theNodeId );
411  virtual SMDSAbs_ElementType GetType() const;
412  virtual void SetType( SMDSAbs_ElementType theType );
413 
414  bool AddToRange( long theEntityId );
415  void GetRangeStr( TCollection_AsciiString& );
416  bool SetRangeStr( const TCollection_AsciiString& );
417 
418  protected:
420 
421  TColStd_SequenceOfInteger myMin;
422  TColStd_SequenceOfInteger myMax;
423  TColStd_MapOfInteger myIds;
424 
426  };
427 
428  typedef boost::shared_ptr<RangeOfIds> RangeOfIdsPtr;
429 
430 
431  /*
432  Class : Comparator
433  Description : Base class for comparators
434  */
436  public:
437  Comparator();
438  virtual ~Comparator();
439  virtual void SetMesh( const SMDS_Mesh* theMesh );
440  virtual void SetMargin(double theValue);
441  virtual void SetNumFunctor(NumericalFunctorPtr theFunct);
442  virtual bool IsSatisfy( long theElementId ) = 0;
443  virtual SMDSAbs_ElementType GetType() const;
444  double GetMargin();
445 
446  protected:
447  double myMargin;
449  };
450  typedef boost::shared_ptr<Comparator> ComparatorPtr;
451 
452 
453  /*
454  Class : LessThan
455  Description : Comparator "<"
456  */
457  class SMESHCONTROLS_EXPORT LessThan: public virtual Comparator{
458  public:
459  virtual bool IsSatisfy( long theElementId );
460  };
461 
462 
463  /*
464  Class : MoreThan
465  Description : Comparator ">"
466  */
467  class SMESHCONTROLS_EXPORT MoreThan: public virtual Comparator{
468  public:
469  virtual bool IsSatisfy( long theElementId );
470  };
471 
472 
473  /*
474  Class : EqualTo
475  Description : Comparator "="
476  */
477  class SMESHCONTROLS_EXPORT EqualTo: public virtual Comparator{
478  public:
479  EqualTo();
480  virtual bool IsSatisfy( long theElementId );
481  virtual void SetTolerance( double theTol );
482  virtual double GetTolerance();
483 
484  private:
485  double myToler;
486  };
487  typedef boost::shared_ptr<EqualTo> EqualToPtr;
488 
489 
490  /*
491  Class : LogicalNOT
492  Description : Logical NOT predicate
493  */
495  public:
496  LogicalNOT();
497  virtual ~LogicalNOT();
498  virtual bool IsSatisfy( long theElementId );
499  virtual void SetMesh( const SMDS_Mesh* theMesh );
500  virtual void SetPredicate(PredicatePtr thePred);
501  virtual SMDSAbs_ElementType GetType() const;
502 
503  private:
505  };
506  typedef boost::shared_ptr<LogicalNOT> LogicalNOTPtr;
507 
508 
509  /*
510  Class : LogicalBinary
511  Description : Base class for binary logical predicate
512  */
514  public:
515  LogicalBinary();
516  virtual ~LogicalBinary();
517  virtual void SetMesh( const SMDS_Mesh* theMesh );
518  virtual void SetPredicate1(PredicatePtr thePred);
519  virtual void SetPredicate2(PredicatePtr thePred);
520  virtual SMDSAbs_ElementType GetType() const;
521 
522  protected:
525  };
526  typedef boost::shared_ptr<LogicalBinary> LogicalBinaryPtr;
527 
528 
529  /*
530  Class : LogicalAND
531  Description : Logical AND
532  */
534  public:
535  virtual bool IsSatisfy( long theElementId );
536  };
537 
538 
539  /*
540  Class : LogicalOR
541  Description : Logical OR
542  */
544  public:
545  virtual bool IsSatisfy( long theElementId );
546  };
547 
548 
549  /*
550  Class : ManifoldPart
551  Description : Predicate for manifold part of mesh
552  */
554  public:
555 
556  /* internal class for algorithm uses */
557  class Link
558  {
559  public:
560  Link( SMDS_MeshNode* theNode1,
561  SMDS_MeshNode* theNode2 );
562  ~Link();
563 
564  bool IsEqual( const ManifoldPart::Link& theLink ) const;
565  bool operator<(const ManifoldPart::Link& x) const;
566 
569  };
570 
571  bool IsEqual( const ManifoldPart::Link& theLink1,
572  const ManifoldPart::Link& theLink2 );
573 
574  typedef std::set<ManifoldPart::Link> TMapOfLink;
575  typedef std::vector<SMDS_MeshFace*> TVectorOfFacePtr;
576  typedef std::vector<ManifoldPart::Link> TVectorOfLink;
577  typedef std::map<SMDS_MeshFace*,int> TDataMapFacePtrInt;
578  typedef std::map<ManifoldPart::Link,SMDS_MeshFace*> TDataMapOfLinkFacePtr;
579 
580  ManifoldPart();
581  ~ManifoldPart();
582  virtual void SetMesh( const SMDS_Mesh* theMesh );
583  // inoke when all parameters already set
584  virtual bool IsSatisfy( long theElementId );
585  virtual SMDSAbs_ElementType GetType() const;
586 
587  void SetAngleTolerance( const double theAngToler );
588  double GetAngleTolerance() const;
589  void SetIsOnlyManifold( const bool theIsOnly );
590  void SetStartElem( const long theStartElemId );
591 
592  private:
593  bool process();
594  bool findConnected( const TDataMapFacePtrInt& theAllFacePtrInt,
595  SMDS_MeshFace* theStartFace,
596  TMapOfLink& theNonManifold,
597  TColStd_MapOfInteger& theResFaces );
598  bool isInPlane( const SMDS_MeshFace* theFace1,
599  const SMDS_MeshFace* theFace2 );
600  void expandBoundary( TMapOfLink& theMapOfBoundary,
601  TVectorOfLink& theSeqOfBoundary,
602  TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
603  TMapOfLink& theNonManifold,
604  SMDS_MeshFace* theNextFace ) const;
605 
606  void getFacesByLink( const Link& theLink,
607  TVectorOfFacePtr& theFaces ) const;
608 
609  private:
611  TColStd_MapOfInteger myMapIds;
612  TColStd_MapOfInteger myMapBadGeomIds;
615  double myAngToler;
618 
619  };
620  typedef boost::shared_ptr<ManifoldPart> ManifoldPartPtr;
621 
622 
623  /*
624  Class : ElementsOnSurface
625  Description : Predicate elements that lying on indicated surface
626  (plane or cylinder)
627  */
629  public:
632  virtual void SetMesh( const SMDS_Mesh* theMesh );
633  virtual bool IsSatisfy( long theElementId );
634  virtual SMDSAbs_ElementType GetType() const;
635 
636  void SetTolerance( const double theToler );
637  double GetTolerance() const;
638  void SetSurface( const TopoDS_Shape& theShape,
639  const SMDSAbs_ElementType theType );
640  void SetUseBoundaries( bool theUse );
641  bool GetUseBoundaries() const { return myUseBoundaries; }
642 
643  private:
644  void process();
645  void process( const SMDS_MeshElement* theElem );
646  bool isOnSurface( const SMDS_MeshNode* theNode );
647 
648  private:
650  TColStd_MapOfInteger myIds;
652  //Handle(Geom_Surface) mySurf;
653  TopoDS_Face mySurf;
654  double myToler;
656  GeomAPI_ProjectPointOnSurf myProjector;
657  };
658 
659  typedef boost::shared_ptr<ElementsOnSurface> ElementsOnSurfacePtr;
660 
661 
662  /*
663  Class : ElementsOnShape
664  Description : Predicate elements that lying on indicated shape
665  (1D, 2D or 3D)
666  */
668  {
669  public:
670  ElementsOnShape();
671  ~ElementsOnShape();
672 
673  virtual void SetMesh (const SMDS_Mesh* theMesh);
674  virtual bool IsSatisfy (long theElementId);
675  virtual SMDSAbs_ElementType GetType() const;
676 
677  void SetTolerance (const double theToler);
678  double GetTolerance() const;
679  void SetAllNodes (bool theAllNodes);
680  bool GetAllNodes() const { return myAllNodesFlag; }
681  void SetShape (const TopoDS_Shape& theShape,
682  const SMDSAbs_ElementType theType);
683 
684  private:
685  void addShape (const TopoDS_Shape& theShape);
686  void process();
687  void process (const SMDS_MeshElement* theElem);
688 
689  private:
691  TColStd_MapOfInteger myIds;
693  TopoDS_Shape myShape;
694  double myToler;
696 
697  TopTools_MapOfShape myShapesMap;
698  TopAbs_ShapeEnum myCurShapeType; // type of current sub-shape
699  BRepClass3d_SolidClassifier myCurSC; // current SOLID
700  GeomAPI_ProjectPointOnSurf myCurProjFace; // current FACE
701  TopoDS_Face myCurFace; // current FACE
702  GeomAPI_ProjectPointOnCurve myCurProjEdge; // current EDGE
703  gp_Pnt myCurPnt; // current VERTEX
704  };
705 
706  typedef boost::shared_ptr<ElementsOnShape> ElementsOnShapePtr;
707 
708 
709  /*
710  Class : FreeFaces
711  Description : Predicate for free faces
712  */
713  class SMESHCONTROLS_EXPORT FreeFaces: public virtual Predicate{
714  public:
715  FreeFaces();
716  virtual void SetMesh( const SMDS_Mesh* theMesh );
717  virtual bool IsSatisfy( long theElementId );
718  virtual SMDSAbs_ElementType GetType() const;
719 
720  private:
722  };
723 
724  /*
725  Class : LinearOrQuadratic
726  Description : Predicate for free faces
727  */
729  public:
731  virtual void SetMesh( const SMDS_Mesh* theMesh );
732  virtual bool IsSatisfy( long theElementId );
733  void SetType( SMDSAbs_ElementType theType );
734  virtual SMDSAbs_ElementType GetType() const;
735 
736  private:
739  };
740  typedef boost::shared_ptr<LinearOrQuadratic> LinearOrQuadraticPtr;
741 
742  /*
743  Class : GroupColor
744  Description : Functor for check color of group to whic mesh element belongs to
745  */
747  public:
748  GroupColor();
749  virtual void SetMesh( const SMDS_Mesh* theMesh );
750  virtual bool IsSatisfy( long theElementId );
751  void SetType( SMDSAbs_ElementType theType );
752  virtual SMDSAbs_ElementType GetType() const;
753  void SetColorStr( const TCollection_AsciiString& );
754  void GetColorStr( TCollection_AsciiString& ) const;
755 
756  private:
757  typedef std::set< long > TIDs;
758 
759  Quantity_Color myColor;
762  };
763  typedef boost::shared_ptr<GroupColor> GroupColorPtr;
764 
765  /*
766  Class : ElemGeomType
767  Description : Predicate to check element geometry type
768  */
770  public:
771  ElemGeomType();
772  virtual void SetMesh( const SMDS_Mesh* theMesh );
773  virtual bool IsSatisfy( long theElementId );
774  void SetType( SMDSAbs_ElementType theType );
775  virtual SMDSAbs_ElementType GetType() const;
776  void SetGeomType( SMDSAbs_GeometryType theType );
777  virtual SMDSAbs_GeometryType GetGeomType() const;
778 
779  private:
783  };
784  typedef boost::shared_ptr<ElemGeomType> ElemGeomTypePtr;
785 
786  /*
787  FILTER
788  */
790  public:
791  Filter();
792  virtual ~Filter();
793  virtual void SetPredicate(PredicatePtr thePred);
794 
795  typedef std::vector<long> TIdSequence;
796 
797  virtual
798  void
799  GetElementsId( const SMDS_Mesh* theMesh,
800  TIdSequence& theSequence );
801 
802  static
803  void
804  GetElementsId( const SMDS_Mesh* theMesh,
805  PredicatePtr thePredicate,
806  TIdSequence& theSequence );
807 
808  protected:
810  };
811  };
812 };
813 
814 
815 #endif
TColStd_SequenceOfInteger myMin
TColStd_MapOfInteger myMapBadGeomIds
std::vector< SMDS_MeshFace * > TVectorOfFacePtr
virtual double GetValue(const TSequenceOfXYZ &thePoints)
std::vector< ManifoldPart::Link > TVectorOfLink
boost::shared_ptr< EqualTo > EqualToPtr
const SMDS_MeshElement * myCurrElement
boost::shared_ptr< ManifoldPart > ManifoldPartPtr
std::set< ManifoldPart::Link > TMapOfLink
boost::shared_ptr< NumericalFunctor > NumericalFunctorPtr
boost::shared_ptr< Comparator > ComparatorPtr
GeomAPI_ProjectPointOnCurve myCurProjEdge
boost::shared_ptr< MultiConnection2D > MultiConnection2DPtr
BRepClass3d_SolidClassifier myCurSC
std::map< SMDS_MeshFace *, int > TDataMapFacePtrInt
GeomAPI_ProjectPointOnSurf myProjector
SMDSAbs_ElementType
Type (node, edge, face or volume) of elements.
std::vector< long > TIdSequence
std::map< ManifoldPart::Link, SMDS_MeshFace * > TDataMapOfLinkFacePtr
TColStd_SequenceOfInteger myMax
boost::shared_ptr< ElemGeomType > ElemGeomTypePtr
boost::shared_ptr< LinearOrQuadratic > LinearOrQuadraticPtr
std::vector< gp_XYZ >::size_type size_type
boost::shared_ptr< FreeEdges > FreeEdgesPtr
boost::shared_ptr< LogicalBinary > LogicalBinaryPtr
boost::shared_ptr< Length2D > Length2DPtr
Standard_Boolean IsEqual(SMDS_MeshElementPtr theOne, SMDS_MeshElementPtr theTwo)
boost::shared_ptr< RangeOfIds > RangeOfIdsPtr
boost::shared_ptr< Predicate > PredicatePtr
Base class for elements.
boost::shared_ptr< LogicalNOT > LogicalNOTPtr
GeomAPI_ProjectPointOnSurf myCurProjFace
boost::shared_ptr< ElementsOnShape > ElementsOnShapePtr
TDataMapFacePtrInt myAllFacePtrIntDMap
#define SMESHCONTROLS_EXPORT
boost::shared_ptr< GroupColor > GroupColorPtr
SMDSAbs_GeometryType
boost::shared_ptr< ElementsOnSurface > ElementsOnSurfacePtr