Salome HOME
Fix for issue #1000
[modules/shaper.git] / src / GeomAlgoAPI / GeomAlgoAPI_Prism.cpp
1 // Copyright (C) 2014-20xx CEA/DEN, EDF R&D
2
3 // File:        GeomAlgoAPI_Prism.cpp
4 // Created:     5 May 2015
5 // Author:      Dmitry Bobylev
6
7 #include <GeomAlgoAPI_Prism.h>
8
9 #include <GeomAPI_Face.h>
10 #include <GeomAPI_Pln.h>
11 #include <GeomAPI_Pnt.h>
12 #include <GeomAPI_ShapeExplorer.h>
13 #include <GeomAPI_XYZ.h>
14 #include <GeomAlgoAPI_DFLoader.h>
15 #include <GeomAlgoAPI_FaceBuilder.h>
16 #include <GeomAlgoAPI_MakeShapeList.h>
17 #include <GeomAlgoAPI_ShapeTools.h>
18
19 #include <Bnd_Box.hxx>
20 #include <BRep_Builder.hxx>
21 #include <BRep_Tool.hxx>
22 #include <BRepAlgoAPI_Cut.hxx>
23 #include <BRepBndLib.hxx>
24 #include <BRepBuilderAPI_MakeEdge.hxx>
25 #include <BRepBuilderAPI_MakeShape.hxx>
26 #include <BRepBuilderAPI_MakeWire.hxx>
27 #include <BRepCheck_Analyzer.hxx>
28 #include <BRepExtrema_ExtCF.hxx>
29 #include <BRepFeat_MakePrism.hxx>
30 #include <BRepGProp.hxx>
31 #include <BRepOffsetAPI_MakePipe.hxx>
32 #include <BRepTools.hxx>
33 #include <Geom_Plane.hxx>
34 #include <gp_Pln.hxx>
35 #include <GProp_GProps.hxx>
36 #include <IntAna_IntConicQuad.hxx>
37 #include <IntAna_Quadric.hxx>
38 #include <TCollection_AsciiString.hxx>
39 #include <TopExp_Explorer.hxx>
40 #include <TopoDS.hxx>
41 #include <TopoDS_Edge.hxx>
42 #include <TopoDS_Shell.hxx>
43 #include <TopoDS_Solid.hxx>
44 #include <TopoDS_Wire.hxx>
45 #include <TopTools_ListIteratorOfListOfShape.hxx>
46
47 //=================================================================================================
48 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(std::shared_ptr<GeomAPI_Shape> theBasis,
49                                      double                         theToSize,
50                                      double                         theFromSize)
51 : myDone(false)
52 {
53   build(theBasis, std::shared_ptr<GeomAPI_Shape>(), theToSize, std::shared_ptr<GeomAPI_Shape>(), theFromSize);
54 }
55
56 //=================================================================================================
57 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(std::shared_ptr<GeomAPI_Shape> theBasis,
58                                      std::shared_ptr<GeomAPI_Shape> theToShape,
59                                      double                         theToSize,
60                                      std::shared_ptr<GeomAPI_Shape> theFromShape,
61                                      double                         theFromSize)
62 : myDone(false)
63 {
64   build(theBasis, theToShape, theToSize, theFromShape, theFromSize);
65 }
66
67 //=================================================================================================
68 void GeomAlgoAPI_Prism::build(const std::shared_ptr<GeomAPI_Shape>& theBasis,
69                               const std::shared_ptr<GeomAPI_Shape>& theToShape,
70                               double                                theToSize,
71                               const std::shared_ptr<GeomAPI_Shape>& theFromShape,
72                               double                                theFromSize)
73 {
74   if(!theBasis ||
75     (((!theFromShape && !theToShape) || (theFromShape && theToShape && theFromShape->isEqual(theToShape)))
76     && (theFromSize == -theToSize))) {
77     return;
78   }
79
80   // If bounding faces was not set creating them.
81   std::shared_ptr<GeomAPI_Face> aBaseFace;
82   if(theBasis->shapeType() == GeomAPI_Shape::FACE) {
83     aBaseFace = std::shared_ptr<GeomAPI_Face>(new GeomAPI_Face(theBasis));
84   } else if(theBasis->shapeType() == GeomAPI_Shape::SHELL){
85     GeomAPI_ShapeExplorer anExp(theBasis, GeomAPI_Shape::FACE);
86     if(anExp.more()) {
87       std::shared_ptr<GeomAPI_Shape> aFaceOnShell = anExp.current();
88       aBaseFace = std::shared_ptr<GeomAPI_Face>(new GeomAPI_Face(aFaceOnShell));
89     }
90   }
91   if(!aBaseFace.get()) {
92     return;
93   }
94   std::shared_ptr<GeomAPI_Pln>   aBasePln = aBaseFace->getPlane();
95   std::shared_ptr<GeomAPI_Dir>   aBaseDir = aBasePln->direction();
96   std::shared_ptr<GeomAPI_Pnt>   aBaseLoc = aBasePln->location();
97   std::shared_ptr<GeomAPI_Shape> aBasePlane = GeomAlgoAPI_FaceBuilder::planarFace(aBaseLoc, aBaseDir);
98
99   std::shared_ptr<GeomAPI_Shape> aBoundingFromShape = theFromShape ? theFromShape : aBasePlane;
100   std::shared_ptr<GeomAPI_Shape> aBoundingToShape   = theToShape   ? theToShape   : aBasePlane;
101
102   // Moving bounding faces according to "from" and "to" sizes.
103   std::shared_ptr<GeomAPI_Face> aFromFace(new GeomAPI_Face(aBoundingFromShape));
104   std::shared_ptr<GeomAPI_Pln>  aFromPln = aFromFace->getPlane();
105   std::shared_ptr<GeomAPI_Pnt>  aFromLoc = aFromPln->location();
106   std::shared_ptr<GeomAPI_Dir>  aFromDir = aFromPln->direction();
107
108   std::shared_ptr<GeomAPI_Face> aToFace(new GeomAPI_Face(aBoundingToShape));
109   std::shared_ptr<GeomAPI_Pln>  aToPln = aToFace->getPlane();
110   std::shared_ptr<GeomAPI_Pnt>  aToLoc = aToPln->location();
111   std::shared_ptr<GeomAPI_Dir>  aToDir = aToPln->direction();
112
113   bool aSign = aFromLoc->xyz()->dot(aBaseDir->xyz()) > aToLoc->xyz()->dot(aBaseDir->xyz());
114
115   std::shared_ptr<GeomAPI_Pnt> aFromPnt(new GeomAPI_Pnt(aFromLoc->xyz()->added(aBaseDir->xyz()->multiplied(
116                                                         aSign ? theFromSize : -theFromSize))));
117   aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
118
119   std::shared_ptr<GeomAPI_Pnt> aToPnt(new GeomAPI_Pnt(aToLoc->xyz()->added(aBaseDir->xyz()->multiplied(
120                                                       aSign ? -theToSize : theToSize))));
121   aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
122
123   if(theBasis->shapeType() == GeomAPI_Shape::FACE) {
124     TopoDS_Face aBasis = TopoDS::Face(aBaseFace->impl<TopoDS_Shape>());
125     const gp_Dir& aNormal = aBaseDir->impl<gp_Dir>();
126     BRepFeat_MakePrism* aBuilder = new BRepFeat_MakePrism(aBasis, aBasis, aBasis, aNormal, 2, Standard_True);
127
128     if(aBuilder) {
129       const TopoDS_Shape& aFromShape = aBoundingFromShape->impl<TopoDS_Shape>();
130       const TopoDS_Shape& aToShape   = aBoundingToShape->impl<TopoDS_Shape>();
131       aBuilder->Perform(aFromShape, aToShape);
132       myDone = aBuilder->IsDone() == Standard_True;
133       if(myDone){
134         TopoDS_Shape aResult = aBuilder->Shape();
135         TopExp_Explorer anExp(aResult, TopAbs_SOLID);
136         if(!anExp.More()) {
137           return;
138         }
139         if(aResult.ShapeType() == TopAbs_COMPOUND) {
140           aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
141         }
142         // fill data map to keep correct orientation of sub-shapes
143         myMap.reset(new GeomAPI_DataMapOfShapeShape);
144         for (TopExp_Explorer Exp(aResult,TopAbs_FACE); Exp.More(); Exp.Next()) {
145           std::shared_ptr<GeomAPI_Shape> aCurrentShape(new GeomAPI_Shape());
146           aCurrentShape->setImpl(new TopoDS_Shape(Exp.Current()));
147           myMap->bind(aCurrentShape, aCurrentShape);
148         }
149         myShape.reset(new GeomAPI_Shape);
150         myShape->setImpl(new TopoDS_Shape(aResult));
151         std::shared_ptr<GeomAPI_Shape> aFrom(new GeomAPI_Shape());
152         aFrom->setImpl(new TopoDS_Shape(aBuilder->Modified(aFromShape).First()));
153         myFromFaces.push_back(aFrom);
154         std::shared_ptr<GeomAPI_Shape> aTo(new GeomAPI_Shape());
155         aTo->setImpl(new TopoDS_Shape(aBuilder->Modified(aToShape).First()));
156         myToFaces.push_back(aTo);
157         myMkShape.reset(new GeomAlgoAPI_MakeShape(aBuilder));
158       }
159     }
160   } else {
161     // Getting bounding box for base shape.
162     const TopoDS_Shape& aBasisShape = theBasis->impl<TopoDS_Shape>();
163     Bnd_Box aBndBox;
164     BRepBndLib::Add(aBasisShape, aBndBox);
165     Standard_Real aXArr[2] = {aBndBox.CornerMin().X(), aBndBox.CornerMax().X()};
166     Standard_Real aYArr[2] = {aBndBox.CornerMin().Y(), aBndBox.CornerMax().Y()};
167     Standard_Real aZArr[2] = {aBndBox.CornerMin().Z(), aBndBox.CornerMax().Z()};
168     gp_Pnt aPoints[8];
169     int aNum = 0;
170     for(int i = 0; i < 2; i++) {
171       for(int j = 0; j < 2; j++) {
172         for(int k = 0; k < 2; k++) {
173           aPoints[aNum] = gp_Pnt(aXArr[i], aYArr[j], aZArr[k]);
174           aNum++;
175         }
176       }
177     }
178
179     // Project points to bounding planes. Search max distance to them.
180     IntAna_Quadric aBndToQuadric(gp_Pln(aToPnt->impl<gp_Pnt>(), aToDir->impl<gp_Dir>()));
181     IntAna_Quadric aBndFromQuadric(gp_Pln(aFromPnt->impl<gp_Pnt>(), aFromDir->impl<gp_Dir>()));
182     Standard_Real aMaxToDist = 0, aMaxFromDist = 0;
183     gp_Vec aNormal(aBaseDir->impl<gp_Dir>());
184     for(int i = 0; i < 8; i++) {
185       gp_Lin aLine(aPoints[i], aNormal);
186       IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
187       IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
188       if(aToIntAna.NbPoints() == 0 || aFromIntAna.NbPoints() == 0) {
189         return;
190       }
191       const gp_Pnt& aPntOnToFace = aToIntAna.Point(1);
192       const gp_Pnt& aPntOnFromFace = aFromIntAna.Point(1);
193       if(aPoints[i].Distance(aPntOnToFace) > aMaxToDist) {
194         aMaxToDist = aPoints[i].Distance(aPntOnToFace);
195       }
196       if(aPoints[i].Distance(aPntOnFromFace) > aMaxFromDist) {
197         aMaxFromDist = aPoints[i].Distance(aPntOnFromFace);
198       }
199     }
200     // We added 1 just to be sure that pipe is long enough for boolean operation.
201     Standard_Real aPipeLength = aMaxToDist + aMaxFromDist + 1;
202
203     // Making wire for pipe.
204     std::shared_ptr<GeomAPI_Pnt> aCentreOfMass = GeomAlgoAPI_ShapeTools::centreOfMass(theBasis);
205     const gp_Pnt aCentrePnt = aCentreOfMass->impl<gp_Pnt>();
206     TopoDS_Face aFace = TopoDS::Face(aBaseFace->impl<TopoDS_Shape>());
207     gp_Pnt aPipeStartPnt = aCentrePnt.Translated(aNormal.Scaled(aPipeLength));
208     gp_Pnt aPipeEndPnt = aCentrePnt.Translated(aNormal.Scaled(-aPipeLength));
209     TopoDS_Edge aPipeEdge = BRepBuilderAPI_MakeEdge(aPipeStartPnt, aPipeEndPnt);
210     TopoDS_Wire aPipeWire = BRepBuilderAPI_MakeWire(aPipeEdge).Wire();
211
212     // Making pipe.
213     ListOfMakeShape aListOfMakeShape;
214     BRepOffsetAPI_MakePipe* aPipeBuilder = new BRepOffsetAPI_MakePipe(aPipeWire, aBasisShape);
215     if(!aPipeBuilder) {
216       return;
217     }
218     std::shared_ptr<GeomAPI_Shape> aWire(new GeomAPI_Shape);
219     std::shared_ptr<GeomAPI_Shape> aBShape(new GeomAPI_Shape);
220     aWire->setImpl(new TopoDS_Shape(aPipeWire));
221     aBShape->setImpl(new TopoDS_Shape(aBasisShape));
222     aListOfMakeShape.push_back(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aPipeBuilder, aWire, aBShape)));
223     TopoDS_Shape aResult = aPipeBuilder->Shape();
224
225     // Orienting bounding planes.
226     gp_Lin aLine(aCentrePnt, aNormal);
227     IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
228     IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
229     Standard_Real aToParameter = aToIntAna.ParamOnConic(1);
230     Standard_Real aFromParameter = aFromIntAna.ParamOnConic(1);
231     if(aToParameter > aFromParameter) {
232       gp_Vec aVec = aToDir->impl<gp_Dir>();
233       if((aVec * aNormal) > 0) {
234         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
235         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
236       }
237       aVec = aFromDir->impl<gp_Dir>();
238       if((aVec * aNormal) < 0) {
239         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
240         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
241       }
242     } else {
243       gp_Vec aVec = aToDir->impl<gp_Dir>();
244       if((aVec * aNormal) < 0) {
245         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
246         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
247       }
248       aVec = aFromDir->impl<gp_Dir>();
249       if((aVec * aNormal) > 0) {
250         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
251         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
252       }
253     }
254
255     // Making solids from bounding planes.
256     TopoDS_Shell aToShell, aFromShell;
257     TopoDS_Solid aToSolid, aFromSolid;
258     const TopoDS_Shape& aToShape   = aBoundingToShape->impl<TopoDS_Shape>();
259     const TopoDS_Shape& aFromShape = aBoundingFromShape->impl<TopoDS_Shape>();
260     BRep_Builder aBoundingBuilder;
261     aBoundingBuilder.MakeShell(aToShell);
262     aBoundingBuilder.MakeShell(aFromShell);
263     aBoundingBuilder.Add(aToShell, aToShape);
264     aBoundingBuilder.Add(aFromShell, aFromShape);
265     aBoundingBuilder.MakeSolid(aToSolid);
266     aBoundingBuilder.MakeSolid(aFromSolid);
267     aBoundingBuilder.Add(aToSolid, aToShell);
268     aBoundingBuilder.Add(aFromSolid, aFromShell);
269
270     // Cutting with to plane.
271     BRepAlgoAPI_Cut* aToCutBuilder = new BRepAlgoAPI_Cut(aResult, aToSolid);
272     aToCutBuilder->Build();
273     if(!aToCutBuilder->IsDone()) {
274       return;
275     }
276     aListOfMakeShape.push_back(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aToCutBuilder)));
277     const TopTools_ListOfShape& aToShapes = aToCutBuilder->Modified(aToShape);
278     for(TopTools_ListIteratorOfListOfShape anIt(aToShapes); anIt.More(); anIt.Next()) {
279       std::shared_ptr<GeomAPI_Shape> aShape(new GeomAPI_Shape());
280       aShape->setImpl(new TopoDS_Shape(anIt.Value()));
281       myToFaces.push_back(aShape);
282     }
283     aResult = aToCutBuilder->Shape();
284
285     // Cutting with from plane.
286     BRepAlgoAPI_Cut* aFromCutBuilder = new BRepAlgoAPI_Cut(aResult, aFromSolid);
287     aFromCutBuilder->Build();
288     if(!aFromCutBuilder->IsDone()) {
289       return;
290     }
291     aListOfMakeShape.push_back(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aFromCutBuilder)));
292     const TopTools_ListOfShape& aFromShapes = aFromCutBuilder->Modified(aFromShape);
293     for(TopTools_ListIteratorOfListOfShape anIt(aFromShapes); anIt.More(); anIt.Next()) {
294       std::shared_ptr<GeomAPI_Shape> aShape(new GeomAPI_Shape());
295       aShape->setImpl(new TopoDS_Shape(anIt.Value()));
296       myFromFaces.push_back(aShape);
297     }
298     aResult = aFromCutBuilder->Shape();
299
300     TopExp_Explorer anExp(aResult, TopAbs_SOLID);
301     if(!anExp.More()) {
302       return;
303     }
304     if(aResult.ShapeType() == TopAbs_COMPOUND) {
305       aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
306     }
307     if(aResult.ShapeType() == TopAbs_COMPOUND) {
308       std::shared_ptr<GeomAPI_Shape> aCompound(new GeomAPI_Shape);
309       aCompound->setImpl(new TopoDS_Shape(aResult));
310       ListOfShape aCompSolids, aFreeSolids;
311       GeomAlgoAPI_ShapeTools::combineShapes(aCompound, GeomAPI_Shape::COMPSOLID, aCompSolids, aFreeSolids);
312       if(aCompSolids.size() == 1 && aFreeSolids.size() == 0) {
313         aResult = aCompSolids.front()->impl<TopoDS_Shape>();
314       } else if (aCompSolids.size() > 1 || (aCompSolids.size() >= 1 && aFreeSolids.size() >= 1)) {
315         TopoDS_Compound aResultComp;
316         TopoDS_Builder aBuilder;
317         aBuilder.MakeCompound(aResultComp);
318         for(ListOfShape::const_iterator anIter = aCompSolids.cbegin(); anIter != aCompSolids.cend(); anIter++) {
319           aBuilder.Add(aResultComp, (*anIter)->impl<TopoDS_Shape>());
320         }
321         for(ListOfShape::const_iterator anIter = aFreeSolids.cbegin(); anIter != aFreeSolids.cend(); anIter++) {
322           aBuilder.Add(aResultComp, (*anIter)->impl<TopoDS_Shape>());
323         }
324         aResult = aResultComp;
325       }
326     }
327
328     // Fill data map to keep correct orientation of sub-shapes.
329     myMap = std::shared_ptr<GeomAPI_DataMapOfShapeShape>(new GeomAPI_DataMapOfShapeShape);
330     for (TopExp_Explorer Exp(aResult,TopAbs_FACE); Exp.More(); Exp.Next()) {
331       std::shared_ptr<GeomAPI_Shape> aCurrentShape(new GeomAPI_Shape());
332       aCurrentShape->setImpl(new TopoDS_Shape(Exp.Current()));
333       myMap->bind(aCurrentShape, aCurrentShape);
334     }
335     myShape = std::shared_ptr<GeomAPI_Shape>(new GeomAPI_Shape);
336     myShape->setImpl(new TopoDS_Shape(aResult));
337     myMkShape = std::shared_ptr<GeomAlgoAPI_MakeShapeList>(new GeomAlgoAPI_MakeShapeList(aListOfMakeShape));
338     myDone = true;
339   }
340 }
341
342 //=================================================================================================
343 bool GeomAlgoAPI_Prism::isDone() const
344 {
345   return myDone;
346 }
347
348 //=================================================================================================
349 bool GeomAlgoAPI_Prism::isValid() const
350 {
351   BRepCheck_Analyzer aChecker(myShape->impl<TopoDS_Shape>());
352   return (aChecker.IsValid() == Standard_True);
353 }
354
355 //=================================================================================================
356 bool GeomAlgoAPI_Prism::hasVolume() const
357 {
358   bool hasVolume(false);
359   if(isValid()) {
360     const TopoDS_Shape& aRShape = myShape->impl<TopoDS_Shape>();
361     GProp_GProps aGProp;
362     BRepGProp::VolumeProperties(aRShape, aGProp);
363     if(aGProp.Mass() > Precision::Confusion())
364       hasVolume = true;
365   }
366   return hasVolume;
367 }
368
369 //=================================================================================================
370 std::shared_ptr<GeomAPI_Shape> GeomAlgoAPI_Prism::shape() const
371 {
372   return myShape;
373 }
374
375 //=================================================================================================
376 const ListOfShape& GeomAlgoAPI_Prism::fromFaces() const
377 {
378   return myFromFaces;
379 }
380
381 //=================================================================================================
382 const ListOfShape& GeomAlgoAPI_Prism::toFaces() const
383 {
384   return myToFaces;
385 }
386
387 //=================================================================================================
388 std::shared_ptr<GeomAPI_DataMapOfShapeShape> GeomAlgoAPI_Prism::mapOfShapes() const
389 {
390   return myMap;
391 }
392
393 //=================================================================================================
394 std::shared_ptr<GeomAlgoAPI_MakeShape> GeomAlgoAPI_Prism::makeShape() const
395 {
396   return myMkShape;
397 }