]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_LocalLength.cxx
Salome HOME
cosmétique
[modules/smesh.git] / src / StdMeshers / StdMeshers_LocalLength.cxx
1 // Copyright (C) 2007-2020  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, 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 //  SMESH SMESH : implementation of SMESH idl descriptions
24 //  File   : StdMeshers_LocalLength.cxx
25 //           Moved here from SMESH_LocalLength.cxx
26 //  Author : Paul RASCLE, EDF
27 //  Module : SMESH
28 //
29 #include "StdMeshers_LocalLength.hxx"
30
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_Algo.hxx"
33
34 #include "utilities.h"
35
36 #include <BRep_Tool.hxx>
37 #include <GCPnts_AbscissaPoint.hxx>
38 #include <GeomAdaptor_Curve.hxx>
39 #include <Geom_Curve.hxx>
40 #include <TopExp.hxx>
41 #include <TopLoc_Location.hxx>
42 #include <TopTools_IndexedMapOfShape.hxx>
43 #include <TopoDS.hxx>
44 #include <TopoDS_Edge.hxx>
45 #include <Precision.hxx>
46
47 using namespace std;
48
49 //=============================================================================
50 /*!
51  *  
52  */
53 //=============================================================================
54
55 StdMeshers_LocalLength::StdMeshers_LocalLength(int hypId, SMESH_Gen * gen)
56   :SMESH_Hypothesis(hypId, gen)
57 {
58   _length = 1.;
59   _precision = Precision::Confusion();
60   _name = "LocalLength";
61   _param_algo_dim = 1; // is used by SMESH_Regular_1D
62 }
63
64 //=============================================================================
65 /*!
66  *  
67  */
68 //=============================================================================
69
70 StdMeshers_LocalLength::~StdMeshers_LocalLength()
71 {
72 }
73
74 //=============================================================================
75 /*!
76  *  
77  */
78 //=============================================================================
79
80 void StdMeshers_LocalLength::SetLength(double length)
81 {
82   double oldLength = _length;
83   if (length <= 0)
84     throw SALOME_Exception(LOCALIZED("length must be positive"));
85   _length = length;
86   const double precision = 1e-7;
87   if (fabs(oldLength - _length) > precision)
88     NotifySubMeshesHypothesisModification();
89 }
90
91 //=============================================================================
92 /*!
93  *  
94  */
95 //=============================================================================
96
97 double StdMeshers_LocalLength::GetLength() const
98 {
99   return _length;
100 }
101
102 //=============================================================================
103 /*!
104  *  
105  */
106 //=============================================================================
107 void StdMeshers_LocalLength::SetPrecision (double thePrecision)
108 {
109   double oldPrecision = _precision;
110   if (_precision < 0)
111     throw SALOME_Exception(LOCALIZED("precision cannot be negative"));
112   _precision = thePrecision;
113   const double precision = 1e-8;
114   if (fabs(oldPrecision - _precision) > precision)
115     NotifySubMeshesHypothesisModification();
116 }
117
118 //=============================================================================
119 /*!
120  *  
121  */
122 //=============================================================================
123 double StdMeshers_LocalLength::GetPrecision() const
124 {
125   return _precision;
126 }
127
128 //=============================================================================
129 /*!
130  *  
131  */
132 //=============================================================================
133
134 ostream & StdMeshers_LocalLength::SaveTo(ostream & save)
135 {
136   save << this->_length << " " << this->_precision;
137   return save;
138 }
139
140 //=============================================================================
141 /*!
142  *  
143  */
144 //=============================================================================
145
146 istream & StdMeshers_LocalLength::LoadFrom(istream & load)
147 {
148   bool isOK = true;
149   double a;
150
151   isOK = static_cast<bool>(load >> a);
152   if (isOK)
153     this->_length = a;
154   else
155     load.clear(ios::badbit | load.rdstate());
156
157   isOK = static_cast<bool>(load >> a);
158   if (isOK)
159     this->_precision = a;
160   else
161   {
162     load.clear(ios::badbit | load.rdstate());
163     // old format, without precision
164     _precision = 0.;
165   }
166
167   return load;
168 }
169
170 //=============================================================================
171 /*!
172  *  
173  */
174 //=============================================================================
175
176 ostream & operator <<(ostream & save, StdMeshers_LocalLength & hyp)
177 {
178   return hyp.SaveTo( save );
179 }
180
181 //=============================================================================
182 /*!
183  *  
184  */
185 //=============================================================================
186
187 istream & operator >>(istream & load, StdMeshers_LocalLength & hyp)
188 {
189   return hyp.LoadFrom( load );
190 }
191
192 //================================================================================
193 /*!
194  * \brief Initialize segment length by the mesh built on the geometry
195  * \param theMesh - the built mesh
196  * \param theShape - the geometry of interest
197  * \retval bool - true if parameter values have been successfully defined
198  */
199 //================================================================================
200
201 bool StdMeshers_LocalLength::SetParametersByMesh(const SMESH_Mesh*   theMesh,
202                                                  const TopoDS_Shape& theShape)
203 {
204   if ( !theMesh || theShape.IsNull() )
205     return false;
206
207   _length = 0.;
208
209   Standard_Real UMin, UMax;
210   TopLoc_Location L;
211
212   int nbEdges = 0;
213   TopTools_IndexedMapOfShape edgeMap;
214   TopExp::MapShapes( theShape, TopAbs_EDGE, edgeMap );
215   for ( int iE = 1; iE <= edgeMap.Extent(); ++iE )
216   {
217     const TopoDS_Edge& edge = TopoDS::Edge( edgeMap( iE ));
218     Handle(Geom_Curve) C = BRep_Tool::Curve( edge, L, UMin, UMax );
219     if ( C.IsNull() )
220       continue;
221     GeomAdaptor_Curve AdaptCurve(C, UMin, UMax);
222
223     vector< double > params;
224     SMESHDS_Mesh* aMeshDS = const_cast< SMESH_Mesh* >( theMesh )->GetMeshDS();
225     if ( SMESH_Algo::GetNodeParamOnEdge( aMeshDS, edge, params ))
226     {
227       for ( size_t i = 1; i < params.size(); ++i )
228         _length += GCPnts_AbscissaPoint::Length( AdaptCurve, params[ i-1 ], params[ i ]);
229       nbEdges += params.size() - 1;
230     }
231   }
232   if ( nbEdges )
233     _length /= nbEdges;
234
235   _precision = Precision::Confusion();
236
237   return nbEdges;
238 }
239 //================================================================================
240 /*!
241  * \brief Initialize my parameter values by default parameters.
242  *  \retval bool - true if parameter values have been successfully defined
243  */
244 //================================================================================
245
246 bool StdMeshers_LocalLength::SetParametersByDefaults(const TDefaults&  dflts,
247                                                      const SMESH_Mesh* /*theMesh*/)
248 {
249   return ( _length = dflts._elemLength );
250 }
251