SALOME - SMESH
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 
71  class SMESHCONTROLS_EXPORT TSequenceOfXYZ: public std::vector<gp_XYZ>
72  {
73  public:
75 
76  TSequenceOfXYZ(size_type n);
77 
78  TSequenceOfXYZ(size_type n, const value_type& t);
79 
80  TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ);
81 
82  template <class InputIterator>
83  TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd);
84 
85  TSequenceOfXYZ& operator=(const TSequenceOfXYZ& theSequenceOfXYZ);
86 
87  reference operator()(size_type n);
88 
89  const_reference operator()(size_type n) const;
90 
91  private:
92  reference operator[](size_type n);
93 
94  const_reference operator[](size_type n) const;
95  };
96 
97  /*
98  Class : Functor
99  Description : Root of all Functors
100  */
102  {
103  public:
105  virtual void SetMesh( const SMDS_Mesh* theMesh ) = 0;
106  virtual SMDSAbs_ElementType GetType() const = 0;
107  };
108 
109  /*
110  Class : NumericalFunctor
111  Description : Root of all Functors returning numeric value
112  */
114  public:
116  virtual void SetMesh( const SMDS_Mesh* theMesh );
117  virtual double GetValue( long theElementId );
118  virtual double GetValue(const TSequenceOfXYZ& thePoints) { return -1.0;};
119  virtual SMDSAbs_ElementType GetType() const = 0;
120  virtual double GetBadRate( double Value, int nbNodes ) const = 0;
121  long GetPrecision() const;
122  void SetPrecision( const long thePrecision );
123 
124  bool GetPoints(const int theId,
125  TSequenceOfXYZ& theRes) const;
126  static bool GetPoints(const SMDS_MeshElement* theElem,
127  TSequenceOfXYZ& theRes);
128  protected:
132  };
133 
134 
135  /*
136  Class : Volume
137  Description : Functor calculating volume of 3D mesh element
138  */
140  public:
141  virtual double GetValue( long theElementId );
142  //virtual double GetValue( const TSequenceOfXYZ& thePoints );
143  virtual double GetBadRate( double Value, int nbNodes ) const;
144  virtual SMDSAbs_ElementType GetType() const;
145  };
146 
147 
148  /*
149  Class : SMESH_MinimumAngle
150  Description : Functor for calculation of minimum angle
151  */
153  public:
154  virtual double GetValue( const TSequenceOfXYZ& thePoints );
155  virtual double GetBadRate( double Value, int nbNodes ) const;
156  virtual SMDSAbs_ElementType GetType() const;
157  };
158 
159 
160  /*
161  Class : AspectRatio
162  Description : Functor for calculating aspect ratio
163  */
165  public:
166  virtual double GetValue( const TSequenceOfXYZ& thePoints );
167  virtual double GetBadRate( double Value, int nbNodes ) const;
168  virtual SMDSAbs_ElementType GetType() const;
169  };
170 
171 
172  /*
173  Class : AspectRatio3D
174  Description : Functor for calculating aspect ratio of 3D elems.
175  */
177  public:
178  virtual double GetValue( const TSequenceOfXYZ& thePoints );
179  virtual double GetBadRate( double Value, int nbNodes ) const;
180  virtual SMDSAbs_ElementType GetType() const;
181  };
182 
183 
184  /*
185  Class : Warping
186  Description : Functor for calculating warping
187  */
189  public:
190  virtual double GetValue( const TSequenceOfXYZ& thePoints );
191  virtual double GetBadRate( double Value, int nbNodes ) const;
192  virtual SMDSAbs_ElementType GetType() const;
193 
194  private:
195  double ComputeA( const gp_XYZ&, const gp_XYZ&, const gp_XYZ&, const gp_XYZ& ) const;
196  };
197 
198 
199  /*
200  Class : Taper
201  Description : Functor for calculating taper
202  */
204  public:
205  virtual double GetValue( const TSequenceOfXYZ& thePoints );
206  virtual double GetBadRate( double Value, int nbNodes ) const;
207  virtual SMDSAbs_ElementType GetType() const;
208  };
209 
210 
211  /*
212  Class : Skew
213  Description : Functor for calculating skew in degrees
214  */
216  public:
217  virtual double GetValue( const TSequenceOfXYZ& thePoints );
218  virtual double GetBadRate( double Value, int nbNodes ) const;
219  virtual SMDSAbs_ElementType GetType() const;
220  };
221 
222 
223  /*
224  Class : Area
225  Description : Functor for calculating area
226  */
228  public:
229  virtual double GetValue( const TSequenceOfXYZ& thePoints );
230  virtual double GetBadRate( double Value, int nbNodes ) const;
231  virtual SMDSAbs_ElementType GetType() const;
232  };
233 
234 
235  /*
236  Class : Length
237  Description : Functor for calculating length of edge
238  */
240  public:
241  virtual double GetValue( const TSequenceOfXYZ& thePoints );
242  virtual double GetBadRate( double Value, int nbNodes ) const;
243  virtual SMDSAbs_ElementType GetType() const;
244  };
245 
246  /*
247  Class : Length2D
248  Description : Functor for calculating length of edge
249  */
251  public:
252  virtual double GetValue( long theElementId );
253  virtual double GetBadRate( double Value, int nbNodes ) const;
254  virtual SMDSAbs_ElementType GetType() const;
255  struct Value{
256  double myLength;
257  long myPntId[2];
258  Value(double theLength, long thePntId1, long thePntId2);
259  bool operator<(const Value& x) const;
260  };
261  typedef std::set<Value> TValues;
262  void GetValues(TValues& theValues);
263  };
264  typedef boost::shared_ptr<Length2D> Length2DPtr;
265 
266  /*
267  Class : MultiConnection
268  Description : Functor for calculating number of faces conneted to the edge
269  */
271  public:
272  virtual double GetValue( long theElementId );
273  virtual double GetValue( const TSequenceOfXYZ& thePoints );
274  virtual double GetBadRate( double Value, int nbNodes ) const;
275  virtual SMDSAbs_ElementType GetType() const;
276  };
277 
278  /*
279  Class : MultiConnection2D
280  Description : Functor for calculating number of faces conneted to the edge
281  */
283  public:
284  virtual double GetValue( long theElementId );
285  virtual double GetValue( const TSequenceOfXYZ& thePoints );
286  virtual double GetBadRate( double Value, int nbNodes ) const;
287  virtual SMDSAbs_ElementType GetType() const;
288  struct Value{
289  long myPntId[2];
290  Value(long thePntId1, long thePntId2);
291  bool operator<(const Value& x) const;
292  };
293  typedef std::map<Value,int> MValues;
294 
295  void GetValues(MValues& theValues);
296  };
297  typedef boost::shared_ptr<MultiConnection2D> MultiConnection2DPtr;
298  /*
299  PREDICATES
300  */
301  /*
302  Class : Predicate
303  Description : Base class for all predicates
304  */
305  class SMESHCONTROLS_EXPORT Predicate: public virtual Functor{
306  public:
307  virtual bool IsSatisfy( long theElementId ) = 0;
308  virtual SMDSAbs_ElementType GetType() const = 0;
309  };
310 
311 
312  /*
313  Class : FreeBorders
314  Description : Predicate for free borders
315  */
317  public:
318  FreeBorders();
319  virtual void SetMesh( const SMDS_Mesh* theMesh );
320  virtual bool IsSatisfy( long theElementId );
321  virtual SMDSAbs_ElementType GetType() const;
322 
323  protected:
325  };
326 
327 
328  /*
329  Class : BadOrientedVolume
330  Description : Predicate bad oriented volumes
331  */
333  public:
335  virtual void SetMesh( const SMDS_Mesh* theMesh );
336  virtual bool IsSatisfy( long theElementId );
337  virtual SMDSAbs_ElementType GetType() const;
338 
339  protected:
341  };
342 
343 
344  /*
345  Class : FreeEdges
346  Description : Predicate for free Edges
347  */
348  class SMESHCONTROLS_EXPORT FreeEdges: public virtual Predicate{
349  public:
350  FreeEdges();
351  virtual void SetMesh( const SMDS_Mesh* theMesh );
352  virtual bool IsSatisfy( long theElementId );
353  virtual SMDSAbs_ElementType GetType() const;
354  static bool IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId );
355  typedef long TElemId;
356  struct Border{
358  TElemId myPntId[2];
359  Border(long theElemId, long thePntId1, long thePntId2);
360  bool operator<(const Border& x) const;
361  };
362  typedef std::set<Border> TBorders;
363  void GetBoreders(TBorders& theBorders);
364 
365  protected:
367  };
368  typedef boost::shared_ptr<FreeEdges> FreeEdgesPtr;
369 
370 
371  /*
372  Class : FreeNodes
373  Description : Predicate for free nodes
374  */
375  class SMESHCONTROLS_EXPORT FreeNodes: public virtual Predicate{
376  public:
377  FreeNodes();
378  virtual void SetMesh( const SMDS_Mesh* theMesh );
379  virtual bool IsSatisfy( long theNodeId );
380  virtual SMDSAbs_ElementType GetType() const;
381 
382  protected:
384  };
385 
386 
387  /*
388  Class : RangeOfIds
389  Description : Predicate for Range of Ids.
390  Range may be specified with two ways.
391  1. Using AddToRange method
392  2. With SetRangeStr method. Parameter of this method is a string
393  like as "1,2,3,50-60,63,67,70-"
394  */
396  {
397  public:
398  RangeOfIds();
399  virtual void SetMesh( const SMDS_Mesh* theMesh );
400  virtual bool IsSatisfy( long theNodeId );
401  virtual SMDSAbs_ElementType GetType() const;
402  virtual void SetType( SMDSAbs_ElementType theType );
403 
404  bool AddToRange( long theEntityId );
405  void GetRangeStr( TCollection_AsciiString& );
406  bool SetRangeStr( const TCollection_AsciiString& );
407 
408  protected:
410 
411  TColStd_SequenceOfInteger myMin;
412  TColStd_SequenceOfInteger myMax;
413  TColStd_MapOfInteger myIds;
414 
416  };
417 
418  typedef boost::shared_ptr<RangeOfIds> RangeOfIdsPtr;
419 
420 
421  /*
422  Class : Comparator
423  Description : Base class for comparators
424  */
426  public:
427  Comparator();
428  virtual ~Comparator();
429  virtual void SetMesh( const SMDS_Mesh* theMesh );
430  virtual void SetMargin(double theValue);
431  virtual void SetNumFunctor(NumericalFunctorPtr theFunct);
432  virtual bool IsSatisfy( long theElementId ) = 0;
433  virtual SMDSAbs_ElementType GetType() const;
434  double GetMargin();
435 
436  protected:
437  double myMargin;
439  };
440  typedef boost::shared_ptr<Comparator> ComparatorPtr;
441 
442 
443  /*
444  Class : LessThan
445  Description : Comparator "<"
446  */
447  class SMESHCONTROLS_EXPORT LessThan: public virtual Comparator{
448  public:
449  virtual bool IsSatisfy( long theElementId );
450  };
451 
452 
453  /*
454  Class : MoreThan
455  Description : Comparator ">"
456  */
457  class SMESHCONTROLS_EXPORT MoreThan: public virtual Comparator{
458  public:
459  virtual bool IsSatisfy( long theElementId );
460  };
461 
462 
463  /*
464  Class : EqualTo
465  Description : Comparator "="
466  */
467  class SMESHCONTROLS_EXPORT EqualTo: public virtual Comparator{
468  public:
469  EqualTo();
470  virtual bool IsSatisfy( long theElementId );
471  virtual void SetTolerance( double theTol );
472  virtual double GetTolerance();
473 
474  private:
475  double myToler;
476  };
477  typedef boost::shared_ptr<EqualTo> EqualToPtr;
478 
479 
480  /*
481  Class : LogicalNOT
482  Description : Logical NOT predicate
483  */
485  public:
486  LogicalNOT();
487  virtual ~LogicalNOT();
488  virtual bool IsSatisfy( long theElementId );
489  virtual void SetMesh( const SMDS_Mesh* theMesh );
490  virtual void SetPredicate(PredicatePtr thePred);
491  virtual SMDSAbs_ElementType GetType() const;
492 
493  private:
495  };
496  typedef boost::shared_ptr<LogicalNOT> LogicalNOTPtr;
497 
498 
499  /*
500  Class : LogicalBinary
501  Description : Base class for binary logical predicate
502  */
504  public:
505  LogicalBinary();
506  virtual ~LogicalBinary();
507  virtual void SetMesh( const SMDS_Mesh* theMesh );
508  virtual void SetPredicate1(PredicatePtr thePred);
509  virtual void SetPredicate2(PredicatePtr thePred);
510  virtual SMDSAbs_ElementType GetType() const;
511 
512  protected:
515  };
516  typedef boost::shared_ptr<LogicalBinary> LogicalBinaryPtr;
517 
518 
519  /*
520  Class : LogicalAND
521  Description : Logical AND
522  */
524  public:
525  virtual bool IsSatisfy( long theElementId );
526  };
527 
528 
529  /*
530  Class : LogicalOR
531  Description : Logical OR
532  */
534  public:
535  virtual bool IsSatisfy( long theElementId );
536  };
537 
538 
539  /*
540  Class : ManifoldPart
541  Description : Predicate for manifold part of mesh
542  */
544  public:
545 
546  /* internal class for algorithm uses */
547  class Link
548  {
549  public:
550  Link( SMDS_MeshNode* theNode1,
551  SMDS_MeshNode* theNode2 );
552  ~Link();
553 
554  bool IsEqual( const ManifoldPart::Link& theLink ) const;
555  bool operator<(const ManifoldPart::Link& x) const;
556 
559  };
560 
561  bool IsEqual( const ManifoldPart::Link& theLink1,
562  const ManifoldPart::Link& theLink2 );
563 
564  typedef std::set<ManifoldPart::Link> TMapOfLink;
565  typedef std::vector<SMDS_MeshFace*> TVectorOfFacePtr;
566  typedef std::vector<ManifoldPart::Link> TVectorOfLink;
567  typedef std::map<SMDS_MeshFace*,int> TDataMapFacePtrInt;
568  typedef std::map<ManifoldPart::Link,SMDS_MeshFace*> TDataMapOfLinkFacePtr;
569 
570  ManifoldPart();
571  ~ManifoldPart();
572  virtual void SetMesh( const SMDS_Mesh* theMesh );
573  // inoke when all parameters already set
574  virtual bool IsSatisfy( long theElementId );
575  virtual SMDSAbs_ElementType GetType() const;
576 
577  void SetAngleTolerance( const double theAngToler );
578  double GetAngleTolerance() const;
579  void SetIsOnlyManifold( const bool theIsOnly );
580  void SetStartElem( const long theStartElemId );
581 
582  private:
583  bool process();
584  bool findConnected( const TDataMapFacePtrInt& theAllFacePtrInt,
585  SMDS_MeshFace* theStartFace,
586  TMapOfLink& theNonManifold,
587  TColStd_MapOfInteger& theResFaces );
588  bool isInPlane( const SMDS_MeshFace* theFace1,
589  const SMDS_MeshFace* theFace2 );
590  void expandBoundary( TMapOfLink& theMapOfBoundary,
591  TVectorOfLink& theSeqOfBoundary,
592  TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
593  TMapOfLink& theNonManifold,
594  SMDS_MeshFace* theNextFace ) const;
595 
596  void getFacesByLink( const Link& theLink,
597  TVectorOfFacePtr& theFaces ) const;
598 
599  private:
601  TColStd_MapOfInteger myMapIds;
602  TColStd_MapOfInteger myMapBadGeomIds;
605  double myAngToler;
608 
609  };
610  typedef boost::shared_ptr<ManifoldPart> ManifoldPartPtr;
611 
612 
613  /*
614  Class : ElementsOnSurface
615  Description : Predicate elements that lying on indicated surface
616  (plane or cylinder)
617  */
619  public:
622  virtual void SetMesh( const SMDS_Mesh* theMesh );
623  virtual bool IsSatisfy( long theElementId );
624  virtual SMDSAbs_ElementType GetType() const;
625 
626  void SetTolerance( const double theToler );
627  double GetTolerance() const;
628  void SetSurface( const TopoDS_Shape& theShape,
629  const SMDSAbs_ElementType theType );
630  void SetUseBoundaries( bool theUse );
631  bool GetUseBoundaries() const { return myUseBoundaries; }
632 
633  private:
634  void process();
635  void process( const SMDS_MeshElement* theElem );
636  bool isOnSurface( const SMDS_MeshNode* theNode );
637 
638  private:
640  TColStd_MapOfInteger myIds;
642  //Handle(Geom_Surface) mySurf;
643  TopoDS_Face mySurf;
644  double myToler;
646  GeomAPI_ProjectPointOnSurf myProjector;
647  };
648 
649  typedef boost::shared_ptr<ElementsOnSurface> ElementsOnSurfacePtr;
650 
651 
652  /*
653  Class : ElementsOnShape
654  Description : Predicate elements that lying on indicated shape
655  (1D, 2D or 3D)
656  */
658  {
659  public:
660  ElementsOnShape();
661  ~ElementsOnShape();
662 
663  virtual void SetMesh (const SMDS_Mesh* theMesh);
664  virtual bool IsSatisfy (long theElementId);
665  virtual SMDSAbs_ElementType GetType() const;
666 
667  void SetTolerance (const double theToler);
668  double GetTolerance() const;
669  void SetAllNodes (bool theAllNodes);
670  bool GetAllNodes() const { return myAllNodesFlag; }
671  void SetShape (const TopoDS_Shape& theShape,
672  const SMDSAbs_ElementType theType);
673 
674  private:
675  void addShape (const TopoDS_Shape& theShape);
676  void process();
677  void process (const SMDS_MeshElement* theElem);
678 
679  private:
681  TColStd_MapOfInteger myIds;
684  double myToler;
686 
687  TopTools_MapOfShape myShapesMap;
688  TopAbs_ShapeEnum myCurShapeType; // type of current sub-shape
689  BRepClass3d_SolidClassifier myCurSC; // current SOLID
690  GeomAPI_ProjectPointOnSurf myCurProjFace; // current FACE
691  TopoDS_Face myCurFace; // current FACE
692  GeomAPI_ProjectPointOnCurve myCurProjEdge; // current EDGE
693  gp_Pnt myCurPnt; // current VERTEX
694  };
695 
696  typedef boost::shared_ptr<ElementsOnShape> ElementsOnShapePtr;
697 
698 
699  /*
700  Class : FreeFaces
701  Description : Predicate for free faces
702  */
703  class SMESHCONTROLS_EXPORT FreeFaces: public virtual Predicate{
704  public:
705  FreeFaces();
706  virtual void SetMesh( const SMDS_Mesh* theMesh );
707  virtual bool IsSatisfy( long theElementId );
708  virtual SMDSAbs_ElementType GetType() const;
709 
710  private:
712  };
713 
714  /*
715  Class : LinearOrQuadratic
716  Description : Predicate for free faces
717  */
719  public:
721  virtual void SetMesh( const SMDS_Mesh* theMesh );
722  virtual bool IsSatisfy( long theElementId );
723  void SetType( SMDSAbs_ElementType theType );
724  virtual SMDSAbs_ElementType GetType() const;
725 
726  private:
729  };
730  typedef boost::shared_ptr<LinearOrQuadratic> LinearOrQuadraticPtr;
731 
732  /*
733  Class : GroupColor
734  Description : Functor for check color of group to whic mesh element belongs to
735  */
737  public:
738  GroupColor();
739  virtual void SetMesh( const SMDS_Mesh* theMesh );
740  virtual bool IsSatisfy( long theElementId );
741  void SetType( SMDSAbs_ElementType theType );
742  virtual SMDSAbs_ElementType GetType() const;
743  void SetColorStr( const TCollection_AsciiString& );
744  void GetColorStr( TCollection_AsciiString& ) const;
745 
746  private:
747  typedef std::set< long > TIDs;
748 
749  Quantity_Color myColor;
752  };
753  typedef boost::shared_ptr<GroupColor> GroupColorPtr;
754 
755  /*
756  Class : ElemGeomType
757  Description : Predicate to check element geometry type
758  */
760  public:
761  ElemGeomType();
762  virtual void SetMesh( const SMDS_Mesh* theMesh );
763  virtual bool IsSatisfy( long theElementId );
764  void SetType( SMDSAbs_ElementType theType );
765  virtual SMDSAbs_ElementType GetType() const;
766  void SetGeomType( SMDSAbs_GeometryType theType );
767  virtual SMDSAbs_GeometryType GetGeomType() const;
768 
769  private:
773  };
774  typedef boost::shared_ptr<ElemGeomType> ElemGeomTypePtr;
775 
776  /*
777  FILTER
778  */
780  public:
781  Filter();
782  virtual ~Filter();
783  virtual void SetPredicate(PredicatePtr thePred);
784 
785  typedef std::vector<long> TIdSequence;
786 
787  virtual
788  void
789  GetElementsId( const SMDS_Mesh* theMesh,
790  TIdSequence& theSequence );
791 
792  static
793  void
794  GetElementsId( const SMDS_Mesh* theMesh,
795  PredicatePtr thePredicate,
796  TIdSequence& theSequence );
797 
798  protected:
800  };
801  };
802 };
803 
804 
805 #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
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