]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_SegmentLengthAroundVertex.cxx
Salome HOME
speed up NotifySubMeshesHypothesisModification()
[modules/smesh.git] / src / StdMeshers / StdMeshers_SegmentLengthAroundVertex.cxx
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_SegmentLengthAroundVertex.cxx
25 //  Module : SMESH
26 //
27 #include "StdMeshers_SegmentLengthAroundVertex.hxx"
28
29 #include "SMESH_Mesh.hxx"
30 #include "SMESH_Algo.hxx"
31 #include "SMDS_MeshNode.hxx"
32 #include "SMESHDS_Mesh.hxx"
33 #include "SMESHDS_SubMesh.hxx"
34 #include "SMESH_MesherHelper.hxx"
35
36 #include "utilities.h"
37
38 #include <BRepAdaptor_Curve.hxx>
39 #include <GCPnts_AbscissaPoint.hxx>
40 #include <TopTools_IndexedMapOfShape.hxx>
41 #include <TopoDS.hxx>
42 #include <TopoDS_Edge.hxx>
43
44 using namespace std;
45
46 //=============================================================================
47 /*!
48  *  
49  */
50 //=============================================================================
51
52 StdMeshers_SegmentLengthAroundVertex::StdMeshers_SegmentLengthAroundVertex
53                                        (int hypId, int studyId, SMESH_Gen * gen)
54   :SMESH_Hypothesis(hypId, studyId, gen)
55 {
56   _length = 1.;
57   _name = "SegmentLengthAroundVertex";
58   _param_algo_dim = 0; // is used by StdMeshers_SegmentAroundVertex_0D
59 }
60
61 //=============================================================================
62 /*!
63  *  
64  */
65 //=============================================================================
66
67 StdMeshers_SegmentLengthAroundVertex::~StdMeshers_SegmentLengthAroundVertex()
68 {
69 }
70
71 //=============================================================================
72 /*!
73  *  
74  */
75 //=============================================================================
76
77 void StdMeshers_SegmentLengthAroundVertex::SetLength(double length) throw(SALOME_Exception)
78 {
79   if (length <= 0)
80     throw SALOME_Exception(LOCALIZED("length must be positive"));
81   if (_length != length) {
82     _length = length;
83     NotifySubMeshesHypothesisModification();
84   }
85 }
86
87 //=============================================================================
88 /*!
89  *  
90  */
91 //=============================================================================
92
93 double StdMeshers_SegmentLengthAroundVertex::GetLength() const
94 {
95   return _length;
96 }
97
98 //=============================================================================
99 /*!
100  *  
101  */
102 //=============================================================================
103
104 ostream & StdMeshers_SegmentLengthAroundVertex::SaveTo(ostream & save)
105 {
106   save << this->_length;
107   return save;
108 }
109
110 //=============================================================================
111 /*!
112  *  
113  */
114 //=============================================================================
115
116 istream & StdMeshers_SegmentLengthAroundVertex::LoadFrom(istream & load)
117 {
118   bool isOK = true;
119   double a;
120   isOK = (load >> a);
121   if (isOK)
122     this->_length = a;
123   else
124     load.clear(ios::badbit | load.rdstate());
125   return load;
126 }
127
128 //=============================================================================
129 /*!
130  *  
131  */
132 //=============================================================================
133
134 ostream & operator <<(ostream & save, StdMeshers_SegmentLengthAroundVertex & hyp)
135 {
136   return hyp.SaveTo( save );
137 }
138
139 //=============================================================================
140 /*!
141  *  
142  */
143 //=============================================================================
144
145 istream & operator >>(istream & load, StdMeshers_SegmentLengthAroundVertex & hyp)
146 {
147   return hyp.LoadFrom( load );
148 }
149
150 //================================================================================
151 /*!
152  * \brief Initialize segment length by the mesh built on the geometry
153  * \param theMesh - the built mesh
154  * \param theShape - the geometry of interest
155  * \retval bool - true if parameter values have been successfully defined
156  */
157 //================================================================================
158
159 bool StdMeshers_SegmentLengthAroundVertex::SetParametersByMesh(const SMESH_Mesh*   theMesh,
160                                                                const TopoDS_Shape& theShape)
161 {
162   if ( !theMesh || theShape.IsNull() || theShape.ShapeType() != TopAbs_VERTEX )
163     return false;
164
165   SMESH_MeshEditor editor( const_cast<SMESH_Mesh*>( theMesh ) );
166   SMESH_MesherHelper helper( *editor.GetMesh() );
167
168   // get node built on theShape vertex
169   SMESHDS_Mesh* meshDS = editor.GetMeshDS();
170   SMESHDS_SubMesh* smV = meshDS->MeshElements( theShape );
171   if ( !smV || smV->NbNodes() == 0 )
172     return false;
173   const SMDS_MeshNode* vNode = smV->GetNodes()->next();
174
175   // calculate average length of segments sharing vNode
176
177   _length = 0.;
178   int nbSegs = 0;
179
180   SMDS_ElemIteratorPtr segIt = vNode->GetInverseElementIterator(SMDSAbs_Edge);
181   while ( segIt->more() ) {
182     const SMDS_MeshElement* seg = segIt->next();
183     // get geom edge
184     int shapeID = editor.FindShape( seg );
185     if (!shapeID) continue;
186     const TopoDS_Shape& s = meshDS->IndexToShape( shapeID );
187     if ( s.IsNull() || s.ShapeType() != TopAbs_EDGE ) continue;
188     const TopoDS_Edge& edge = TopoDS::Edge( s );
189     // params of edge ends
190     double u0 = helper.GetNodeU( edge, seg->GetNode(0) );
191     double u1 = helper.GetNodeU( edge, seg->GetNode(1) );
192     // length
193     BRepAdaptor_Curve AdaptCurve( edge );
194     _length += GCPnts_AbscissaPoint::Length( AdaptCurve, u0, u1);
195     nbSegs++;
196   }
197   
198   if ( nbSegs > 1 )
199     _length /= nbSegs;
200
201   return nbSegs;
202 }
203
204 //================================================================================
205 /*!
206  * \brief Initialize my parameter values by default parameters.
207  *  \retval bool - true if parameter values have been successfully defined
208  */
209 //================================================================================
210
211 bool StdMeshers_SegmentLengthAroundVertex::SetParametersByDefaults(const TDefaults&,
212                                                                    const SMESH_Mesh*)
213 {
214   return false;
215 }
216