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