]> SALOME platform Git repositories - modules/shaper.git/blob - src/GeomAlgoAPI/GeomAlgoAPI_Prism.cpp
Salome HOME
Issue #1343 Extrusion fixes
[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_ShapeTools.h>
17
18 #include <Bnd_Box.hxx>
19 #include <BRep_Builder.hxx>
20 #include <BRepAlgoAPI_Cut.hxx>
21 #include <BRepBndLib.hxx>
22 #include <BRepBuilderAPI_FindPlane.hxx>
23 #include <BRepBuilderAPI_Transform.hxx>
24 #include <BRepPrimAPI_MakePrism.hxx>
25 #include <Geom_Plane.hxx>
26 #include <gp_Pln.hxx>
27 #include <IntAna_IntConicQuad.hxx>
28 #include <IntAna_Quadric.hxx>
29 #include <TopExp_Explorer.hxx>
30 #include <TopoDS_Shell.hxx>
31 #include <TopoDS_Solid.hxx>
32 #include <TopTools_ListIteratorOfListOfShape.hxx>
33
34 #include <BRepTools.hxx>
35
36 //=================================================================================================
37 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr theBaseShape,
38                                      const double       theToSize,
39                                      const double       theFromSize)
40 {
41   build(theBaseShape, std::shared_ptr<GeomAPI_Dir>(), GeomShapePtr(), theToSize, GeomShapePtr(), theFromSize);
42 }
43
44 //=================================================================================================
45 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr                 theBaseShape,
46                                      const std::shared_ptr<GeomAPI_Dir> theDirection,
47                                      const double                       theToSize,
48                                      const double                       theFromSize)
49 {
50   build(theBaseShape, theDirection, GeomShapePtr(), theToSize, GeomShapePtr(), theFromSize);
51 }
52
53 //=================================================================================================
54 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr theBaseShape,
55                                      const GeomShapePtr theToShape,
56                                      const double       theToSize,
57                                      const GeomShapePtr theFromShape,
58                                      const double       theFromSize)
59 {
60   build(theBaseShape, std::shared_ptr<GeomAPI_Dir>(), theToShape, theToSize, theFromShape, theFromSize);
61 }
62
63 //=================================================================================================
64 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr                 theBaseShape,
65                                      const std::shared_ptr<GeomAPI_Dir> theDirection,
66                                      const GeomShapePtr                 theToShape,
67                                      const double                       theToSize,
68                                      const GeomShapePtr                 theFromShape,
69                                      const double                       theFromSize)
70 {
71   build(theBaseShape, theDirection, theToShape, theToSize, theFromShape, theFromSize);
72 }
73
74 //=================================================================================================
75 void GeomAlgoAPI_Prism::build(const GeomShapePtr&                theBaseShape,
76                               const std::shared_ptr<GeomAPI_Dir> theDirection,
77                               const GeomShapePtr&                theToShape,
78                               const double                       theToSize,
79                               const GeomShapePtr&                theFromShape,
80                               const double                       theFromSize)
81 {
82   if(!theBaseShape.get() ||
83     (((!theFromShape.get() && !theToShape.get()) || (theFromShape.get() && theToShape.get() && theFromShape->isEqual(theToShape)))
84     && (theFromSize == -theToSize))) {
85     return;
86   }
87
88   // Getting base shape.
89   const TopoDS_Shape& aBaseShape = theBaseShape->impl<TopoDS_Shape>();
90   TopAbs_ShapeEnum aShapeTypeToExp;
91   switch(aBaseShape.ShapeType()) {
92     case TopAbs_VERTEX:
93       aShapeTypeToExp = TopAbs_VERTEX;
94       break;
95     case TopAbs_EDGE:
96     case TopAbs_WIRE:
97       aShapeTypeToExp = TopAbs_EDGE;
98       break;
99     case TopAbs_FACE:
100     case TopAbs_SHELL:
101       aShapeTypeToExp = TopAbs_FACE;
102       break;
103   }
104
105   // Getting direction.
106   gp_Vec aDirVec;
107   std::shared_ptr<GeomAPI_Pnt> aBaseLoc;
108   std::shared_ptr<GeomAPI_Dir> aBaseDir;
109   GeomShapePtr aBasePlane;
110   const bool isBoundingShapesSet = theFromShape.get() || theToShape.get();
111   BRepBuilderAPI_FindPlane aFindPlane(aBaseShape);
112   if(aBaseShape.ShapeType() == TopAbs_VERTEX || aBaseShape.ShapeType() == TopAbs_EDGE ||
113      aFindPlane.Found() != Standard_True) {
114     // Direction and both bounding planes should be set or empty.
115     if(!theDirection.get() || (isBoundingShapesSet && (!theToShape.get() || !theFromShape.get()))) {
116       return;
117     }
118
119     aBaseDir = theDirection;
120     aDirVec = theDirection->impl<gp_Dir>();
121   } else {
122     Handle(Geom_Plane) aPlane = aFindPlane.Plane();
123     gp_Pnt aLoc = aPlane->Axis().Location();
124     aDirVec = aPlane->Axis().Direction();
125
126     aBaseLoc.reset(new GeomAPI_Pnt(aLoc.X(), aLoc.Y(), aLoc.Z()));
127     aBaseDir.reset(new GeomAPI_Dir(aDirVec.X(), aDirVec.Y(), aDirVec.Z()));
128     aBasePlane = GeomAlgoAPI_FaceBuilder::planarFace(aBaseLoc, aBaseDir);
129   }
130
131   TopoDS_Shape aResult;
132   if(!isBoundingShapesSet) {
133     // Moving base shape.
134     gp_Trsf aTrsf;
135     aTrsf.SetTranslation(aDirVec * -theFromSize);
136     BRepBuilderAPI_Transform* aTransformBuilder = new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
137     if(!aTransformBuilder) {
138       return;
139     }
140     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aTransformBuilder)));
141     if(!aTransformBuilder->IsDone()) {
142       return;
143     }
144     TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
145
146     // Making prism.
147     BRepPrimAPI_MakePrism* aPrismBuilder = new BRepPrimAPI_MakePrism(aMovedBase, aDirVec * (theFromSize + theToSize));
148     if(!aPrismBuilder) {
149       return;
150     }
151     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aPrismBuilder)));
152     if(!aPrismBuilder->IsDone()) {
153       return;
154     }
155     aResult = aPrismBuilder->Shape();
156
157     // Setting naming.
158     for(TopExp_Explorer anExp(aMovedBase, aShapeTypeToExp); anExp.More(); anExp.Next()) {
159       const TopoDS_Shape& aFace = anExp.Current();
160       GeomShapePtr aFromShape(new GeomAPI_Shape), aToShape(new GeomAPI_Shape);
161       aFromShape->setImpl(new TopoDS_Shape(aPrismBuilder->FirstShape(aFace)));
162       aToShape->setImpl(new TopoDS_Shape(aPrismBuilder->LastShape(aFace)));
163       this->addFromShape(aFromShape);
164       this->addToShape(aToShape);
165     }
166   } else {
167     GeomShapePtr aBoundingFromShape = theFromShape ? theFromShape : aBasePlane;
168     GeomShapePtr aBoundingToShape   = theToShape   ? theToShape   : aBasePlane;
169
170     // Moving prism bounding faces according to "from" and "to" sizes.
171     std::shared_ptr<GeomAPI_Face> aFromFace(new GeomAPI_Face(aBoundingFromShape));
172     std::shared_ptr<GeomAPI_Pln>  aFromPln = aFromFace->getPlane();
173     std::shared_ptr<GeomAPI_Pnt>  aFromLoc = aFromPln->location();
174     std::shared_ptr<GeomAPI_Dir>  aFromDir = aFromPln->direction();
175
176     std::shared_ptr<GeomAPI_Face> aToFace(new GeomAPI_Face(aBoundingToShape));
177     std::shared_ptr<GeomAPI_Pln>  aToPln = aToFace->getPlane();
178     std::shared_ptr<GeomAPI_Pnt>  aToLoc = aToPln->location();
179     std::shared_ptr<GeomAPI_Dir>  aToDir = aToPln->direction();
180
181     bool aSign = aFromLoc->xyz()->dot(aBaseDir->xyz()) > aToLoc->xyz()->dot(aBaseDir->xyz());
182
183     std::shared_ptr<GeomAPI_Pnt> aFromPnt(new GeomAPI_Pnt(aFromLoc->xyz()->added(aBaseDir->xyz()->multiplied(
184                                                           aSign ? theFromSize : -theFromSize))));
185     aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
186
187     std::shared_ptr<GeomAPI_Pnt> aToPnt(new GeomAPI_Pnt(aToLoc->xyz()->added(aBaseDir->xyz()->multiplied(
188                                                         aSign ? -theToSize : theToSize))));
189     aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
190
191     // Getting bounding box for base shape.
192     Bnd_Box aBndBox;
193     BRepBndLib::Add(aBaseShape, aBndBox);
194     Standard_Real aXArr[2] = {aBndBox.CornerMin().X(), aBndBox.CornerMax().X()};
195     Standard_Real aYArr[2] = {aBndBox.CornerMin().Y(), aBndBox.CornerMax().Y()};
196     Standard_Real aZArr[2] = {aBndBox.CornerMin().Z(), aBndBox.CornerMax().Z()};
197     gp_Pnt aPoints[8];
198     int aNum = 0;
199     for(int i = 0; i < 2; i++) {
200       for(int j = 0; j < 2; j++) {
201         for(int k = 0; k < 2; k++) {
202           aPoints[aNum] = gp_Pnt(aXArr[i], aYArr[j], aZArr[k]);
203           aNum++;
204         }
205       }
206     }
207
208     // Project points to bounding planes. Search max distance to them.
209     IntAna_Quadric aBndToQuadric(gp_Pln(aToPnt->impl<gp_Pnt>(), aToDir->impl<gp_Dir>()));
210     IntAna_Quadric aBndFromQuadric(gp_Pln(aFromPnt->impl<gp_Pnt>(), aFromDir->impl<gp_Dir>()));
211     Standard_Real aMaxToDist = 0, aMaxFromDist = 0;
212     for(int i = 0; i < 8; i++) {
213       gp_Lin aLine(aPoints[i], aDirVec);
214       IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
215       IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
216       if(aToIntAna.NbPoints() == 0 || aFromIntAna.NbPoints() == 0) {
217         return;
218       }
219       const gp_Pnt& aPntOnToFace = aToIntAna.Point(1);
220       const gp_Pnt& aPntOnFromFace = aFromIntAna.Point(1);
221       if(aPoints[i].Distance(aPntOnToFace) > aMaxToDist) {
222         aMaxToDist = aPoints[i].Distance(aPntOnToFace);
223       }
224       if(aPoints[i].Distance(aPntOnFromFace) > aMaxFromDist) {
225         aMaxFromDist = aPoints[i].Distance(aPntOnFromFace);
226       }
227     }
228
229     // We added 1 just to be sure that prism is long enough for boolean operation.
230     double aPrismLength = aMaxToDist + aMaxFromDist + 1;
231
232     // Moving base shape.
233     gp_Trsf aTrsf;
234     aTrsf.SetTranslation(aDirVec * -aPrismLength);
235     BRepBuilderAPI_Transform* aTransformBuilder = new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
236     if(!aTransformBuilder) {
237       return;
238     }
239     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aTransformBuilder)));
240     if(!aTransformBuilder->IsDone()) {
241       return;
242     }
243     TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
244
245     // Making prism.
246     BRepPrimAPI_MakePrism* aPrismBuilder = new BRepPrimAPI_MakePrism(aMovedBase, aDirVec * 2 * aPrismLength);
247     if(!aPrismBuilder) {
248       return;
249     }
250     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aPrismBuilder)));
251     if(!aPrismBuilder->IsDone()) {
252       return;
253     }
254     aResult = aPrismBuilder->Shape();
255
256     // Orienting bounding planes.
257     std::shared_ptr<GeomAPI_Pnt> aCentreOfMass = GeomAlgoAPI_ShapeTools::centreOfMass(theBaseShape);
258     const gp_Pnt& aCentrePnt = aCentreOfMass->impl<gp_Pnt>();
259     gp_Lin aLine(aCentrePnt, aDirVec);
260     IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
261     IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
262     Standard_Real aToParameter = aToIntAna.ParamOnConic(1);
263     Standard_Real aFromParameter = aFromIntAna.ParamOnConic(1);
264     if(aToParameter > aFromParameter) {
265       gp_Vec aVec = aToDir->impl<gp_Dir>();
266       if((aVec * aDirVec) > 0) {
267         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
268         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
269       }
270       aVec = aFromDir->impl<gp_Dir>();
271       if((aVec * aDirVec) < 0) {
272         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
273         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
274       }
275     } else {
276       gp_Vec aVec = aToDir->impl<gp_Dir>();
277       if((aVec * aDirVec) < 0) {
278         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
279         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
280       }
281       aVec = aFromDir->impl<gp_Dir>();
282       if((aVec * aDirVec) > 0) {
283         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
284         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
285       }
286     }
287
288     // Making solids from bounding planes.
289     TopoDS_Shell aToShell, aFromShell;
290     TopoDS_Solid aToSolid, aFromSolid;
291     const TopoDS_Shape& aToShape   = aBoundingToShape->impl<TopoDS_Shape>();
292     const TopoDS_Shape& aFromShape = aBoundingFromShape->impl<TopoDS_Shape>();
293     BRep_Builder aBoundingBuilder;
294     aBoundingBuilder.MakeShell(aToShell);
295     aBoundingBuilder.Add(aToShell, aToShape);
296     aBoundingBuilder.MakeShell(aFromShell);
297     aBoundingBuilder.Add(aFromShell, aFromShape);
298     aBoundingBuilder.MakeSolid(aToSolid);
299     aBoundingBuilder.Add(aToSolid, aToShell);
300     aBoundingBuilder.MakeSolid(aFromSolid);
301     aBoundingBuilder.Add(aFromSolid, aFromShell);
302
303     // Cutting with to plane.
304     BRepAlgoAPI_Cut* aToCutBuilder = new BRepAlgoAPI_Cut(aResult, aToSolid);
305     aToCutBuilder->Build();
306     if(!aToCutBuilder->IsDone()) {
307       return;
308     }
309     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aToCutBuilder)));
310     const TopTools_ListOfShape& aToShapes = aToCutBuilder->Modified(aToShape);
311     for(TopTools_ListIteratorOfListOfShape anIt(aToShapes); anIt.More(); anIt.Next()) {
312       GeomShapePtr aShape(new GeomAPI_Shape());
313       aShape->setImpl(new TopoDS_Shape(anIt.Value()));
314       this->addToShape(aShape);
315     }
316     aResult = aToCutBuilder->Shape();
317     if(aResult.ShapeType() == TopAbs_COMPOUND) {
318       aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
319     }
320
321     // Cutting with from plane.
322     BRepAlgoAPI_Cut* aFromCutBuilder = new BRepAlgoAPI_Cut(aResult, aFromSolid);
323     aFromCutBuilder->Build();
324     if(!aFromCutBuilder->IsDone()) {
325       return;
326     }
327     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aFromCutBuilder)));
328     const TopTools_ListOfShape& aFromShapes = aFromCutBuilder->Modified(aFromShape);
329     for(TopTools_ListIteratorOfListOfShape anIt(aFromShapes); anIt.More(); anIt.Next()) {
330       GeomShapePtr aShape(new GeomAPI_Shape());
331       aShape->setImpl(new TopoDS_Shape(anIt.Value()));
332       this->addFromShape(aShape);
333     }
334     aResult = aFromCutBuilder->Shape();
335
336     TopoDS_Iterator anIt(aResult);
337     if(!anIt.More()) {
338       return;
339     }
340     if(aResult.ShapeType() == TopAbs_COMPOUND) {
341       aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
342     }
343     if(aResult.ShapeType() == TopAbs_COMPOUND) {
344       GeomShapePtr aCompound(new GeomAPI_Shape);
345       aCompound->setImpl(new TopoDS_Shape(aResult));
346       ListOfShape aCompSolids, aFreeSolids;
347       GeomAlgoAPI_ShapeTools::combineShapes(aCompound, GeomAPI_Shape::COMPSOLID, aCompSolids, aFreeSolids);
348       if(aCompSolids.size() == 1 && aFreeSolids.size() == 0) {
349         aResult = aCompSolids.front()->impl<TopoDS_Shape>();
350       } else if (aCompSolids.size() > 1 || (aCompSolids.size() >= 1 && aFreeSolids.size() >= 1)) {
351         TopoDS_Compound aResultComp;
352         TopoDS_Builder aBuilder;
353         aBuilder.MakeCompound(aResultComp);
354         for(ListOfShape::const_iterator anIter = aCompSolids.cbegin(); anIter != aCompSolids.cend(); anIter++) {
355           aBuilder.Add(aResultComp, (*anIter)->impl<TopoDS_Shape>());
356         }
357         for(ListOfShape::const_iterator anIter = aFreeSolids.cbegin(); anIter != aFreeSolids.cend(); anIter++) {
358           aBuilder.Add(aResultComp, (*anIter)->impl<TopoDS_Shape>());
359         }
360         aResult = aResultComp;
361       }
362     }
363   }
364
365   // Setting result.
366   if(aResult.IsNull()) {
367     return;
368   }
369   GeomShapePtr aShape(new GeomAPI_Shape());
370   aShape->setImpl(new TopoDS_Shape(aResult));
371   this->setShape(aShape);
372   this->setDone(true);
373 }