]> SALOME platform Git repositories - modules/geom.git/blob - src/GEOMImpl/GEOMImpl_MeasureDriver.cxx
Salome HOME
54989d8e2c2cd65090120b1523bbc87e98f757e6
[modules/geom.git] / src / GEOMImpl / GEOMImpl_MeasureDriver.cxx
1 // Copyright (C) 2007-2013  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 #include <Standard_Stream.hxx>
24
25 #include <GEOMImpl_MeasureDriver.hxx>
26 #include <GEOMImpl_IMeasure.hxx>
27 #include <GEOMImpl_Types.hxx>
28
29 #include <GEOM_Function.hxx>
30
31 #include <GEOMUtils.hxx>
32
33 #include <BRep_Tool.hxx>
34 #include <BRepTools.hxx>
35 #include <BRepGProp.hxx>
36 #include <BRepBuilderAPI_Copy.hxx>
37 #include <BRepBuilderAPI_MakeEdge.hxx>
38 #include <BRepBuilderAPI_MakeVertex.hxx>
39 #include <BRepPrimAPI_MakeBox.hxx>
40
41 #include <TopAbs.hxx>
42 #include <TopoDS.hxx>
43 #include <TopoDS_Shape.hxx>
44 #include <TopoDS_Edge.hxx>
45 #include <TopoDS_Wire.hxx>
46 #include <TopExp.hxx>
47 #include <TopExp_Explorer.hxx>
48 #include <TopTools_IndexedMapOfShape.hxx>
49
50 #include <GProp_GProps.hxx>
51 #include <GeomLProp_SLProps.hxx>
52 #include <Geom_Surface.hxx>
53
54 #include <Bnd_Box.hxx>
55 #include <BRepBndLib.hxx>
56 #include <BRepAdaptor_Surface.hxx>
57 #include <ShapeAnalysis_Surface.hxx>
58
59 #include <gp_Pnt.hxx>
60 #include <Precision.hxx>
61 #include <Standard_NullObject.hxx>
62 #include <StdFail_NotDone.hxx>
63
64 //=======================================================================
65 //function : GetID
66 //purpose  :
67 //======================================================================= 
68 const Standard_GUID& GEOMImpl_MeasureDriver::GetID()
69 {
70   static Standard_GUID aMeasureDriver("FF1BBB65-5D14-4df2-980B-3A668264EA16");
71   return aMeasureDriver; 
72 }
73
74
75 //=======================================================================
76 //function : GEOMImpl_MeasureDriver
77 //purpose  : 
78 //=======================================================================
79 GEOMImpl_MeasureDriver::GEOMImpl_MeasureDriver() 
80 {
81 }
82
83 //=======================================================================
84 //function : Execute
85 //purpose  :
86 //======================================================================= 
87 Standard_Integer GEOMImpl_MeasureDriver::Execute(TFunction_Logbook& log) const
88 {
89   if (Label().IsNull()) return 0;    
90   Handle(GEOM_Function) aFunction = GEOM_Function::GetFunction(Label());
91
92   GEOMImpl_IMeasure aCI (aFunction);
93   Standard_Integer aType = aFunction->GetType();
94
95   TopoDS_Shape aShape;
96
97   if (aType == CDG_MEASURE)
98   {
99     Handle(GEOM_Function) aRefBase = aCI.GetBase();
100     TopoDS_Shape aShapeBase = aRefBase->GetValue();
101     if (aShapeBase.IsNull())
102       Standard_NullObject::Raise("Shape for centre of mass calculation is null");
103
104     gp_Ax3 aPos = GEOMUtils::GetPosition(aShapeBase);
105     gp_Pnt aCenterMass = aPos.Location();
106     aShape = BRepBuilderAPI_MakeVertex(aCenterMass).Shape();
107   }
108   else if (aType == BND_BOX_MEASURE)
109   {
110     Handle(GEOM_Function) aRefBase = aCI.GetBase();
111     TopoDS_Shape aShapeBase = aRefBase->GetValue();
112     if (aShapeBase.IsNull())
113       Standard_NullObject::Raise("Shape for bounding box calculation is null");
114
115     BRepBuilderAPI_Copy aCopyTool (aShapeBase);
116     if (!aCopyTool.IsDone())
117       Standard_NullObject::Raise("Shape for bounding box calculation is bad");
118
119     aShapeBase = aCopyTool.Shape();
120
121     // remove triangulation to obtain more exact boundaries
122     BRepTools::Clean(aShapeBase);
123
124     Bnd_Box B;
125     BRepBndLib::Add(aShapeBase, B);
126
127     Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
128     B.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
129
130     if (Xmax - Xmin < Precision::Confusion()) Xmax += Precision::Confusion();
131     if (Ymax - Ymin < Precision::Confusion()) Ymax += Precision::Confusion();
132     if (Zmax - Zmin < Precision::Confusion()) Zmax += Precision::Confusion();
133
134     gp_Pnt P1 (Xmin, Ymin, Zmin);
135     gp_Pnt P2 (Xmax, Ymax, Zmax);
136
137     BRepPrimAPI_MakeBox MB (P1, P2);
138     MB.Build();
139     if (!MB.IsDone()) StdFail_NotDone::Raise("Bounding box cannot be computed from the given shape");
140     aShape = MB.Shape();
141   }
142   else if (aType == VERTEX_BY_INDEX)
143   {
144     Handle(GEOM_Function) aRefBase = aCI.GetBase();
145     TopoDS_Shape aShapeBase = aRefBase->GetValue();
146     if (aShapeBase.IsNull()) {
147       Standard_NullObject::Raise("Shape for centre of mass calculation is null");
148     }
149
150     int index = aCI.GetIndex();
151     gp_Pnt aVertex;
152
153     if (aShapeBase.ShapeType() == TopAbs_VERTEX) {
154       if ( index != 1 )
155         Standard_NullObject::Raise("Vertex index is out of range");
156       else
157         aVertex = BRep_Tool::Pnt(TopoDS::Vertex(aShapeBase));
158     } else if (aShapeBase.ShapeType() == TopAbs_EDGE) {
159       TopoDS_Vertex aV1, aV2;
160       TopoDS_Edge anEdgeE = TopoDS::Edge(aShapeBase);
161       
162       TopExp::Vertices(anEdgeE, aV1, aV2);
163       gp_Pnt aP1 = BRep_Tool::Pnt(aV1);
164       gp_Pnt aP2 = BRep_Tool::Pnt(aV2);
165
166       if (index < 0 || index > 1)
167         Standard_NullObject::Raise("Vertex index is out of range");
168
169       if ( ( anEdgeE.Orientation() == TopAbs_FORWARD && index == 0 ) ||
170            ( anEdgeE.Orientation() == TopAbs_REVERSED && index == 1 ) )
171         aVertex = aP1;
172       else
173       aVertex = aP2;
174     } else if (aShapeBase.ShapeType() == TopAbs_WIRE) {
175       TopTools_IndexedMapOfShape anEdgeShapes;
176       TopTools_IndexedMapOfShape aVertexShapes;
177       TopoDS_Vertex aV1, aV2;
178       TopoDS_Wire aWire = TopoDS::Wire(aShapeBase);
179       TopExp_Explorer exp (aWire, TopAbs_EDGE);
180       for (; exp.More(); exp.Next()) {
181         anEdgeShapes.Add(exp.Current());
182         TopoDS_Edge E = TopoDS::Edge(exp.Current());
183         TopExp::Vertices(E, aV1, aV2);
184         if ( aVertexShapes.Extent() == 0)
185           aVertexShapes.Add(aV1);
186         if ( !aV1.IsSame( aVertexShapes(aVertexShapes.Extent()) ) )
187           aVertexShapes.Add(aV1);
188         if ( !aV2.IsSame( aVertexShapes(aVertexShapes.Extent()) ) )
189           aVertexShapes.Add(aV2);
190       }
191
192       if (index < 0 || index > aVertexShapes.Extent())
193         Standard_NullObject::Raise("Vertex index is out of range");
194
195       if (aWire.Orientation() == TopAbs_FORWARD)
196         aVertex = BRep_Tool::Pnt(TopoDS::Vertex(aVertexShapes(index+1)));
197       else
198         aVertex = BRep_Tool::Pnt(TopoDS::Vertex(aVertexShapes(aVertexShapes.Extent() - index)));
199     } else {
200       Standard_NullObject::Raise("Shape for vertex calculation is not an edge or wire");
201     }
202
203     aShape = BRepBuilderAPI_MakeVertex(aVertex).Shape();
204   }
205   else if (aType == VECTOR_FACE_NORMALE)
206   {
207     // Face
208     Handle(GEOM_Function) aRefBase = aCI.GetBase();
209     TopoDS_Shape aShapeBase = aRefBase->GetValue();
210     if (aShapeBase.IsNull()) {
211       Standard_NullObject::Raise("Face for normale calculation is null");
212     }
213     if (aShapeBase.ShapeType() != TopAbs_FACE) {
214       Standard_NullObject::Raise("Shape for normale calculation is not a face");
215     }
216     TopoDS_Face aFace = TopoDS::Face(aShapeBase);
217
218     // Point
219     gp_Pnt p1 (0,0,0);
220
221     Handle(GEOM_Function) aPntFunc = aCI.GetPoint();
222     if (!aPntFunc.IsNull())
223     {
224       TopoDS_Shape anOptPnt = aPntFunc->GetValue();
225       if (anOptPnt.IsNull())
226         Standard_NullObject::Raise("Invalid shape given for point argument");
227       p1 = BRep_Tool::Pnt(TopoDS::Vertex(anOptPnt));
228     }
229     else
230     {
231       gp_Ax3 aPos = GEOMUtils::GetPosition(aFace);
232       p1 = aPos.Location();
233     }
234
235     // Point parameters on surface
236     Handle(Geom_Surface) aSurf = BRep_Tool::Surface(aFace);
237     Handle(ShapeAnalysis_Surface) aSurfAna = new ShapeAnalysis_Surface (aSurf);
238     gp_Pnt2d pUV = aSurfAna->ValueOfUV(p1, Precision::Confusion());
239
240     // Normal direction
241     gp_Vec Vec1,Vec2;
242     BRepAdaptor_Surface SF (aFace);
243     SF.D1(pUV.X(), pUV.Y(), p1, Vec1, Vec2);
244     if (Vec1.Magnitude() < Precision::Confusion()) {
245       gp_Vec tmpV;
246       gp_Pnt tmpP;
247       SF.D1(pUV.X(), pUV.Y()-0.1, tmpP, Vec1, tmpV);
248     }
249     else if (Vec2.Magnitude() < Precision::Confusion()) {
250       gp_Vec tmpV;
251       gp_Pnt tmpP;
252       SF.D1(pUV.X()-0.1, pUV.Y(), tmpP, tmpV, Vec2);
253     }
254
255     gp_Vec V = Vec1.Crossed(Vec2);
256     Standard_Real mod = V.Magnitude();
257     if (mod < Precision::Confusion())
258       Standard_NullObject::Raise("Normal vector of a face has null magnitude");
259
260     // Set length of normal vector to average radius of curvature
261     Standard_Real radius = 0.0;
262     GeomLProp_SLProps aProperties (aSurf, pUV.X(), pUV.Y(), 2, Precision::Confusion());
263     if (aProperties.IsCurvatureDefined()) {
264       Standard_Real radius1 = Abs(aProperties.MinCurvature());
265       Standard_Real radius2 = Abs(aProperties.MaxCurvature());
266       if (Abs(radius1) > Precision::Confusion()) {
267         radius = 1.0 / radius1;
268         if (Abs(radius2) > Precision::Confusion()) {
269           radius = (radius + 1.0 / radius2) / 2.0;
270         }
271       }
272       else {
273         if (Abs(radius2) > Precision::Confusion()) {
274           radius = 1.0 / radius2;
275         }
276       }
277     }
278
279     // Set length of normal vector to average dimension of the face
280     // (only if average radius of curvature is not appropriate)
281     if (radius < Precision::Confusion()) {
282         Bnd_Box B;
283         Standard_Real Xmin, Xmax, Ymin, Ymax, Zmin, Zmax;
284         BRepBndLib::Add(aFace, B);
285         B.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
286         radius = ((Xmax - Xmin) + (Ymax - Ymin) + (Zmax - Zmin)) / 3.0;
287     }
288
289     if (radius < Precision::Confusion())
290       radius = 1.0;
291
292     V *= radius / mod;
293
294     // consider the face orientation
295     if (aFace.Orientation() == TopAbs_REVERSED ||
296         aFace.Orientation() == TopAbs_INTERNAL) {
297       V = - V;
298     }
299
300     // Edge
301     gp_Pnt p2 = p1.Translated(V);
302     BRepBuilderAPI_MakeEdge aBuilder (p1, p2);
303     if (!aBuilder.IsDone())
304       Standard_NullObject::Raise("Vector construction failed");
305     aShape = aBuilder.Shape();
306   }
307   else {
308   }
309
310   if (aShape.IsNull()) return 0;
311
312   aFunction->SetValue(aShape);
313
314   log.SetTouched(Label()); 
315
316   return 1;
317 }
318
319
320 //=======================================================================
321 //function :  GEOMImpl_MeasureDriver_Type_
322 //purpose  :
323 //======================================================================= 
324 Standard_EXPORT Handle_Standard_Type& GEOMImpl_MeasureDriver_Type_()
325 {
326
327   static Handle_Standard_Type aType1 = STANDARD_TYPE(TFunction_Driver);
328   if ( aType1.IsNull()) aType1 = STANDARD_TYPE(TFunction_Driver);
329   static Handle_Standard_Type aType2 = STANDARD_TYPE(MMgt_TShared);
330   if ( aType2.IsNull()) aType2 = STANDARD_TYPE(MMgt_TShared); 
331   static Handle_Standard_Type aType3 = STANDARD_TYPE(Standard_Transient);
332   if ( aType3.IsNull()) aType3 = STANDARD_TYPE(Standard_Transient);
333  
334
335   static Handle_Standard_Transient _Ancestors[]= {aType1,aType2,aType3,NULL};
336   static Handle_Standard_Type _aType = new Standard_Type("GEOMImpl_MeasureDriver",
337                                                          sizeof(GEOMImpl_MeasureDriver),
338                                                          1,
339                                                          (Standard_Address)_Ancestors,
340                                                          (Standard_Address)NULL);
341
342   return _aType;
343 }
344
345 //=======================================================================
346 //function : DownCast
347 //purpose  :
348 //======================================================================= 
349 const Handle(GEOMImpl_MeasureDriver) Handle(GEOMImpl_MeasureDriver)::DownCast(const Handle(Standard_Transient)& AnObject)
350 {
351   Handle(GEOMImpl_MeasureDriver) _anOtherObject;
352
353   if (!AnObject.IsNull()) {
354      if (AnObject->IsKind(STANDARD_TYPE(GEOMImpl_MeasureDriver))) {
355        _anOtherObject = Handle(GEOMImpl_MeasureDriver)((Handle(GEOMImpl_MeasureDriver)&)AnObject);
356      }
357   }
358
359   return _anOtherObject ;
360 }