Salome HOME
Issue #1367: Fill feature.
[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 <BRepTools.hxx>
25 #include <Geom_Curve.hxx>
26 #include <Geom2d_Curve.hxx>
27 #include <BRepLib_CheckCurveOnSurface.hxx>
28 #include <BRepPrimAPI_MakePrism.hxx>
29 #include <Geom_Plane.hxx>
30 #include <Geom_RectangularTrimmedSurface.hxx>
31 #include <gp_Pln.hxx>
32 #include <IntAna_IntConicQuad.hxx>
33 #include <IntAna_Quadric.hxx>
34 #include <IntTools_Context.hxx>
35 #include <TopExp_Explorer.hxx>
36 #include <TopoDS.hxx>
37 #include <TopoDS_Edge.hxx>
38 #include <TopoDS_Shell.hxx>
39 #include <TopoDS_Solid.hxx>
40 #include <TopTools_ListIteratorOfListOfShape.hxx>
41
42 //=================================================================================================
43 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr theBaseShape,
44                                      const double       theToSize,
45                                      const double       theFromSize)
46 {
47   build(theBaseShape, std::shared_ptr<GeomAPI_Dir>(), GeomShapePtr(), theToSize, GeomShapePtr(), theFromSize);
48 }
49
50 //=================================================================================================
51 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr                 theBaseShape,
52                                      const std::shared_ptr<GeomAPI_Dir> theDirection,
53                                      const double                       theToSize,
54                                      const double                       theFromSize)
55 {
56   build(theBaseShape, theDirection, GeomShapePtr(), theToSize, GeomShapePtr(), theFromSize);
57 }
58
59 //=================================================================================================
60 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr theBaseShape,
61                                      const GeomShapePtr theToShape,
62                                      const double       theToSize,
63                                      const GeomShapePtr theFromShape,
64                                      const double       theFromSize)
65 {
66   build(theBaseShape, std::shared_ptr<GeomAPI_Dir>(), theToShape, theToSize, theFromShape, theFromSize);
67 }
68
69 //=================================================================================================
70 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr                 theBaseShape,
71                                      const std::shared_ptr<GeomAPI_Dir> theDirection,
72                                      const GeomShapePtr                 theToShape,
73                                      const double                       theToSize,
74                                      const GeomShapePtr                 theFromShape,
75                                      const double                       theFromSize)
76 {
77   build(theBaseShape, theDirection, theToShape, theToSize, theFromShape, theFromSize);
78 }
79
80 //=================================================================================================
81 void GeomAlgoAPI_Prism::build(const GeomShapePtr&                theBaseShape,
82                               const std::shared_ptr<GeomAPI_Dir> theDirection,
83                               const GeomShapePtr&                theToShape,
84                               const double                       theToSize,
85                               const GeomShapePtr&                theFromShape,
86                               const double                       theFromSize)
87 {
88   if(!theBaseShape.get() ||
89     (((!theFromShape.get() && !theToShape.get()) || (theFromShape.get() && theToShape.get() && theFromShape->isEqual(theToShape)))
90     && (theFromSize == -theToSize))) {
91     return;
92   }
93
94   // Getting base shape.
95   const TopoDS_Shape& aBaseShape = theBaseShape->impl<TopoDS_Shape>();
96   TopAbs_ShapeEnum aShapeTypeToExp;
97   switch(aBaseShape.ShapeType()) {
98     case TopAbs_VERTEX:
99       aShapeTypeToExp = TopAbs_VERTEX;
100       break;
101     case TopAbs_EDGE:
102     case TopAbs_WIRE:
103       aShapeTypeToExp = TopAbs_EDGE;
104       break;
105     case TopAbs_FACE:
106     case TopAbs_SHELL:
107       aShapeTypeToExp = TopAbs_FACE;
108       break;
109     default:
110       return;
111   }
112
113   // Getting direction.
114   gp_Vec aDirVec;
115   std::shared_ptr<GeomAPI_Pnt> aBaseLoc;
116   std::shared_ptr<GeomAPI_Dir> aBaseDir;
117   GeomShapePtr aBasePlane;
118   const bool isBoundingShapesSet = theFromShape.get() || theToShape.get();
119   BRepBuilderAPI_FindPlane aFindPlane(aBaseShape);
120   if(theDirection.get()) {
121     aBaseDir = theDirection;
122     aDirVec = theDirection->impl<gp_Dir>();
123   } else {
124     if(aBaseShape.ShapeType() == TopAbs_VERTEX
125         || aBaseShape.ShapeType() == TopAbs_EDGE
126         || aFindPlane.Found() == Standard_False) {
127       return;
128     }
129
130     Handle(Geom_Plane) aPlane;
131     if(aBaseShape.ShapeType() == TopAbs_FACE || aBaseShape.ShapeType() == TopAbs_SHELL) {
132       TopExp_Explorer anExp(aBaseShape, TopAbs_FACE);
133       const TopoDS_Shape& aFace = anExp.Current();
134       Handle(Geom_Surface) aSurface = BRep_Tool::Surface(TopoDS::Face(aFace));
135       if(aSurface->DynamicType() == STANDARD_TYPE(Geom_RectangularTrimmedSurface)) {
136         Handle(Geom_RectangularTrimmedSurface) aTrimSurface =
137           Handle(Geom_RectangularTrimmedSurface)::DownCast(aSurface);
138         aSurface = aTrimSurface->BasisSurface();
139       }
140       if(aSurface->DynamicType() != STANDARD_TYPE(Geom_Plane)) {
141         return;
142       }
143       aPlane = Handle(Geom_Plane)::DownCast(aSurface);
144     } else {
145       aPlane = aFindPlane.Plane();
146     }
147     gp_Pnt aLoc = aPlane->Axis().Location();
148     aDirVec = aPlane->Axis().Direction();
149     aBaseLoc.reset(new GeomAPI_Pnt(aLoc.X(), aLoc.Y(), aLoc.Z()));
150     aBaseDir.reset(new GeomAPI_Dir(aDirVec.X(), aDirVec.Y(), aDirVec.Z()));
151   }
152   if(!aBaseLoc.get()) {
153     gp_Pnt aLoc;
154     gp_XYZ aDirXYZ = aDirVec.XYZ();
155     Standard_Real aMinParam = Precision::Infinite();
156     for(TopExp_Explorer anExp(aBaseShape, TopAbs_VERTEX); anExp.More(); anExp.Next()) {
157       const TopoDS_Shape& aVertex = anExp.Current();
158       gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(aVertex));
159       double aParam = aDirXYZ.Dot(aPnt.XYZ());
160       if(aParam < aMinParam) {
161         aMinParam = aParam;
162         aLoc = aPnt;
163       }
164     }
165     aBaseLoc.reset(new GeomAPI_Pnt(aLoc.X(), aLoc.Y(), aLoc.Z()));
166   }
167   aBasePlane = GeomAlgoAPI_FaceBuilder::planarFace(aBaseLoc, aBaseDir);
168
169   TopoDS_Shape aResult;
170   if(!isBoundingShapesSet) {
171     // Moving base shape.
172     gp_Trsf aTrsf;
173     aTrsf.SetTranslation(aDirVec * -theFromSize);
174     BRepBuilderAPI_Transform* aTransformBuilder = new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
175     if(!aTransformBuilder) {
176       return;
177     }
178     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aTransformBuilder)));
179     if(!aTransformBuilder->IsDone()) {
180       return;
181     }
182     TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
183
184     // Making prism.
185     BRepPrimAPI_MakePrism* aPrismBuilder = new BRepPrimAPI_MakePrism(aMovedBase, aDirVec * (theFromSize + theToSize));
186     if(!aPrismBuilder) {
187       return;
188     }
189     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aPrismBuilder)));
190     if(!aPrismBuilder->IsDone()) {
191       return;
192     }
193     aResult = aPrismBuilder->Shape();
194
195     // Setting naming.
196     for(TopExp_Explorer anExp(aMovedBase, aShapeTypeToExp); anExp.More(); anExp.Next()) {
197       const TopoDS_Shape& aShape = anExp.Current();
198       GeomShapePtr aFromShape(new GeomAPI_Shape), aToShape(new GeomAPI_Shape);
199       aFromShape->setImpl(new TopoDS_Shape(aPrismBuilder->FirstShape(aShape)));
200       aToShape->setImpl(new TopoDS_Shape(aPrismBuilder->LastShape(aShape)));
201       this->addFromShape(aFromShape);
202       this->addToShape(aToShape);
203     }
204   } else {
205     GeomShapePtr aBoundingFromShape = theFromShape ? theFromShape : aBasePlane;
206     GeomShapePtr aBoundingToShape   = theToShape   ? theToShape   : aBasePlane;
207
208     // Moving prism bounding faces according to "from" and "to" sizes.
209     std::shared_ptr<GeomAPI_Pln> aFromPln = GeomAPI_Face(aBoundingFromShape).getPlane();
210     std::shared_ptr<GeomAPI_Pnt> aFromLoc = aFromPln->location();
211     std::shared_ptr<GeomAPI_Dir> aFromDir = aFromPln->direction();
212
213     std::shared_ptr<GeomAPI_Pln> aToPln = GeomAPI_Face(aBoundingToShape).getPlane();
214     std::shared_ptr<GeomAPI_Pnt> aToLoc = aToPln->location();
215     std::shared_ptr<GeomAPI_Dir> aToDir = aToPln->direction();
216
217     bool aSign = aFromLoc->xyz()->dot(aBaseDir->xyz()) > aToLoc->xyz()->dot(aBaseDir->xyz());
218
219     std::shared_ptr<GeomAPI_Pnt> aFromPnt(new GeomAPI_Pnt(aFromLoc->xyz()->added(aBaseDir->xyz()->multiplied(
220                                                           aSign ? theFromSize : -theFromSize))));
221     aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
222
223     std::shared_ptr<GeomAPI_Pnt> aToPnt(new GeomAPI_Pnt(aToLoc->xyz()->added(aBaseDir->xyz()->multiplied(
224                                                         aSign ? -theToSize : theToSize))));
225     aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
226
227     // Getting bounding box for base shape.
228     Bnd_Box aBndBox;
229     BRepBndLib::Add(aBaseShape, aBndBox);
230     Standard_Real aXArr[2] = {aBndBox.CornerMin().X(), aBndBox.CornerMax().X()};
231     Standard_Real aYArr[2] = {aBndBox.CornerMin().Y(), aBndBox.CornerMax().Y()};
232     Standard_Real aZArr[2] = {aBndBox.CornerMin().Z(), aBndBox.CornerMax().Z()};
233     gp_Pnt aPoints[8];
234     int aNum = 0;
235     for(int i = 0; i < 2; i++) {
236       for(int j = 0; j < 2; j++) {
237         for(int k = 0; k < 2; k++) {
238           aPoints[aNum] = gp_Pnt(aXArr[i], aYArr[j], aZArr[k]);
239           aNum++;
240         }
241       }
242     }
243
244     // Project points to bounding planes. Search max distance to them.
245     IntAna_Quadric aBndToQuadric(gp_Pln(aToPnt->impl<gp_Pnt>(), aToDir->impl<gp_Dir>()));
246     IntAna_Quadric aBndFromQuadric(gp_Pln(aFromPnt->impl<gp_Pnt>(), aFromDir->impl<gp_Dir>()));
247     Standard_Real aMaxToDist = 0, aMaxFromDist = 0;
248     for(int i = 0; i < 8; i++) {
249       gp_Lin aLine(aPoints[i], aDirVec);
250       IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
251       IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
252       if(aToIntAna.NbPoints() == 0 || aFromIntAna.NbPoints() == 0) {
253         return;
254       }
255       const gp_Pnt& aPntOnToFace = aToIntAna.Point(1);
256       const gp_Pnt& aPntOnFromFace = aFromIntAna.Point(1);
257       if(aPoints[i].Distance(aPntOnToFace) > aMaxToDist) {
258         aMaxToDist = aPoints[i].Distance(aPntOnToFace);
259       }
260       if(aPoints[i].Distance(aPntOnFromFace) > aMaxFromDist) {
261         aMaxFromDist = aPoints[i].Distance(aPntOnFromFace);
262       }
263     }
264
265     // We added 1 just to be sure that prism is long enough for boolean operation.
266     double aPrismLength = aMaxToDist + aMaxFromDist + 1;
267
268     // Moving base shape.
269     gp_Trsf aTrsf;
270     aTrsf.SetTranslation(aDirVec * -aPrismLength);
271     BRepBuilderAPI_Transform* aTransformBuilder = new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
272     if(!aTransformBuilder) {
273       return;
274     }
275     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aTransformBuilder)));
276     if(!aTransformBuilder->IsDone()) {
277       return;
278     }
279     TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
280
281     // Making prism.
282     BRepPrimAPI_MakePrism* aPrismBuilder = new BRepPrimAPI_MakePrism(aMovedBase, aDirVec * 2 * aPrismLength);
283     if(!aPrismBuilder) {
284       return;
285     }
286     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aPrismBuilder)));
287     if(!aPrismBuilder->IsDone()) {
288       return;
289     }
290     aResult = aPrismBuilder->Shape();
291
292     // Orienting bounding planes.
293     std::shared_ptr<GeomAPI_Pnt> aCentreOfMass = GeomAlgoAPI_ShapeTools::centreOfMass(theBaseShape);
294     const gp_Pnt& aCentrePnt = aCentreOfMass->impl<gp_Pnt>();
295     gp_Lin aLine(aCentrePnt, aDirVec);
296     IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
297     IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
298     Standard_Real aToParameter = aToIntAna.ParamOnConic(1);
299     Standard_Real aFromParameter = aFromIntAna.ParamOnConic(1);
300     if(aToParameter > aFromParameter) {
301       gp_Vec aVec = aToDir->impl<gp_Dir>();
302       if((aVec * aDirVec) > 0) {
303         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
304         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
305       }
306       aVec = aFromDir->impl<gp_Dir>();
307       if((aVec * aDirVec) < 0) {
308         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
309         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
310       }
311     } else {
312       gp_Vec aVec = aToDir->impl<gp_Dir>();
313       if((aVec * aDirVec) < 0) {
314         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
315         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
316       }
317       aVec = aFromDir->impl<gp_Dir>();
318       if((aVec * aDirVec) > 0) {
319         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
320         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
321       }
322     }
323
324     // Making solids from bounding planes.
325     TopoDS_Shell aToShell, aFromShell;
326     TopoDS_Solid aToSolid, aFromSolid;
327     const TopoDS_Shape& aToShape   = aBoundingToShape->impl<TopoDS_Shape>();
328     const TopoDS_Shape& aFromShape = aBoundingFromShape->impl<TopoDS_Shape>();
329     TopoDS_Face aToFace   = TopoDS::Face(aToShape);
330     TopoDS_Face aFromFace = TopoDS::Face(aFromShape);
331     BRep_Builder aBoundingBuilder;
332     aBoundingBuilder.MakeShell(aToShell);
333     aBoundingBuilder.Add(aToShell, aToShape);
334     aBoundingBuilder.MakeShell(aFromShell);
335     aBoundingBuilder.Add(aFromShell, aFromShape);
336     aBoundingBuilder.MakeSolid(aToSolid);
337     aBoundingBuilder.Add(aToSolid, aToShell);
338     aBoundingBuilder.MakeSolid(aFromSolid);
339     aBoundingBuilder.Add(aFromSolid, aFromShell);
340
341     // Cutting with to plane.
342     BRepAlgoAPI_Cut* aToCutBuilder = new BRepAlgoAPI_Cut(aResult, aToSolid);
343     aToCutBuilder->Build();
344     if(!aToCutBuilder->IsDone()) {
345       return;
346     }
347     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aToCutBuilder)));
348     aResult = aToCutBuilder->Shape();
349     if(aResult.ShapeType() == TopAbs_COMPOUND) {
350       aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
351     }
352     if(aShapeTypeToExp == TopAbs_FACE) {
353       const TopTools_ListOfShape& aToShapes = aToCutBuilder->Modified(aToShape);
354       for(TopTools_ListIteratorOfListOfShape anIt(aToShapes); anIt.More(); anIt.Next()) {
355         GeomShapePtr aGeomSh(new GeomAPI_Shape());
356         aGeomSh->setImpl(new TopoDS_Shape(anIt.Value()));
357         this->addToShape(aGeomSh);
358       }
359     }
360
361     // Cutting with from plane.
362     BRepAlgoAPI_Cut* aFromCutBuilder = new BRepAlgoAPI_Cut(aResult, aFromSolid);
363     aFromCutBuilder->Build();
364     if(!aFromCutBuilder->IsDone()) {
365       return;
366     }
367     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(new GeomAlgoAPI_MakeShape(aFromCutBuilder)));
368     aResult = aFromCutBuilder->Shape();
369     TopoDS_Iterator aCheckIt(aResult);
370     if(!aCheckIt.More()) {
371       return;
372     }
373     if(aResult.ShapeType() == TopAbs_COMPOUND) {
374       aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
375     }
376     if(aShapeTypeToExp == TopAbs_FACE) {
377       const TopTools_ListOfShape& aFromShapes = aFromCutBuilder->Modified(aFromShape);
378       for(TopTools_ListIteratorOfListOfShape anIt(aFromShapes); anIt.More(); anIt.Next()) {
379         GeomShapePtr aGeomSh(new GeomAPI_Shape());
380         aGeomSh->setImpl(new TopoDS_Shape(anIt.Value()));
381         this->addFromShape(aGeomSh);
382       }
383     }
384
385     // Naming for extrusion from vertex, edge.
386     for(TopExp_Explorer anExp(aResult, aShapeTypeToExp); anExp.More(); anExp.Next()) {
387       const TopoDS_Shape& aShape = anExp.Current();
388       GeomShapePtr aGeomSh(new GeomAPI_Shape());
389       if(aShapeTypeToExp == TopAbs_VERTEX) {
390         gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(aShape));
391         IntTools_Context anIntTools;
392         if(anIntTools.IsValidPointForFace(aPnt, aToFace, Precision::Confusion()) == Standard_True) {
393           aGeomSh->setImpl(new TopoDS_Shape(aShape));
394           this->addToShape(aGeomSh);
395         }
396         if(anIntTools.IsValidPointForFace(aPnt, aFromFace, Precision::Confusion()) == Standard_True) {
397           aGeomSh->setImpl(new TopoDS_Shape(aShape));
398           this->addFromShape(aGeomSh);
399         }
400       } else if(aShapeTypeToExp == TopAbs_EDGE) {
401         TopoDS_Edge anEdge = TopoDS::Edge(aShape);
402         BRepLib_CheckCurveOnSurface anEdgeCheck(anEdge, aToFace);
403         anEdgeCheck.Perform();
404         if(anEdgeCheck.MaxDistance() < Precision::Confusion()) {
405           aGeomSh->setImpl(new TopoDS_Shape(aShape));
406           this->addToShape(aGeomSh);
407         }
408         anEdgeCheck.Init(anEdge, aFromFace);
409         anEdgeCheck.Perform();
410         if(anEdgeCheck.MaxDistance() < Precision::Confusion()) {
411           aGeomSh->setImpl(new TopoDS_Shape(aShape));
412           this->addFromShape(aGeomSh);
413         }
414       } else {
415         break;
416       }
417     }
418
419     if(aResult.ShapeType() == TopAbs_COMPOUND) {
420       std::shared_ptr<GeomAPI_Shape> aGeomShape(new GeomAPI_Shape);
421       aGeomShape->setImpl(new TopoDS_Shape(aResult));
422       ListOfShape aCompSolids, aFreeSolids;
423       aGeomShape = GeomAlgoAPI_ShapeTools::combineShapes(aGeomShape,
424                                                          GeomAPI_Shape::COMPSOLID,
425                                                          aCompSolids,
426                                                          aFreeSolids);
427       aResult = aGeomShape->impl<TopoDS_Shape>();
428     }
429   }
430
431   // Setting result.
432   if(aResult.IsNull()) {
433     return;
434   }
435   aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
436   GeomShapePtr aGeomSh(new GeomAPI_Shape());
437   aGeomSh->setImpl(new TopoDS_Shape(aResult));
438   this->setShape(aGeomSh);
439   this->setDone(true);
440 }