]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_NumberOfSegments.hxx
Salome HOME
speed up NotifySubMeshesHypothesisModification()
[modules/smesh.git] / src / StdMeshers / StdMeshers_NumberOfSegments.hxx
1 // Copyright (C) 2007-2011  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
23 //  SMESH SMESH : implementaion of SMESH idl descriptions
24 //  File   : StdMeshers_NumberOfSegments.hxx
25 //           Moved here from SMESH_NumberOfSegments.hxx
26 //  Author : Paul RASCLE, EDF
27 //  Module : SMESH
28 //
29 #ifndef _SMESH_NUMBEROFSEGMENTS_HXX_
30 #define _SMESH_NUMBEROFSEGMENTS_HXX_
31
32 #include "SMESH_StdMeshers.hxx"
33
34 #include "SMESH_Hypothesis.hxx"
35 #include "Utils_SALOME_Exception.hxx"
36 #include <vector>
37
38 /*!
39  * \brief This class represents hypothesis for 1d algorithm
40  * 
41  * It provides parameters for subdivision an edge by various
42  * distribution types, considering the given number of resulting segments
43  */
44 class STDMESHERS_EXPORT StdMeshers_NumberOfSegments:
45   public SMESH_Hypothesis
46 {
47 public:
48   StdMeshers_NumberOfSegments(int hypId, int studyId, SMESH_Gen* gen);
49   virtual ~StdMeshers_NumberOfSegments();
50
51   // Builds point distribution according to passed function
52   const std::vector<double>& BuildDistributionExpr( const char*, int, int ) throw ( SALOME_Exception );
53   const std::vector<double>& BuildDistributionTab( const std::vector<double>&, int, int ) throw ( SALOME_Exception );
54
55   /*!
56    * \brief Set the number of segments
57     * \param segmentsNumber - must be greater than zero
58    */
59   void SetNumberOfSegments(int segmentsNumber)
60     throw (SALOME_Exception);
61
62   /*!
63    * \brief Get the number of segments
64    */
65   int GetNumberOfSegments() const;
66
67   /*!
68    * \brief This enumeration presents available types of distribution
69    */
70   enum DistrType
71   {
72     DT_Regular, //!< equidistant distribution
73     DT_Scale,   //!< scale distribution
74     DT_TabFunc, //!< distribution with density function presented by table
75     DT_ExprFunc //!< distribution with density function presented by expression
76   };
77
78   /*!
79    * \brief Set distribution type
80    */
81   void SetDistrType(DistrType typ)
82     throw (SALOME_Exception);
83
84   /*!
85    * \brief Get distribution type
86    */
87   DistrType GetDistrType() const;
88
89   /*!
90    * \brief Set scale factor for scale distribution
91    * \param scaleFactor - positive value different from 1
92    * 
93    * Throws SALOME_Exception if distribution type is not DT_Scale,
94    * or scaleFactor is not a positive value different from 1
95    */
96   virtual void SetScaleFactor(double scaleFactor)
97     throw (SALOME_Exception);
98
99   /*!
100    * \brief Get scale factor for scale distribution
101    * 
102    * Throws SALOME_Exception if distribution type is not DT_Scale
103    */
104   double GetScaleFactor() const
105     throw (SALOME_Exception);
106
107   /*!
108    * \brief Set table function for distribution DT_TabFunc
109     * \param table - this vector contains the pairs (parameter, value)
110    * following each by other, so the number of elements in the vector
111    * must be even. The parameters must be in range [0,1] and sorted in
112    * increase order. The values of function must be positive.
113    * 
114    * Throws SALOME_Exception if distribution type is not DT_TabFunc
115    */
116   void SetTableFunction(const std::vector<double>& table)
117     throw (SALOME_Exception);
118
119   /*!
120    * \brief Get table function for distribution DT_TabFunc
121    * 
122    * Throws SALOME_Exception if distribution type is not DT_TabFunc
123    */
124   const std::vector<double>& GetTableFunction() const
125     throw (SALOME_Exception);
126
127   /*!
128    * \brief Set expression function for distribution DT_ExprFunc
129     * \param expr - string containing the expression of the function
130     *               f(t), e.g. "sin(t)"
131    * 
132    * Throws SALOME_Exception if distribution type is not DT_ExprFunc
133    */
134   void SetExpressionFunction( const char* expr)
135     throw (SALOME_Exception);
136
137   /*!
138    * \brief Get expression function for distribution DT_ExprFunc
139    * 
140    * Throws SALOME_Exception if distribution type is not DT_ExprFunc
141    */
142   const char* GetExpressionFunction() const
143     throw (SALOME_Exception);
144
145   /*!
146    * \brief Set conversion mode. When it is 0, it means "exponent mode":
147    * the function of distribution of density is used as an exponent of 10, i,e, 10^f(t).
148    * When it is 1, it means "cut negative mode". The function of distribution is used as
149    * F(t), where F(t0)=f(t0), if f(t0)>=0, otherwise F(t0) = 0.
150    * This mode is sensible only when function distribution is used (DT_TabFunc or DT_ExprFunc)
151    * 
152    * Throws SALOME_Exception if distribution type is not functional
153    */
154   void SetConversionMode( int conv )
155     throw (SALOME_Exception);
156
157   /*!
158    * \brief Returns conversion mode
159    * 
160    * Throws SALOME_Exception if distribution type is not functional
161    */
162   int ConversionMode() const
163     throw (SALOME_Exception);
164
165   void SetReversedEdges( std::vector<int>& ids);
166
167   void SetObjectEntry( const char* entry ) { _objEntry = entry; }
168
169   const char* GetObjectEntry() { return _objEntry.c_str(); }
170
171   const std::vector<int>& GetReversedEdges() const { return _edgeIDs; }
172
173   /*!
174    * \brief Initialize number of segments by the mesh built on the geometry
175    * \param theMesh - the built mesh
176    * \param theShape - the geometry of interest
177    * \retval bool - true if parameter values have been successfully defined
178    */
179   virtual bool SetParametersByMesh(const SMESH_Mesh* theMesh, const TopoDS_Shape& theShape);
180
181   /*!
182    * \brief Initialize my parameter values by default parameters.
183    *  \retval bool - true if parameter values have been successfully defined
184    */
185   virtual bool SetParametersByDefaults(const TDefaults& dflts, const SMESH_Mesh* theMesh=0);
186
187   virtual std::ostream & SaveTo(std::ostream & save);
188   virtual std::istream & LoadFrom(std::istream & load);
189   friend std::ostream& operator << (std::ostream & save, StdMeshers_NumberOfSegments & hyp);
190   friend std::istream& operator >> (std::istream & load, StdMeshers_NumberOfSegments & hyp);
191
192 protected:
193   int                 _numberOfSegments; //!< an edge will be split on to this number of segments
194   DistrType           _distrType;        //!< the type of distribution of density function
195   double              _scaleFactor;      //!< the scale parameter for DT_Scale
196   std::vector<double> _table, _distr;    //!< the table for DT_TabFunc, a sequence of pairs of numbers
197   std::string         _func;             //!< the expression of the function for DT_ExprFunc
198   int                 _convMode;         //!< flag of conversion mode: 0=exponent, 1=cut negative
199   std::vector<int>    _edgeIDs;          //!< list of reversed edges ids
200   std::string         _objEntry;          //!< Entry of the main object to reverse edges
201 };
202
203 #endif