Salome HOME
2576423e74ceb37811b60bbd5cd556b347a1cf96
[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 aBaseVec;
148   std::shared_ptr<GeomAPI_Pnt> aBaseLoc;
149   std::shared_ptr<GeomAPI_Dir> aBaseDir;
150   BRepBuilderAPI_FindPlane aFindPlane(aBaseShape);
151   if(aFindPlane.Found() == Standard_True)
152   {
153     Handle(Geom_Plane) aPlane;
154     if(aBaseShape.ShapeType() == TopAbs_FACE || aBaseShape.ShapeType() == TopAbs_SHELL) {
155       TopExp_Explorer anExp(aBaseShape, TopAbs_FACE);
156       const TopoDS_Shape& aFace = anExp.Current();
157       Handle(Geom_Surface) aSurface = BRep_Tool::Surface(TopoDS::Face(aFace));
158       if(aSurface->DynamicType() == STANDARD_TYPE(Geom_RectangularTrimmedSurface)) {
159         Handle(Geom_RectangularTrimmedSurface) aTrimSurface =
160           Handle(Geom_RectangularTrimmedSurface)::DownCast(aSurface);
161         aSurface = aTrimSurface->BasisSurface();
162       }
163       if(aSurface->DynamicType() != STANDARD_TYPE(Geom_Plane)) {
164         return;
165       }
166       aPlane = Handle(Geom_Plane)::DownCast(aSurface);
167     } else {
168       aPlane = aFindPlane.Plane();
169     }
170     gp_Pnt aLoc = aPlane->Axis().Location();
171     aBaseVec = aPlane->Axis().Direction();
172     aBaseLoc.reset(new GeomAPI_Pnt(aLoc.X(), aLoc.Y(), aLoc.Z()));
173     aBaseDir.reset(new GeomAPI_Dir(aBaseVec.X(), aBaseVec.Y(), aBaseVec.Z()));
174   }
175   else if (theDirection.get())
176   {
177     aBaseDir = theDirection;
178     aBaseVec = theDirection->impl<gp_Dir>();
179   }
180   else
181   {
182     return;
183   }
184
185   if(!aBaseLoc.get()) {
186     gp_Pnt aLoc;
187     gp_XYZ aDirXYZ = aBaseVec.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
201   GeomShapePtr aBasePlane = GeomAlgoAPI_FaceBuilder::planarFace(aBaseLoc, aBaseDir);
202
203   gp_Vec anExtVec;
204   std::shared_ptr<GeomAPI_Dir> anExtDir;
205   if (theDirection.get())
206   {
207     anExtDir = theDirection;
208     anExtVec = theDirection->impl<gp_Dir>();
209   }
210   else
211   {
212     anExtDir = aBaseDir;
213     anExtVec = aBaseDir->impl<gp_Dir>();
214   }
215
216
217   TopoDS_Shape aResult;
218   const bool isBoundingShapesSet = theFromShape.get() || theToShape.get();
219   if(!isBoundingShapesSet) {
220     // Moving base shape.
221     gp_Trsf aTrsf;
222     aTrsf.SetTranslation(anExtVec * -theFromSize);
223     BRepBuilderAPI_Transform* aTransformBuilder =
224       new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
225     if(!aTransformBuilder) {
226       return;
227     }
228     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
229       new GeomAlgoAPI_MakeShape(aTransformBuilder)));
230     if(!aTransformBuilder->IsDone()) {
231       return;
232     }
233     TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
234
235     // Making prism.
236     BRepPrimAPI_MakePrism* aPrismBuilder =
237       new BRepPrimAPI_MakePrism(aMovedBase, anExtVec * (theFromSize + theToSize));
238     if(!aPrismBuilder) {
239       return;
240     }
241     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
242       new GeomAlgoAPI_MakeShape(aPrismBuilder)));
243     if(!aPrismBuilder->IsDone()) {
244       return;
245     }
246     aResult = aPrismBuilder->Shape();
247
248     // Setting naming.
249     if(aShapeTypeToExp == TopAbs_COMPOUND) {
250       storeGenerationHistory(this, aMovedBase, TopAbs_EDGE, aPrismBuilder);
251       storeGenerationHistory(this, aMovedBase, TopAbs_FACE, aPrismBuilder);
252     } else {
253       storeGenerationHistory(this, aMovedBase, aShapeTypeToExp, aPrismBuilder);
254     }
255   } else {
256     GeomShapePtr aBoundingFromShape = theFromShape ? theFromShape : aBasePlane;
257     GeomShapePtr aBoundingToShape   = theToShape   ? theToShape   : aBasePlane;
258
259     // Moving prism bounding faces according to "from" and "to" sizes.
260     std::shared_ptr<GeomAPI_Pln> aFromPln = GeomAPI_Face(aBoundingFromShape).getPlane();
261     std::shared_ptr<GeomAPI_Pnt> aFromLoc = aFromPln->location();
262     std::shared_ptr<GeomAPI_Dir> aFromDir = aFromPln->direction();
263
264     std::shared_ptr<GeomAPI_Pln> aToPln = GeomAPI_Face(aBoundingToShape).getPlane();
265     std::shared_ptr<GeomAPI_Pnt> aToLoc = aToPln->location();
266     std::shared_ptr<GeomAPI_Dir> aToDir = aToPln->direction();
267
268     bool aSign = aFromLoc->xyz()->dot(anExtDir->xyz()) > aToLoc->xyz()->dot(anExtDir->xyz());
269
270     std::shared_ptr<GeomAPI_Pnt> aFromPnt(
271       new GeomAPI_Pnt(aFromLoc->xyz()->added(anExtDir->xyz()->multiplied(
272                       aSign ? theFromSize : -theFromSize))));
273     aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
274
275     std::shared_ptr<GeomAPI_Pnt> aToPnt(
276       new GeomAPI_Pnt(aToLoc->xyz()->added(anExtDir->xyz()->multiplied(
277                       aSign ? -theToSize : theToSize))));
278     aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
279
280     // Getting bounding box for base shape.
281     Bnd_Box aBndBox;
282     BRepBndLib::Add(aBaseShape, aBndBox);
283     Standard_Real aXArr[2] = {aBndBox.CornerMin().X(), aBndBox.CornerMax().X()};
284     Standard_Real aYArr[2] = {aBndBox.CornerMin().Y(), aBndBox.CornerMax().Y()};
285     Standard_Real aZArr[2] = {aBndBox.CornerMin().Z(), aBndBox.CornerMax().Z()};
286     gp_Pnt aPoints[8];
287     int aNum = 0;
288     for(int i = 0; i < 2; i++) {
289       for(int j = 0; j < 2; j++) {
290         for(int k = 0; k < 2; k++) {
291           aPoints[aNum] = gp_Pnt(aXArr[i], aYArr[j], aZArr[k]);
292           aNum++;
293         }
294       }
295     }
296
297     // Project points to bounding planes. Search max distance to them.
298     IntAna_Quadric aBndToQuadric(gp_Pln(aToPnt->impl<gp_Pnt>(), aToDir->impl<gp_Dir>()));
299     IntAna_Quadric aBndFromQuadric(gp_Pln(aFromPnt->impl<gp_Pnt>(), aFromDir->impl<gp_Dir>()));
300     Standard_Real aMaxToDist = 0, aMaxFromDist = 0;
301     for(int i = 0; i < 8; i++) {
302       gp_Lin aLine(aPoints[i], anExtVec);
303       IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
304       IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
305       if(aToIntAna.NbPoints() == 0 || aFromIntAna.NbPoints() == 0) {
306         return;
307       }
308       const gp_Pnt& aPntOnToFace = aToIntAna.Point(1);
309       const gp_Pnt& aPntOnFromFace = aFromIntAna.Point(1);
310       if(aPoints[i].Distance(aPntOnToFace) > aMaxToDist) {
311         aMaxToDist = aPoints[i].Distance(aPntOnToFace);
312       }
313       if(aPoints[i].Distance(aPntOnFromFace) > aMaxFromDist) {
314         aMaxFromDist = aPoints[i].Distance(aPntOnFromFace);
315       }
316     }
317
318     // We added 1 just to be sure that prism is long enough for boolean operation.
319     double aPrismLength = aMaxToDist + aMaxFromDist + 1;
320
321     // Moving base shape.
322     gp_Trsf aTrsf;
323     aTrsf.SetTranslation(anExtVec * -aPrismLength);
324     BRepBuilderAPI_Transform* aTransformBuilder = new BRepBuilderAPI_Transform(aBaseShape, aTrsf);
325     if(!aTransformBuilder) {
326       return;
327     }
328     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
329       new GeomAlgoAPI_MakeShape(aTransformBuilder)));
330     if(!aTransformBuilder->IsDone()) {
331       return;
332     }
333     TopoDS_Shape aMovedBase = aTransformBuilder->Shape();
334
335     // Making prism.
336     BRepPrimAPI_MakePrism* aPrismBuilder =
337       new BRepPrimAPI_MakePrism(aMovedBase, anExtVec * 2 * aPrismLength);
338     if(!aPrismBuilder) {
339       return;
340     }
341     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
342       new GeomAlgoAPI_MakeShape(aPrismBuilder)));
343     if(!aPrismBuilder->IsDone()) {
344       return;
345     }
346     aResult = aPrismBuilder->Shape();
347
348     // Orienting bounding planes.
349     std::shared_ptr<GeomAPI_Pnt> aCentreOfMass = GeomAlgoAPI_ShapeTools::centreOfMass(theBaseShape);
350     const gp_Pnt& aCentrePnt = aCentreOfMass->impl<gp_Pnt>();
351     gp_Lin aLine(aCentrePnt, anExtVec);
352     IntAna_IntConicQuad aToIntAna(aLine, aBndToQuadric);
353     IntAna_IntConicQuad aFromIntAna(aLine, aBndFromQuadric);
354     Standard_Real aToParameter = aToIntAna.ParamOnConic(1);
355     Standard_Real aFromParameter = aFromIntAna.ParamOnConic(1);
356     if(aToParameter > aFromParameter) {
357       gp_Vec aVec = aToDir->impl<gp_Dir>();
358       if((aVec * anExtVec) > 0) {
359         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
360         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
361       }
362       aVec = aFromDir->impl<gp_Dir>();
363       if((aVec * anExtVec) < 0) {
364         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
365         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
366       }
367     } else {
368       gp_Vec aVec = aToDir->impl<gp_Dir>();
369       if((aVec * anExtVec) < 0) {
370         aToDir->setImpl(new gp_Dir(aVec.Reversed()));
371         aBoundingToShape = GeomAlgoAPI_FaceBuilder::planarFace(aToPnt, aToDir);
372       }
373       aVec = aFromDir->impl<gp_Dir>();
374       if((aVec * anExtVec) > 0) {
375         aFromDir->setImpl(new gp_Dir(aVec.Reversed()));
376         aBoundingFromShape = GeomAlgoAPI_FaceBuilder::planarFace(aFromPnt, aFromDir);
377       }
378     }
379
380     // Making solids from bounding planes.
381     TopoDS_Shell aToShell, aFromShell;
382     TopoDS_Solid aToSolid, aFromSolid;
383     const TopoDS_Shape& aToShape   = aBoundingToShape->impl<TopoDS_Shape>();
384     const TopoDS_Shape& aFromShape = aBoundingFromShape->impl<TopoDS_Shape>();
385     TopoDS_Face aToFace   = TopoDS::Face(aToShape);
386     TopoDS_Face aFromFace = TopoDS::Face(aFromShape);
387     BRep_Builder aBoundingBuilder;
388     aBoundingBuilder.MakeShell(aToShell);
389     aBoundingBuilder.Add(aToShell, aToShape);
390     aBoundingBuilder.MakeShell(aFromShell);
391     aBoundingBuilder.Add(aFromShell, aFromShape);
392     aBoundingBuilder.MakeSolid(aToSolid);
393     aBoundingBuilder.Add(aToSolid, aToShell);
394     aBoundingBuilder.MakeSolid(aFromSolid);
395     aBoundingBuilder.Add(aFromSolid, aFromShell);
396
397     // Cutting with to plane.
398     BRepAlgoAPI_Cut* aToCutBuilder = new BRepAlgoAPI_Cut(aResult, aToSolid);
399     aToCutBuilder->Build();
400     if(!aToCutBuilder->IsDone()) {
401       return;
402     }
403     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
404       new GeomAlgoAPI_MakeShape(aToCutBuilder)));
405     aResult = aToCutBuilder->Shape();
406     if(aResult.ShapeType() == TopAbs_COMPOUND) {
407       aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
408     }
409     if(aShapeTypeToExp == TopAbs_FACE || aShapeTypeToExp == TopAbs_COMPOUND) {
410       const TopTools_ListOfShape& aToShapes = aToCutBuilder->Modified(aToShape);
411       for(TopTools_ListIteratorOfListOfShape anIt(aToShapes); anIt.More(); anIt.Next()) {
412         GeomShapePtr aGeomSh(new GeomAPI_Shape());
413         aGeomSh->setImpl(new TopoDS_Shape(anIt.Value()));
414         this->addToShape(aGeomSh);
415       }
416     }
417
418     // Cutting with from plane.
419     BRepAlgoAPI_Cut* aFromCutBuilder = new BRepAlgoAPI_Cut(aResult, aFromSolid);
420     aFromCutBuilder->Build();
421     if(!aFromCutBuilder->IsDone()) {
422       return;
423     }
424     this->appendAlgo(std::shared_ptr<GeomAlgoAPI_MakeShape>(
425       new GeomAlgoAPI_MakeShape(aFromCutBuilder)));
426     aResult = aFromCutBuilder->Shape();
427     TopoDS_Iterator aCheckIt(aResult);
428     if(!aCheckIt.More()) {
429       return;
430     }
431     if(aResult.ShapeType() == TopAbs_COMPOUND) {
432       aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
433     }
434     if(aShapeTypeToExp == TopAbs_FACE || aShapeTypeToExp == TopAbs_COMPOUND) {
435       const TopTools_ListOfShape& aFromShapes = aFromCutBuilder->Modified(aFromShape);
436       for(TopTools_ListIteratorOfListOfShape anIt(aFromShapes); anIt.More(); anIt.Next()) {
437         GeomShapePtr aGeomSh(new GeomAPI_Shape());
438         aGeomSh->setImpl(new TopoDS_Shape(anIt.Value()));
439         this->addFromShape(aGeomSh);
440       }
441     }
442
443     // Naming for extrusion from vertex, edge.
444     if(aShapeTypeToExp == TopAbs_COMPOUND) {
445       storeGenerationHistory(this, aResult, TopAbs_EDGE, aToFace, aFromFace);
446       storeGenerationHistory(this, aResult, TopAbs_FACE, aToFace, aFromFace);
447     } else {
448       storeGenerationHistory(this, aResult, aShapeTypeToExp, aToFace, aFromFace);
449     }
450
451     if(aResult.ShapeType() == TopAbs_COMPOUND) {
452       std::shared_ptr<GeomAPI_Shape> aGeomShape(new GeomAPI_Shape);
453       aGeomShape->setImpl(new TopoDS_Shape(aResult));
454       ListOfShape aCompSolids, aFreeSolids;
455       aGeomShape = GeomAlgoAPI_ShapeTools::combineShapes(aGeomShape,
456                                                          GeomAPI_Shape::COMPSOLID,
457                                                          aCompSolids,
458                                                          aFreeSolids);
459       aResult = aGeomShape->impl<TopoDS_Shape>();
460     }
461   }
462
463   // Setting result.
464   if(aResult.IsNull()) {
465     return;
466   }
467   aResult = GeomAlgoAPI_DFLoader::refineResult(aResult);
468   GeomShapePtr aGeomSh(new GeomAPI_Shape());
469   aGeomSh->setImpl(new TopoDS_Shape(aResult));
470   this->setShape(aGeomSh);
471   this->setDone(true);
472 }
473
474 // Auxilary functions:
475 //==================================================================================================
476 void storeGenerationHistory(GeomAlgoAPI_Prism* thePrismAlgo,
477                             const TopoDS_Shape& theBase,
478                             const TopAbs_ShapeEnum theType,
479                             BRepPrimAPI_MakePrism* thePrismBuilder)
480 {
481   for(TopExp_Explorer anExp(theBase, theType); anExp.More(); anExp.Next()) {
482     const TopoDS_Shape& aShape = anExp.Current();
483     GeomShapePtr aFromShape(new GeomAPI_Shape), aToShape(new GeomAPI_Shape);
484     aFromShape->setImpl(new TopoDS_Shape(thePrismBuilder->FirstShape(aShape)));
485     aToShape->setImpl(new TopoDS_Shape(thePrismBuilder->LastShape(aShape)));
486     thePrismAlgo->addFromShape(aFromShape);
487     thePrismAlgo->addToShape(aToShape);
488   }
489 }
490
491 //==================================================================================================
492 void storeGenerationHistory(GeomAlgoAPI_Prism* thePrismAlgo,
493                             const TopoDS_Shape& theResult,
494                             const TopAbs_ShapeEnum theType,
495                             const TopoDS_Face& theToFace,
496                             const TopoDS_Face& theFromFace)
497 {
498   for(TopExp_Explorer anExp(theResult, theType); anExp.More(); anExp.Next()) {
499     const TopoDS_Shape& aShape = anExp.Current();
500     GeomShapePtr aGeomSh(new GeomAPI_Shape());
501     if(theType == TopAbs_VERTEX) {
502       gp_Pnt aPnt = BRep_Tool::Pnt(TopoDS::Vertex(aShape));
503       IntTools_Context anIntTools;
504       if(anIntTools.IsValidPointForFace(aPnt,
505           theToFace, Precision::Confusion()) == Standard_True) {
506         aGeomSh->setImpl(new TopoDS_Shape(aShape));
507         thePrismAlgo->addToShape(aGeomSh);
508       }
509       if(anIntTools.IsValidPointForFace(aPnt,
510           theFromFace, Precision::Confusion()) == Standard_True) {
511         aGeomSh->setImpl(new TopoDS_Shape(aShape));
512         thePrismAlgo->addFromShape(aGeomSh);
513       }
514     } else if(theType == TopAbs_EDGE) {
515       TopoDS_Edge anEdge = TopoDS::Edge(aShape);
516       BRepLib_CheckCurveOnSurface anEdgeCheck(anEdge, theToFace);
517       anEdgeCheck.Perform();
518       if(anEdgeCheck.MaxDistance() < Precision::Confusion()) {
519         aGeomSh->setImpl(new TopoDS_Shape(aShape));
520         thePrismAlgo->addToShape(aGeomSh);
521       }
522       anEdgeCheck.Init(anEdge, theFromFace);
523       anEdgeCheck.Perform();
524       if(anEdgeCheck.MaxDistance() < Precision::Confusion()) {
525         aGeomSh->setImpl(new TopoDS_Shape(aShape));
526         thePrismAlgo->addFromShape(aGeomSh);
527       }
528     } else {
529       break;
530     }
531   }
532 }