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