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