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