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