]> SALOME platform Git repositories - modules/shaper.git/blob - src/GeomAlgoAPI/GeomAlgoAPI_Prism.cpp
Salome HOME
eae39590229d4909691e4df9cdf40fc9f8b3093f
[modules/shaper.git] / src / GeomAlgoAPI / GeomAlgoAPI_Prism.cpp
1 // Copyright (C) 2014-2024  CEA, EDF
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 #include "GeomAlgoAPI_Prism.h"
21
22 #include <GeomAPI_Ax1.h>
23 #include <GeomAPI_Dir.h>
24 #include <GeomAPI_Face.h>
25 #include <GeomAPI_Pln.h>
26 #include <GeomAPI_Pnt.h>
27 #include <GeomAPI_ShapeExplorer.h>
28 #include <GeomAPI_XYZ.h>
29
30 #include <GeomAlgoAPI_CompoundBuilder.h>
31 #include <GeomAlgoAPI_DFLoader.h>
32 #include <GeomAlgoAPI_FaceBuilder.h>
33 #include <GeomAlgoAPI_Offset.h>
34 #include <GeomAlgoAPI_Partition.h>
35 #include <GeomAlgoAPI_ShapeTools.h>
36 #include <GeomAlgoAPI_Translation.h>
37
38 #include <Bnd_Box.hxx>
39 #include <BRep_Builder.hxx>
40 #include <BRepAlgoAPI_Cut.hxx>
41 #include <BRepBndLib.hxx>
42 #include <BRepBuilderAPI_FindPlane.hxx>
43 #include <BRepBuilderAPI_Transform.hxx>
44 #include <BRepTools.hxx>
45 #include <Geom_Curve.hxx>
46 #include <Geom2d_Curve.hxx>
47 #include <BRepLib_CheckCurveOnSurface.hxx>
48 #include <BRepPrimAPI_MakePrism.hxx>
49 #include <Geom_Plane.hxx>
50 #include <Geom_RectangularTrimmedSurface.hxx>
51 #include <gp_Pln.hxx>
52 #include <IntAna_IntConicQuad.hxx>
53 #include <IntAna_Quadric.hxx>
54 #include <IntTools_Context.hxx>
55 #include <TopExp_Explorer.hxx>
56 #include <TopoDS.hxx>
57 #include <TopoDS_Edge.hxx>
58 #include <TopoDS_Shell.hxx>
59 #include <TopoDS_Solid.hxx>
60 #include <TopTools_ListIteratorOfListOfShape.hxx>
61
62
63 /// Expand planar face to cover the bounding box if theOriginalShape is planar.
64 /// Otherwise, return the same shape;
65 static GeomShapePtr buildPlanarFace(const GeomShapePtr& theOriginalShape,
66                                     const Bnd_Box& theBaseShapeBB);
67
68 /// Build offset for the given shape.
69 /// If the offset algorithm failed, translate the shape along the direction.
70 static GeomShapePtr buildOffset(const GeomShapePtr& theShape,
71                                 const double theOffset,
72                                 const GeomDirPtr theDirection,
73                                 GeomAlgoAPI_MakeShapeList& theMakeShapeList);
74
75 /// Collect base faces of the prism.
76 static void collectPrismBases(const TopoDS_Shape& theBaseShape,
77                               BRepPrimAPI_MakePrism& thePrismAlgo,
78                               ListOfShape& theBoundaries,
79                               const GeomAPI_Shape::ShapeType theTypeToExp);
80
81 /// Collect all solids which contain boundaries but do not contain bases of prism.
82 static GeomShapePtr collectResults(const GeomMakeShapePtr& theOperation,
83                                    const ListOfShape& theBoundaries,
84                                    const ListOfShape& theShapesToExclude,
85                                    const GeomAPI_Shape::ShapeType theTypeToExp);
86
87 static void storeGenerationHistory(GeomAlgoAPI_Prism* thePrismAlgo,
88                                    const TopoDS_Shape& theBase,
89                                    const TopAbs_ShapeEnum theType,
90                                    BRepPrimAPI_MakePrism* thePrismBuilder);
91
92 static void storeGenerationHistory(GeomAlgoAPI_Prism* thePrismAlgo,
93                                    const TopoDS_Shape& theResult,
94                                    const TopAbs_ShapeEnum theType,
95                                    const TopoDS_Face& theToFace,
96                                    const TopoDS_Face& theFromFace);
97
98 static GeomShapePtr toShape(const TopoDS_Shape& theShape)
99 {
100   GeomShapePtr aShape(new GeomAPI_Shape());
101   aShape->setImpl(new TopoDS_Shape(theShape));
102   return aShape;
103 }
104
105 static void changeOrientationIfNeeded(const TopoDS_Shape& theShape, gp_Vec& theNormal)
106 {
107   TopExp_Explorer anExp(theShape, TopAbs_VERTEX);
108   gp_Pnt aPnt0 = BRep_Tool::Pnt(TopoDS::Vertex(anExp.Current()));
109   gp_Dir aDir01;
110   for (anExp.Next(); anExp.More(); anExp.Next()) {
111     gp_Pnt aPnt1 = BRep_Tool::Pnt(TopoDS::Vertex(anExp.Current()));
112     if (aPnt1.SquareDistance(aPnt0) > Precision::SquareConfusion()) {
113       aDir01 = gp_Dir(gp_Vec(aPnt0, aPnt1));
114       break;
115     }
116   }
117   gp_Vec aNormal;
118   for (; anExp.More(); anExp.Next()) {
119     gp_Pnt aPnt2 = BRep_Tool::Pnt(TopoDS::Vertex(anExp.Current()));
120     if (aPnt2.SquareDistance(aPnt0) > Precision::SquareConfusion()) {
121       aNormal = gp_Vec(aDir01) ^ gp_Vec(aPnt0, aPnt2);
122       if (aNormal.SquareMagnitude() > Precision::SquareConfusion())
123         break;
124     }
125   }
126   if (anExp.More() && aNormal.XYZ().Dot(theNormal.XYZ()) < -Precision::Confusion()) {
127     // directions differ, reverse the normal
128     theNormal.Reverse();
129   }
130 }
131
132 //==================================================================================================
133 GeomAlgoAPI_Prism::GeomAlgoAPI_Prism(const GeomShapePtr theBaseShape,
134                                      const GeomDirPtr   theDirection,
135                                      const GeomShapePtr theToShape,
136                                      const double       theToSize,
137                                      const GeomShapePtr theFromShape,
138                                      const double       theFromSize)
139 {
140   if(!theBaseShape.get() ||
141     (((!theFromShape.get() && !theToShape.get()) ||
142     (theFromShape.get() && theToShape.get() && theFromShape->isEqual(theToShape)))
143     && (theFromSize == -theToSize))) {
144     return;
145   }
146
147   // Getting base shape.
148   const TopoDS_Shape& aBaseShape = theBaseShape->impl<TopoDS_Shape>();
149   GeomAPI_Shape::ShapeType aShapeTypeToExp;
150   switch(aBaseShape.ShapeType()) {
151     case TopAbs_VERTEX:
152       aShapeTypeToExp = GeomAPI_Shape::VERTEX;
153       break;
154     case TopAbs_EDGE:
155     case TopAbs_WIRE:
156       aShapeTypeToExp = GeomAPI_Shape::EDGE;
157       break;
158     case TopAbs_FACE:
159     case TopAbs_SHELL:
160       aShapeTypeToExp = GeomAPI_Shape::FACE;
161       break;
162     case TopAbs_COMPOUND:
163       aShapeTypeToExp = GeomAPI_Shape::COMPOUND;
164       break;
165     default:
166       return;
167   }
168
169   // Getting direction.
170   gp_Vec aBaseVec;
171   std::shared_ptr<GeomAPI_Pnt> aBaseLoc;
172   std::shared_ptr<GeomAPI_Dir> aBaseDir;
173   BRepBuilderAPI_FindPlane aFindPlane(aBaseShape);
174   if(aFindPlane.Found() == Standard_True)
175   {
176     bool checkOrientation = false;
177     Handle(Geom_Plane) aPlane;
178     if(aBaseShape.ShapeType() == TopAbs_FACE || aBaseShape.ShapeType() == TopAbs_SHELL) {
179       TopExp_Explorer anExp(aBaseShape, TopAbs_FACE);
180       const TopoDS_Shape& aFace = anExp.Current();
181       Handle(Geom_Surface) aSurface = BRep_Tool::Surface(TopoDS::Face(aFace));
182       if(aSurface->DynamicType() == STANDARD_TYPE(Geom_RectangularTrimmedSurface)) {
183         Handle(Geom_RectangularTrimmedSurface) aTrimSurface =
184           Handle(Geom_RectangularTrimmedSurface)::DownCast(aSurface);
185         aSurface = aTrimSurface->BasisSurface();
186       }
187       if (aSurface->DynamicType() == STANDARD_TYPE(Geom_Plane))
188         aPlane = Handle(Geom_Plane)::DownCast(aSurface);
189     }
190
191     if (aPlane.IsNull()) {
192       aPlane = aFindPlane.Plane();
193       checkOrientation = true;
194     }
195     gp_Pnt aLoc = aPlane->Axis().Location();
196     aBaseVec = aPlane->Axis().Direction();
197
198     if (checkOrientation) {
199       // to stabilize the result of algorithm, if base shape is a wire, compare the orientation
200       // of calculated plane with the normal vector got iterating on vertices
201       changeOrientationIfNeeded(aBaseShape, aBaseVec);
202     }
203
204     aBaseLoc.reset(new GeomAPI_Pnt(aLoc.X(), aLoc.Y(), aLoc.Z()));
205     aBaseDir.reset(new GeomAPI_Dir(aBaseVec.X(), aBaseVec.Y(), aBaseVec.Z()));
206   }
207   else if (theDirection.get())
208   {
209     aBaseDir = theDirection;
210     aBaseVec = theDirection->impl<gp_Dir>();
211   }
212   else
213   {
214     return;
215   }
216
217   if(!aBaseLoc.get()) {
218     gp_Pnt aLoc;
219     gp_XYZ aDirXYZ = aBaseVec.XYZ();
220     Standard_Real aMinParam = Precision::Infinite();
221     for(TopExp_Explorer anExp(aBaseShape, TopAbs_VERTEX); anExp.More(); anExp.Next()) {
222       const TopoDS_Shape& aVertex = anExp.Current();
223       gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(aVertex));
224       double aParam = aDirXYZ.Dot(aPnt.XYZ());
225       if(aParam < aMinParam) {
226         aMinParam = aParam;
227         aLoc = aPnt;
228       }
229     }
230     aBaseLoc.reset(new GeomAPI_Pnt(aLoc.X(), aLoc.Y(), aLoc.Z()));
231   }
232
233   gp_Vec anExtVec;
234   std::shared_ptr<GeomAPI_Dir> anExtDir;
235   if (theDirection.get())
236   {
237     anExtDir = theDirection;
238     anExtVec = theDirection->impl<gp_Dir>();
239   }
240   else
241   {
242     anExtDir = aBaseDir;
243     anExtVec = aBaseDir->impl<gp_Dir>();
244   }
245
246
247   TopoDS_Shape aResult;
248   const bool isBoundingShapesSet = theFromShape.get() || theToShape.get();
249   if(!isBoundingShapesSet) {
250     buildBySizes(theBaseShape, anExtDir, theToSize, theFromSize, aShapeTypeToExp);
251   } else {
252     GeomShapePtr aBasePlane = GeomAlgoAPI_FaceBuilder::squareFace(aBaseLoc, aBaseDir, 100.0);
253
254     GeomShapePtr aBoundingFromShape = theFromShape ? theFromShape : aBasePlane;
255     GeomShapePtr aBoundingToShape   = theToShape   ? theToShape   : aBasePlane;
256
257     bool isFromShapePlanar = aBoundingFromShape->isPlanar();
258     bool isToShapePlanar   = aBoundingToShape->isPlanar();
259
260     // Set signs of offsets if both bounding shapes are planar
261     if (isFromShapePlanar && isToShapePlanar) {
262       std::shared_ptr<GeomAPI_Pln> aFromPln = GeomAPI_Face(aBoundingFromShape).getPlane();
263       std::shared_ptr<GeomAPI_Pln> aToPln = GeomAPI_Face(aBoundingToShape).getPlane();
264       buildByPlanes(theBaseShape, anExtDir,
265                     aToPln, theToSize,
266                     aFromPln, theFromSize,
267                     aShapeTypeToExp);
268     }
269     else {
270       buildByFaces(theBaseShape, anExtDir,
271                    aBoundingToShape, theToSize, isToShapePlanar,
272                    aBoundingFromShape, theFromSize, isFromShapePlanar,
273                    aShapeTypeToExp);
274     }
275   }
276 }
277
278 //==================================================================================================
279 void GeomAlgoAPI_Prism::buildBySizes(const GeomShapePtr             theBaseShape,
280                                      const GeomDirPtr               theDirection,
281                                      const double                   theToSize,
282                                      const double                   theFromSize,
283                                      const GeomAPI_Shape::ShapeType theTypeToExp)
284 {
285   const TopoDS_Shape& aBaseShape = theBaseShape->impl<TopoDS_Shape>();
286   gp_Vec anExtVec = theDirection->impl<gp_Dir>();
287
288   // Moving base shape.
289   gp_Trsf aTrsf;
290   aTrsf.SetTranslation(anExtVec * -theFromSize);
291   BRepBuilderAPI_Transform* aTransformBuilder =
292       new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
293   if (!aTransformBuilder || !aTransformBuilder->IsDone()) {
294     return;
295   }
296   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
297       new GeomAlgoAPI_MakeShape(aTransformBuilder)));
298   TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
299
300   // Making prism.
301   BRepPrimAPI_MakePrism* aPrismBuilder =
302       new BRepPrimAPI_MakePrism(aMovedBase, anExtVec * (theFromSize + theToSize));
303   if (!aPrismBuilder || !aPrismBuilder->IsDone()) {
304     return;
305   }
306   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
307       new GeomAlgoAPI_MakeShape(aPrismBuilder)));
308   TopoDS_Shape aResult = aPrismBuilder->Shape();
309
310   // Setting naming.
311   if(theTypeToExp == GeomAPI_Shape::COMPOUND) {
312     storeGenerationHistory(this, aMovedBase, TopAbs_EDGE, aPrismBuilder);
313     storeGenerationHistory(this, aMovedBase, TopAbs_FACE, aPrismBuilder);
314   } else {
315     storeGenerationHistory(this, aMovedBase, (TopAbs_ShapeEnum)theTypeToExp, aPrismBuilder);
316   }
317
318   // Setting result.
319   if (!aResult.IsNull()) {
320     aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
321     this->setShape(toShape(aResult));
322     this->setDone(true);
323   }
324 }
325
326 //==================================================================================================
327 void GeomAlgoAPI_Prism::buildByPlanes(const GeomShapePtr             theBaseShape,
328                                       const GeomDirPtr               theDirection,
329                                       const GeomPlanePtr             theToPlane,
330                                       const double                   theToSize,
331                                       const GeomPlanePtr             theFromPlane,
332                                       const double                   theFromSize,
333                                       const GeomAPI_Shape::ShapeType theTypeToExp)
334 {
335   const TopoDS_Shape& aBaseShape = theBaseShape->impl<TopoDS_Shape>();
336   gp_Vec anExtVec = theDirection->impl<gp_Dir>();
337
338   // Moving prism bounding faces according to "from" and "to" sizes.
339   std::shared_ptr<GeomAPI_Pnt> aFromLoc = theFromPlane->location();
340   std::shared_ptr<GeomAPI_Dir> aFromDir = theFromPlane->direction();
341
342   std::shared_ptr<GeomAPI_Pnt> aToLoc = theToPlane->location();
343   std::shared_ptr<GeomAPI_Dir> aToDir = theToPlane->direction();
344
345   std::shared_ptr<GeomAPI_XYZ> anExtDir = theDirection->xyz();
346   bool aSign = aFromLoc->xyz()->dot(anExtDir) > aToLoc->xyz()->dot(anExtDir);
347
348   std::shared_ptr<GeomAPI_Pnt> aFromPnt(
349     new GeomAPI_Pnt(aFromLoc->xyz()->added(anExtDir->multiplied(
350                     aSign ? theFromSize : -theFromSize))));
351
352   std::shared_ptr<GeomAPI_Pnt> aToPnt(
353     new GeomAPI_Pnt(aToLoc->xyz()->added(anExtDir->multiplied(
354                     aSign ? -theToSize : theToSize))));
355
356   // Getting bounding box for base shape.
357   Bnd_Box aBndBox;
358   BRepBndLib::Add(aBaseShape, aBndBox);
359   Standard_Real aXArr[2] = {aBndBox.CornerMin().X(), aBndBox.CornerMax().X()};
360   Standard_Real aYArr[2] = {aBndBox.CornerMin().Y(), aBndBox.CornerMax().Y()};
361   Standard_Real aZArr[2] = {aBndBox.CornerMin().Z(), aBndBox.CornerMax().Z()};
362   gp_Pnt aPoints[8];
363   int aNum = 0;
364   for(int i = 0; i < 2; i++) {
365     for(int j = 0; j < 2; j++) {
366       for(int k = 0; k < 2; k++) {
367         aPoints[aNum] = gp_Pnt(aXArr[i], aYArr[j], aZArr[k]);
368         aNum++;
369       }
370     }
371   }
372
373   // Project points to bounding planes. Search max distance to them.
374   IntAna_Quadric aBndToQuadric(gp_Pln(aToPnt->impl<gp_Pnt>(), aToDir->impl<gp_Dir>()));
375   IntAna_Quadric aBndFromQuadric(gp_Pln(aFromPnt->impl<gp_Pnt>(), aFromDir->impl<gp_Dir>()));
376   Standard_Real aMaxToDist = 0, aMaxFromDist = 0;
377   for(int i = 0; i < 8; i++) {
378     gp_Lin aLine(aPoints[i], anExtVec);
379     IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
380     IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
381     if(aToIntAna.NbPoints() == 0 || aFromIntAna.NbPoints() == 0) {
382       return;
383     }
384     const gp_Pnt& aPntOnToFace = aToIntAna.Point(1);
385     const gp_Pnt& aPntOnFromFace = aFromIntAna.Point(1);
386     if(aPoints[i].Distance(aPntOnToFace) > aMaxToDist) {
387       aMaxToDist = aPoints[i].Distance(aPntOnToFace);
388     }
389     if(aPoints[i].Distance(aPntOnFromFace) > aMaxFromDist) {
390       aMaxFromDist = aPoints[i].Distance(aPntOnFromFace);
391     }
392   }
393
394   // We added 1 just to be sure that prism is long enough for boolean operation.
395   double aPrismLength = aMaxToDist + aMaxFromDist + 1;
396
397   // Moving base shape.
398   gp_Trsf aTrsf;
399   aTrsf.SetTranslation(anExtVec * -aPrismLength);
400   BRepBuilderAPI_Transform* aTransformBuilder = new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
401   if(!aTransformBuilder || !aTransformBuilder->IsDone()) {
402     return;
403   }
404   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
405     new GeomAlgoAPI_MakeShape(aTransformBuilder)));
406   TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
407
408   // Making prism.
409   BRepPrimAPI_MakePrism* aPrismBuilder =
410     new BRepPrimAPI_MakePrism(aMovedBase, anExtVec * 2 * aPrismLength);
411   if(!aPrismBuilder || !aPrismBuilder->IsDone()) {
412     return;
413   }
414   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
415     new GeomAlgoAPI_MakeShape(aPrismBuilder)));
416   TopoDS_Shape aResult = aPrismBuilder->Shape();
417
418   BRepBndLib::Add(aResult, aBndBox);
419   aBndBox.Add(aFromPnt->impl<gp_Pnt>());
420   aBndBox.Add(aToPnt->impl<gp_Pnt>());
421   Standard_Real aBndBoxSize = aBndBox.CornerMin().Distance(aBndBox.CornerMax());
422
423   // Orienting bounding planes.
424   std::shared_ptr<GeomAPI_Pnt> aCentreOfMass = GeomAlgoAPI_ShapeTools::centreOfMass(theBaseShape);
425   const gp_Pnt& aCentrePnt = aCentreOfMass->impl<gp_Pnt>();
426   gp_Lin aLine(aCentrePnt, anExtVec);
427   IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
428   IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
429   Standard_Real aToParameter = aToIntAna.ParamOnConic(1);
430   Standard_Real aFromParameter = aFromIntAna.ParamOnConic(1);
431   if(aToParameter > aFromParameter) {
432     gp_Vec aVec = aToDir->impl<gp_Dir>();
433     if((aVec * anExtVec) > 0)
434       aToDir->setImpl(new gp_Dir(aVec.Reversed()));
435     aVec = aFromDir->impl<gp_Dir>();
436     if((aVec * anExtVec) < 0)
437       aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
438   } else {
439     gp_Vec aVec = aToDir->impl<gp_Dir>();
440     if((aVec * anExtVec) < 0)
441       aToDir->setImpl(new gp_Dir(aVec.Reversed()));
442     aVec = aFromDir->impl<gp_Dir>();
443     if((aVec * anExtVec) > 0)
444       aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
445   }
446
447   static const double THE_FACE_SIZE_COEFF = 10.0;
448   GeomShapePtr aBoundingFromShape =
449       GeomAlgoAPI_FaceBuilder::squareFace(aFromPnt, aFromDir, THE_FACE_SIZE_COEFF * aBndBoxSize);
450   GeomShapePtr aBoundingToShape =
451       GeomAlgoAPI_FaceBuilder::squareFace(aToPnt, aToDir, THE_FACE_SIZE_COEFF * aBndBoxSize);
452
453   // bounding planes
454   const TopoDS_Shape& aToShape   = aBoundingToShape->impl<TopoDS_Shape>();
455   const TopoDS_Shape& aFromShape = aBoundingFromShape->impl<TopoDS_Shape>();
456   TopoDS_Face aToFace   = TopoDS::Face(aToShape);
457   TopoDS_Face aFromFace = TopoDS::Face(aFromShape);
458
459   // Solid based on "To" bounding plane
460   gp_Vec aNormal = aToDir->impl<gp_Dir>();
461   BRepPrimAPI_MakePrism* aToPrismBuilder =
462       new BRepPrimAPI_MakePrism(aToShape, aNormal * (-2.0 * aBndBoxSize));
463   if (!aToPrismBuilder || !aToPrismBuilder->IsDone()) {
464     return;
465   }
466   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
467     new GeomAlgoAPI_MakeShape(aToPrismBuilder)));
468   TopoDS_Shape aToSolid = aToPrismBuilder->Shape();
469
470   // Cutting with to plane.
471   BRepAlgoAPI_Cut* aToCutBuilder = new BRepAlgoAPI_Cut(aResult, aToSolid);
472   aToCutBuilder->Build();
473   if(!aToCutBuilder->IsDone()) {
474     return;
475   }
476   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
477     new GeomAlgoAPI_MakeShape(aToCutBuilder)));
478   aResult = aToCutBuilder->Shape();
479   if(aResult.ShapeType() == TopAbs_COMPOUND) {
480     aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
481   }
482   if (theTypeToExp == GeomAPI_Shape::FACE || theTypeToExp == GeomAPI_Shape::COMPOUND) {
483     TopTools_ListOfShape aPrismShapes = aToPrismBuilder->Modified(aToShape);
484     if (aPrismShapes.IsEmpty())
485       aPrismShapes.Append(aToShape);
486     for (TopTools_ListIteratorOfListOfShape anIt1(aPrismShapes); anIt1.More(); anIt1.Next()) {
487       const TopTools_ListOfShape& aToShapes = aToCutBuilder->Modified(anIt1.Value());
488       for (TopTools_ListIteratorOfListOfShape anIt2(aToShapes); anIt2.More(); anIt2.Next()) {
489         GeomShapePtr aGeomSh = toShape(anIt2.Value());
490         fixOrientation(aGeomSh);
491         this->addToShape(aGeomSh);
492       }
493     }
494   }
495
496   // Solid based on "From" bounding plane
497   aNormal = aFromDir->impl<gp_Dir>();
498   BRepPrimAPI_MakePrism* aFromPrismBuilder =
499       new BRepPrimAPI_MakePrism(aFromShape, aNormal * (-2.0 * aBndBoxSize));
500   if (!aFromPrismBuilder || !aFromPrismBuilder->IsDone()) {
501     return;
502   }
503   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
504     new GeomAlgoAPI_MakeShape(aFromPrismBuilder)));
505   TopoDS_Shape aFromSolid = aFromPrismBuilder->Shape();
506
507   // Cutting with from plane.
508   BRepAlgoAPI_Cut* aFromCutBuilder = new BRepAlgoAPI_Cut(aResult, aFromSolid);
509   aFromCutBuilder->Build();
510   if(!aFromCutBuilder->IsDone()) {
511     return;
512   }
513   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
514     new GeomAlgoAPI_MakeShape(aFromCutBuilder)));
515   aResult = aFromCutBuilder->Shape();
516   TopoDS_Iterator aCheckIt(aResult);
517   if(!aCheckIt.More()) {
518     return;
519   }
520   if(aResult.ShapeType() == TopAbs_COMPOUND) {
521     aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
522   }
523   if (theTypeToExp == GeomAPI_Shape::FACE || theTypeToExp == GeomAPI_Shape::COMPOUND) {
524     TopTools_ListOfShape aPrismShapes = aFromPrismBuilder->Modified(aFromShape);
525     if (aPrismShapes.IsEmpty())
526       aPrismShapes.Append(aFromShape);
527     for (TopTools_ListIteratorOfListOfShape anIt1(aPrismShapes); anIt1.More(); anIt1.Next()) {
528       const TopTools_ListOfShape& aFromShapes = aFromCutBuilder->Modified(anIt1.Value());
529       for (TopTools_ListIteratorOfListOfShape anIt2(aFromShapes); anIt2.More(); anIt2.Next()) {
530         GeomShapePtr aGeomSh = toShape(anIt2.Value());
531         fixOrientation(aGeomSh);
532         this->addFromShape(aGeomSh);
533       }
534     }
535   }
536
537   // Naming for extrusion from vertex, edge.
538   if(theTypeToExp == GeomAPI_Shape::COMPOUND) {
539     storeGenerationHistory(this, aResult, TopAbs_EDGE, aToFace, aFromFace);
540     storeGenerationHistory(this, aResult, TopAbs_FACE, aToFace, aFromFace);
541   } else {
542     storeGenerationHistory(this, aResult, (TopAbs_ShapeEnum)theTypeToExp, aToFace, aFromFace);
543   }
544
545   if(aResult.ShapeType() == TopAbs_COMPOUND) {
546     GeomShapePtr aGeomShape = toShape(aResult);
547     ListOfShape aResults;
548     aGeomShape = GeomAlgoAPI_ShapeTools::combineShapes(aGeomShape,
549                                                         GeomAPI_Shape::COMPSOLID,
550                                                         aResults);
551     aResult = aGeomShape->impl<TopoDS_Shape>();
552   }
553
554   // Setting result.
555   if (!aResult.IsNull()) {
556     aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
557     this->setShape(toShape(aResult));
558     this->setDone(true);
559   }
560 }
561
562 //==================================================================================================
563 void GeomAlgoAPI_Prism::buildByFaces(const GeomShapePtr             theBaseShape,
564                                      const GeomDirPtr               theDirection,
565                                      const GeomShapePtr             theToShape,
566                                      const double                   theToSize,
567                                      const bool                     theToIsPlanar,
568                                      const GeomShapePtr             theFromShape,
569                                      const double                   theFromSize,
570                                      const bool                     theFromIsPlanar,
571                                      const GeomAPI_Shape::ShapeType theTypeToExp)
572 {
573   const TopoDS_Shape& aBaseShape = theBaseShape->impl<TopoDS_Shape>();
574   gp_Vec anExtVec = theDirection->impl<gp_Dir>();
575
576   // Moving prism bounding faces according to "from" and "to" sizes.
577   GeomShapePtr aBoundingFromShape = buildOffset(theFromShape, -theFromSize, theDirection, *this);
578   GeomShapePtr aBoundingToShape   = buildOffset(theToShape, theToSize, theDirection, *this);
579
580   // Bounding box for shapes used in prism building.
581   Bnd_Box aBndBox;
582   BRepBndLib::Add(aBaseShape, aBndBox);
583   BRepBndLib::Add(aBoundingFromShape->impl<TopoDS_Shape>(), aBndBox);
584   BRepBndLib::Add(aBoundingToShape->impl<TopoDS_Shape>(), aBndBox);
585   double aPrismLength = 2.0 * aBndBox.CornerMin().Distance(aBndBox.CornerMax());
586
587   // Prism building.
588   gp_Trsf aTrsf;
589   aTrsf.SetTranslation(anExtVec * -aPrismLength);
590   BRepBuilderAPI_Transform* aTransformBuilder = new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
591   if (!aTransformBuilder || !aTransformBuilder->IsDone()) {
592     return;
593   }
594   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
595       new GeomAlgoAPI_MakeShape(aTransformBuilder)));
596   TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
597
598   // Making prism.
599   BRepPrimAPI_MakePrism* aPrismBuilder =
600       new BRepPrimAPI_MakePrism(aMovedBase, anExtVec * 2 * aPrismLength);
601   if (!aPrismBuilder || !aPrismBuilder->IsDone()) {
602     return;
603   }
604   this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
605       new GeomAlgoAPI_MakeShape(aPrismBuilder)));
606
607   GeomShapePtr aResult = toShape(aPrismBuilder->Shape());
608
609   // Prism generatrix
610   ListOfShape aPrismBaseFaces;
611   collectPrismBases(aMovedBase, *aPrismBuilder, aPrismBaseFaces, theTypeToExp);
612
613   // Build planar faces intersecting the prism fully.
614   BRepBndLib::Add(aResult->impl<TopoDS_Shape>(), aBndBox);
615   aBoundingFromShape = buildPlanarFace(aBoundingFromShape, aBndBox);
616   aBoundingToShape   = buildPlanarFace(aBoundingToShape, aBndBox);
617
618   // Perform partition.
619   ListOfShape anObjects, aTools;
620   anObjects.push_back(aResult);
621   aTools.push_back(aBoundingFromShape);
622   aTools.push_back(aBoundingToShape);
623
624   GeomMakeShapePtr aPartition(new GeomAlgoAPI_Partition(anObjects, aTools));
625   if (!aPartition->isDone())
626     return;
627
628   this->appendAlgo(aPartition);
629
630   // Collect pieces of boundary shapes, split by Partition.
631   if (theFromIsPlanar) {
632     ListOfShape anImagesFrom;
633     aPartition->modified(aBoundingFromShape, anImagesFrom);
634     for (ListOfShape::iterator anIt = anImagesFrom.begin(); anIt != anImagesFrom.end(); ++anIt)
635       addFromShape(*anIt);
636   }
637
638   if (theToIsPlanar) {
639     ListOfShape anImagesTo;
640     aPartition->modified(aBoundingToShape, anImagesTo);
641     for (ListOfShape::iterator anIt = anImagesTo.begin(); anIt != anImagesTo.end(); ++anIt)
642       addToShape(*anIt);
643   }
644
645   // Collect results which have both boundaries, selected for extrusion,
646   // but which do not contain top and bottom faces of the prism
647   // (these faces are treated as infinitely distant).
648   aResult = collectResults(aPartition, aTools, aPrismBaseFaces, theTypeToExp);
649   if (aResult && aResult->shapeType() == GeomAPI_Shape::COMPOUND) {
650     ListOfShape aResults;
651     aResult = GeomAlgoAPI_ShapeTools::combineShapes(aResult,
652         theTypeToExp == GeomAPI_Shape::EDGE ? GeomAPI_Shape::SHELL : GeomAPI_Shape::COMPSOLID,
653         aResults);
654
655     if (aResults.size() > 1 &&
656        (GeomAlgoAPI_ShapeTools::hasSharedTopology(aResults, GeomAPI_Shape::EDGE) ||
657         GeomAlgoAPI_ShapeTools::hasSharedTopology(aResults, GeomAPI_Shape::VERTEX))) {
658       // results shuold not have shared topology
659       aResult = GeomShapePtr();
660     }
661   }
662
663   if (aResult) {
664     this->setShape(aResult);
665     this->setDone(true);
666   }
667 }
668
669
670 // Auxilary functions:
671 //==================================================================================================
672 GeomShapePtr buildPlanarFace(const GeomShapePtr& theOriginalShape,
673                              const Bnd_Box& theBaseShapeBB)
674 {
675   GeomPlanePtr aPlane = GeomAPI_Face(theOriginalShape).getPlane();
676   if (!aPlane)
677     return theOriginalShape;
678
679   gp_Pnt aCornerMin = theBaseShapeBB.CornerMin();
680   gp_Pnt aCornerMax = theBaseShapeBB.CornerMax();
681   double aSize = aCornerMin.SquareDistance(aCornerMax);
682
683   gp_Pnt aLocation = aPlane->location()->impl<gp_Pnt>();
684
685   gp_Pnt aCurPnt;
686   for (int x = 0; x < 2; ++x) {
687     aCurPnt.SetX(x == 0 ? aCornerMin.X() : aCornerMax.X());
688     for (int y = 0; y < 2; ++y) {
689       aCurPnt.SetY(y == 0 ? aCornerMin.Y() : aCornerMax.Y());
690       for (int z = 0; z < 2; ++z) {
691         aCurPnt.SetZ(z == 0 ? aCornerMin.Z() : aCornerMax.Z());
692         double aDist = aCurPnt.SquareDistance(aLocation);
693         if (aDist > aSize)
694           aSize = aDist;
695       }
696     }
697   }
698
699   aSize = Sqrt(aSize);
700   return GeomAlgoAPI_FaceBuilder::squareFace(aPlane, 2.0 * aSize);
701 }
702
703 //==================================================================================================
704 GeomShapePtr buildOffset(const GeomShapePtr& theShape,
705                          const double theOffset,
706                          const GeomDirPtr theDirection,
707                          GeomAlgoAPI_MakeShapeList& theMakeShapeList)
708 {
709   if (Abs(theOffset) < Precision::Confusion())
710     return theShape; // no need zero offset
711
712   GeomMakeShapePtr anAlgo(new GeomAlgoAPI_Offset(theShape, theOffset));
713   if (!anAlgo->isDone()) {
714     // offset not done, perform translation
715     std::shared_ptr<GeomAPI_Ax1> anAxis(new GeomAPI_Ax1());
716     anAxis->setDir(theDirection);
717     anAlgo.reset(new GeomAlgoAPI_Translation(theShape, anAxis, theOffset));
718   }
719
720   GeomShapePtr aResult = theShape;
721   if (anAlgo->isDone()) {
722     theMakeShapeList.appendAlgo(anAlgo);
723     aResult = anAlgo->shape();
724   }
725   return aResult;
726 }
727
728 //==================================================================================================
729 void collectPrismBases(const TopoDS_Shape& theBaseShape,
730                        BRepPrimAPI_MakePrism& thePrismAlgo,
731                        ListOfShape& theBoundaries,
732                        const GeomAPI_Shape::ShapeType theTypeToExp)
733 {
734   for (TopExp_Explorer anExp(theBaseShape, (TopAbs_ShapeEnum)theTypeToExp);
735        anExp.More(); anExp.Next()) {
736     theBoundaries.push_back(toShape(thePrismAlgo.FirstShape(anExp.Current())));
737     theBoundaries.push_back(toShape(thePrismAlgo.LastShape(anExp.Current())));
738   }
739 }
740
741 //==================================================================================================
742 typedef std::set<GeomShapePtr, GeomAPI_Shape::Comparator> SetOfShape;
743
744 bool isShapeApplicable(const GeomShapePtr& theSolid,
745                        const std::list<ListOfShape>& theShapesToExist,
746                        const SetOfShape& theShapesToExclude,
747                        const GeomAPI_Shape::ShapeType theTypeToExp)
748 {
749   SetOfShape aFaces;
750   for (GeomAPI_ShapeExplorer aFExp(theSolid, theTypeToExp);
751        aFExp.more(); aFExp.next()) {
752     GeomShapePtr aCurrent = aFExp.current();
753     if (theShapesToExclude.find(aCurrent) != theShapesToExclude.end())
754       return false;
755     aFaces.insert(aCurrent);
756   }
757
758   // check all faces are in solid
759   bool isApplicable = true;
760   for (std::list<ListOfShape>::const_iterator it1 = theShapesToExist.begin();
761        it1 != theShapesToExist.end() && isApplicable; ++it1) {
762     ListOfShape::const_iterator it2 = it1->begin();
763     for (; it2 != it1->end(); ++it2)
764       if (aFaces.find(*it2) != aFaces.end())
765         break;
766     isApplicable = it2 != it1->end();
767   }
768   return isApplicable;
769 }
770
771 void collectModified(const GeomMakeShapePtr& theOperation,
772                      const ListOfShape& theShapes,
773                      std::list<ListOfShape>& theModified)
774 {
775   for (ListOfShape::const_iterator anIt = theShapes.begin();
776        anIt != theShapes.end(); ++anIt) {
777     theModified.push_back(ListOfShape());
778     theOperation->modified(*anIt, theModified.back());
779     theOperation->generated(*anIt, theModified.back());
780     theModified.back().push_back(*anIt);
781   }
782 }
783
784 GeomShapePtr collectResults(const GeomMakeShapePtr& theOperation,
785                             const ListOfShape& theBoundaries,
786                             const ListOfShape& theShapesToExclude,
787                             const GeomAPI_Shape::ShapeType theTypeToExp)
788 {
789   ListOfShape aResults;
790
791   // collect modified shapes
792   std::list<ListOfShape> aModifiedBoundaries;
793   collectModified(theOperation, theBoundaries, aModifiedBoundaries);
794
795   std::list<ListOfShape> aModifiedExclude;
796   collectModified(theOperation, theShapesToExclude, aModifiedExclude);
797   SetOfShape aTabooShapes;
798   for (std::list<ListOfShape>::iterator anIt = aModifiedExclude.begin();
799        anIt != aModifiedExclude.end(); ++anIt)
800     aTabooShapes.insert(anIt->begin(), anIt->end());
801
802   // type of sub-shapes to explode
803   GeomAPI_Shape::ShapeType aSubshapeType;
804   switch (theTypeToExp) {
805   case GeomAPI_Shape::VERTEX:
806     aSubshapeType = GeomAPI_Shape::EDGE;
807     break;
808   case GeomAPI_Shape::EDGE:
809     aSubshapeType = GeomAPI_Shape::FACE;
810     break;
811   case GeomAPI_Shape::FACE:
812     aSubshapeType = GeomAPI_Shape::SOLID;
813     break;
814   default:
815     aSubshapeType = GeomAPI_Shape::COMPOUND;
816   }
817
818   // search applicable solids
819   GeomShapePtr anOperationResult = theOperation->shape();
820   for (GeomAPI_ShapeExplorer anExp(anOperationResult, aSubshapeType);
821         anExp.more(); anExp.next()) {
822     if (isShapeApplicable(anExp.current(), aModifiedBoundaries, aTabooShapes, theTypeToExp))
823       aResults.push_back(anExp.current());
824   }
825
826   GeomShapePtr aResult;
827   if (aResults.size() == 1)
828     aResult = aResults.front();
829   else if (!aResults.empty())
830     aResult = GeomAlgoAPI_CompoundBuilder::compound(aResults);
831   return aResult;
832 }
833
834 //==================================================================================================
835 void storeGenerationHistory(GeomAlgoAPI_Prism* thePrismAlgo,
836                             const TopoDS_Shape& theBase,
837                             const TopAbs_ShapeEnum theType,
838                             BRepPrimAPI_MakePrism* thePrismBuilder)
839 {
840   for(TopExp_Explorer anExp(theBase, theType); anExp.More(); anExp.Next()) {
841     const TopoDS_Shape& aShape = anExp.Current();
842     GeomShapePtr aFromShape(new GeomAPI_Shape), aToShape(new GeomAPI_Shape);
843     aFromShape->setImpl(new TopoDS_Shape(thePrismBuilder->FirstShape(aShape)));
844     aToShape->setImpl(new TopoDS_Shape(thePrismBuilder->LastShape(aShape)));
845     thePrismAlgo->fixOrientation(aFromShape);
846     thePrismAlgo->fixOrientation(aToShape);
847     thePrismAlgo->addFromShape(aFromShape);
848     thePrismAlgo->addToShape(aToShape);
849   }
850 }
851
852 //==================================================================================================
853 void storeGenerationHistory(GeomAlgoAPI_Prism* thePrismAlgo,
854                             const TopoDS_Shape& theResult,
855                             const TopAbs_ShapeEnum theType,
856                             const TopoDS_Face& theToFace,
857                             const TopoDS_Face& theFromFace)
858 {
859   for(TopExp_Explorer anExp(theResult, theType); anExp.More(); anExp.Next()) {
860     const TopoDS_Shape& aShape = anExp.Current();
861     GeomShapePtr aGeomSh(new GeomAPI_Shape());
862     if(theType == TopAbs_VERTEX) {
863       gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(aShape));
864       IntTools_Context anIntTools;
865       if(anIntTools.IsValidPointForFace(aPnt,
866           theToFace, Precision::Confusion()) == Standard_True) {
867         aGeomSh->setImpl(new TopoDS_Shape(aShape));
868         thePrismAlgo->fixOrientation(aGeomSh);
869         thePrismAlgo->addToShape(aGeomSh);
870       }
871       if(anIntTools.IsValidPointForFace(aPnt,
872           theFromFace, Precision::Confusion()) == Standard_True) {
873         aGeomSh->setImpl(new TopoDS_Shape(aShape));
874         thePrismAlgo->fixOrientation(aGeomSh);
875         thePrismAlgo->addFromShape(aGeomSh);
876       }
877     } else if(theType == TopAbs_EDGE) {
878       TopoDS_Edge anEdge = TopoDS::Edge(aShape);
879       BRepLib_CheckCurveOnSurface anEdgeCheck(anEdge, theToFace);
880       anEdgeCheck.Perform();
881       if(anEdgeCheck.MaxDistance() < Precision::Confusion()) {
882         aGeomSh->setImpl(new TopoDS_Shape(aShape));
883         thePrismAlgo->fixOrientation(aGeomSh);
884         thePrismAlgo->addToShape(aGeomSh);
885       }
886       anEdgeCheck.Init(anEdge, theFromFace);
887       anEdgeCheck.Perform();
888       if(anEdgeCheck.MaxDistance() < Precision::Confusion()) {
889         aGeomSh->setImpl(new TopoDS_Shape(aShape));
890         thePrismAlgo->fixOrientation(aGeomSh);
891         thePrismAlgo->addFromShape(aGeomSh);
892       }
893     } else {
894       break;
895     }
896   }
897 }