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