]> SALOME platform Git repositories - modules/shaper.git/blob - src/GeomAlgoAPI/GeomAlgoAPI_ShapeTools.cpp
Salome HOME
Porting Salome to OCCT 7.7.0
[modules/shaper.git] / src / GeomAlgoAPI / GeomAlgoAPI_ShapeTools.cpp
1 // Copyright (C) 2014-2022  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 email : webmaster.salome@opencascade.com
18 //
19
20 #include "GeomAlgoAPI_ShapeTools.h"
21
22 #include "GeomAlgoAPI_SketchBuilder.h"
23
24 #include <Basics_OCCTVersion.hxx>
25
26 #include <GeomAPI_Ax1.h>
27 #include <GeomAPI_Edge.h>
28 #include <GeomAPI_Dir.h>
29 #include <GeomAPI_Face.h>
30 #include <GeomAPI_Pln.h>
31 #include <GeomAPI_Pnt.h>
32 #include <GeomAPI_Wire.h>
33
34 #include <Approx_CurvilinearParameter.hxx>
35
36 #include <Bnd_Box.hxx>
37
38 #include <BRep_Tool.hxx>
39 #include <BRep_Builder.hxx>
40 #include <BRepAlgo.hxx>
41 #include <BRepAlgo_FaceRestrictor.hxx>
42 #include <BRepAdaptor_Curve.hxx>
43 #include <BRepBndLib.hxx>
44 #include <BRepBuilderAPI_Copy.hxx>
45 #include <BRepBuilderAPI_FindPlane.hxx>
46 #include <BRepBuilderAPI_MakeEdge.hxx>
47 #include <BRepBuilderAPI_MakeFace.hxx>
48 #include <BRepBuilderAPI_MakeVertex.hxx>
49 #include <BRepBuilderAPI_MakeWire.hxx>
50 #include <BRepCheck_Analyzer.hxx>
51 #include <BRepExtrema_DistShapeShape.hxx>
52 #include <BRepExtrema_ExtCF.hxx>
53 #include <BRepGProp.hxx>
54 #include <BRepTools.hxx>
55 #include <BRepTools_WireExplorer.hxx>
56 #include <BRepTopAdaptor_FClass2d.hxx>
57 #include <BRepClass_FaceClassifier.hxx>
58 #include <BRepLib_CheckCurveOnSurface.hxx>
59 #include <BRepLProp.hxx>
60
61 #include <BOPAlgo_Builder.hxx>
62
63 #include <Geom2d_Curve.hxx>
64 #include <Geom2d_Curve.hxx>
65
66 #include <Geom_BSplineCurve.hxx>
67 #include <Geom_CylindricalSurface.hxx>
68 #include <Geom_Line.hxx>
69 #include <Geom_Plane.hxx>
70 #include <Geom_RectangularTrimmedSurface.hxx>
71
72 #if OCC_VERSION_LARGE < 0x07070000
73 #include <GeomAdaptor_HCurve.hxx>
74 #else
75 #include <GeomAdaptor_Curve.hxx>
76 #endif
77
78 #include <GeomAPI_ProjectPointOnCurve.hxx>
79 #include <GeomAPI_ShapeIterator.h>
80
81 #include <GeomLib_IsPlanarSurface.hxx>
82 #include <GeomLib_Tool.hxx>
83 #include <GeomAPI_IntCS.hxx>
84
85 #include <gp_Pln.hxx>
86 #include <GProp_GProps.hxx>
87
88 #include <IntAna_IntConicQuad.hxx>
89 #include <IntAna_Quadric.hxx>
90
91 #include <ShapeAnalysis.hxx>
92 #include <ShapeAnalysis_Surface.hxx>
93
94 #include <TopoDS.hxx>
95 #include <TopoDS_Edge.hxx>
96 #include <TopoDS_Face.hxx>
97 #include <TopoDS_Shape.hxx>
98 #include <TopoDS_Shell.hxx>
99 #include <TopoDS_Vertex.hxx>
100 #include <TopoDS_Builder.hxx>
101
102 #include <TopExp.hxx>
103 #include <TopExp_Explorer.hxx>
104
105 #include <TopTools_ListIteratorOfListOfShape.hxx>
106
107 #include <NCollection_Vector.hxx>
108
109 #include <LocalAnalysis_SurfaceContinuity.hxx>
110 #include<array>
111
112 //==================================================================================================
113 static GProp_GProps props(const TopoDS_Shape& theShape)
114 {
115   GProp_GProps aGProps;
116
117   if (theShape.ShapeType() == TopAbs_EDGE || theShape.ShapeType() == TopAbs_WIRE)
118   {
119     BRepGProp::LinearProperties(theShape, aGProps);
120   }
121   else if (theShape.ShapeType() == TopAbs_FACE || theShape.ShapeType() == TopAbs_SHELL)
122   {
123     const Standard_Real anEps = 1.e-6;
124     BRepGProp::SurfaceProperties(theShape, aGProps, anEps);
125   }
126   else if (theShape.ShapeType() == TopAbs_SOLID || theShape.ShapeType() == TopAbs_COMPSOLID)
127   {
128     BRepGProp::VolumeProperties(theShape, aGProps);
129   }
130   else if (theShape.ShapeType() == TopAbs_COMPOUND)
131   {
132     for (TopoDS_Iterator anIt(theShape); anIt.More(); anIt.Next())
133     {
134       aGProps.Add(props(anIt.Value()));
135     }
136   }
137
138   return aGProps;
139 }
140
141 //==================================================================================================
142 double GeomAlgoAPI_ShapeTools::length(const std::shared_ptr<GeomAPI_Shape> theShape)
143 {
144   GProp_GProps aGProps;
145   if (!theShape.get()) {
146     return 0.0;
147   }
148   const TopoDS_Shape& aShape = theShape->impl<TopoDS_Shape>();
149   if (aShape.IsNull()) {
150     return 0.0;
151   }
152
153   BRepGProp::LinearProperties(aShape, aGProps, Standard_True);
154   return  aGProps.Mass();
155 }
156
157 //==================================================================================================
158 double GeomAlgoAPI_ShapeTools::volume(const std::shared_ptr<GeomAPI_Shape> theShape)
159 {
160   if (!theShape.get()) {
161     return 0.0;
162   }
163   const TopoDS_Shape& aShape = theShape->impl<TopoDS_Shape>();
164   if (aShape.IsNull()) {
165     return 0.0;
166   }
167   const Standard_Real anEps = 1.e-6;
168   double aVolume = 0.0;
169   for (TopExp_Explorer anExp(aShape, TopAbs_SOLID); anExp.More(); anExp.Next()) {
170     GProp_GProps aGProps;
171     BRepGProp::VolumeProperties(anExp.Current(), aGProps, anEps);
172     aVolume += aGProps.Mass();
173   }
174   return aVolume;
175 }
176
177 //==================================================================================================
178 double GeomAlgoAPI_ShapeTools::area (const std::shared_ptr<GeomAPI_Shape> theShape)
179 {
180   GProp_GProps aGProps;
181   if (!theShape.get()) {
182     return 0.0;
183   }
184   const TopoDS_Shape& aShape = theShape->impl<TopoDS_Shape>();
185   if (aShape.IsNull()) {
186     return 0.0;
187   }
188   const Standard_Real anEps = 1.e-6;
189
190   BRepGProp::SurfaceProperties(aShape, aGProps, anEps);
191   return aGProps.Mass();
192 }
193
194 //==================================================================================================
195 bool GeomAlgoAPI_ShapeTools::isContinuousFaces(const GeomShapePtr& theFace1,
196                                                const GeomShapePtr& theFace2,
197                                                const GeomPointPtr& thePoint,
198                                                const double & theAngle,
199                                                std::string& theError)
200 {
201
202   #ifdef _DEBUG
203   std::cout << "isContinuousFaces " << std::endl;
204   #endif
205
206   if (!thePoint.get()) {
207       theError = "isContinuousFaces : An invalid argument";
208       return false;
209   }
210   const gp_Pnt& aPoint = thePoint->impl<gp_Pnt>();
211
212   // Getting base shape.
213   if (!theFace1.get()) {
214     theError = "isContinuousFaces : An invalid argument";
215     return false;
216   }
217
218   TopoDS_Shape aShape1 = theFace1->impl<TopoDS_Shape>();
219
220   if (aShape1.IsNull()) {
221     theError = "isContinuousFaces : An invalid argument";
222     return false;
223   }
224
225   // Getting base shape.
226   if (!theFace2.get()) {
227     theError = "isContinuousFaces : An invalid argument";
228     return false;
229   }
230
231   TopoDS_Shape aShape2 = theFace2->impl<TopoDS_Shape>();
232
233   if (aShape2.IsNull()) {
234     theError = "isContinuousFaces : An invalid argument";
235     return false;
236   }
237
238   TopoDS_Face aFace1 = TopoDS::Face(aShape1);
239   if (aFace1.IsNull()) {
240     theError = "isContinuousFaces : An invalid argument";
241     return false;
242   }
243
244   Handle(Geom_Surface) aSurf1 = BRep_Tool::Surface(aFace1);
245   if (aSurf1.IsNull()) {
246     theError = "isContinuousFaces : An invalid surface";
247     return false;
248   }
249
250   ShapeAnalysis_Surface aSAS1(aSurf1);
251   gp_Pnt2d aPointOnFace1 = aSAS1.ValueOfUV(aPoint, Precision::Confusion());
252
253   TopoDS_Face aFace2 = TopoDS::Face(aShape2);
254   if (aFace2.IsNull()) {
255     theError = "isContinuousFaces : An invalid argument";
256     return false;
257   }
258
259   Handle(Geom_Surface) aSurf2 = BRep_Tool::Surface(aFace2);
260   if (aSurf2.IsNull()) {
261     theError = "isContinuousFaces : An invalid surface";
262     return false;
263   }
264
265   ShapeAnalysis_Surface aSAS2(aSurf2);
266   gp_Pnt2d aPointOnFace2= aSAS2.ValueOfUV(aPoint, Precision::Confusion());
267
268   bool aRes = false;
269   try {
270     OCC_CATCH_SIGNALS;
271     LocalAnalysis_SurfaceContinuity aLocAnal(aSurf1,
272                                              aPointOnFace1.X(),
273                                              aPointOnFace1.Y(),
274                                              aSurf2,
275                                              aPointOnFace2.X(),
276                                              aPointOnFace2.Y(),
277                                              GeomAbs_Shape::GeomAbs_G1, // Order
278                                              0.001, // EpsNul
279                                              0.001, // EpsC0
280                                              0.001, // EpsC1
281                                              0.001, // EpsC2
282                                              theAngle * M_PI / 180.0); //EpsG1
283     aRes = aLocAnal.IsG1();
284   }
285   catch (Standard_Failure const& anException) {
286     theError = "LocalAnalysis_SurfaceContinuity error :";
287     theError += anException.GetMessageString();
288   }
289
290   return aRes;
291 }
292
293 //==================================================================================================
294 std::shared_ptr<GeomAPI_Pnt>
295   GeomAlgoAPI_ShapeTools::centreOfMass(const std::shared_ptr<GeomAPI_Shape> theShape)
296 {
297   GProp_GProps aGProps;
298   if (!theShape) {
299     return std::shared_ptr<GeomAPI_Pnt>();
300   }
301   const TopoDS_Shape& aShape = theShape->impl<TopoDS_Shape>();
302   if (aShape.IsNull()) {
303     return std::shared_ptr<GeomAPI_Pnt>();
304   }
305   gp_Pnt aCentre;
306   if (aShape.ShapeType() == TopAbs_VERTEX) {
307     aCentre = BRep_Tool::Pnt(TopoDS::Vertex(aShape));
308   } else {
309     aGProps = props(aShape);
310     aCentre = aGProps.CentreOfMass();
311   }
312
313   return std::shared_ptr<GeomAPI_Pnt>(new GeomAPI_Pnt(aCentre.X(), aCentre.Y(), aCentre.Z()));
314 }
315
316 //==================================================================================================
317 double GeomAlgoAPI_ShapeTools::radius(const std::shared_ptr<GeomAPI_Face>& theCylinder)
318 {
319   double aRadius = -1.0;
320   if (theCylinder->isCylindrical()) {
321     const TopoDS_Shape& aShape = theCylinder->impl<TopoDS_Shape>();
322     Handle(Geom_Surface) aSurf = BRep_Tool::Surface(TopoDS::Face(aShape));
323     Handle(Geom_CylindricalSurface) aCyl = Handle(Geom_CylindricalSurface)::DownCast(aSurf);
324     if (!aCyl.IsNull())
325       aRadius = aCyl->Radius();
326   }
327   return aRadius;
328 }
329
330 //==================================================================================================
331 namespace {
332
333 auto getExtemaDistShape = [](const GeomShapePtr& theShape1,
334     const GeomShapePtr& theShape2) -> BRepExtrema_DistShapeShape
335 {
336   const TopoDS_Shape& aShape1 = theShape1->impl<TopoDS_Shape>();
337   const TopoDS_Shape& aShape2 = theShape2->impl<TopoDS_Shape>();
338
339   BRepExtrema_DistShapeShape aDist(aShape1, aShape2);
340   aDist.Perform();
341   return aDist;
342 };
343 }
344
345 double GeomAlgoAPI_ShapeTools::minimalDistance(const GeomShapePtr& theShape1,
346                                                const GeomShapePtr& theShape2)
347 {
348   BRepExtrema_DistShapeShape aDist = getExtemaDistShape(theShape1, theShape2);
349   return aDist.IsDone() ? aDist.Value() : Precision::Infinite();
350 }
351
352 double GeomAlgoAPI_ShapeTools::minimalDistance(const GeomShapePtr& theShape1,
353                                                const GeomShapePtr& theShape2,
354                                                std::array<double, 3> & fromShape1To2)
355 {
356   BRepExtrema_DistShapeShape aDist = getExtemaDistShape(theShape1, theShape2);
357   const auto & pt1 = aDist.PointOnShape1(1);
358   const auto & pt2 = aDist.PointOnShape2(1) ;
359   fromShape1To2[0] = pt2.X() - pt1.X();
360   fromShape1To2[1] = pt2.Y() - pt1.Y();
361   fromShape1To2[2] = pt2.Z() - pt1.Z();
362   return aDist.IsDone() ? aDist.Value() : Precision::Infinite();
363 }
364
365 //==================================================================================================
366 std::shared_ptr<GeomAPI_Shape> GeomAlgoAPI_ShapeTools::combineShapes(
367   const std::shared_ptr<GeomAPI_Shape> theCompound,
368   const GeomAPI_Shape::ShapeType theType,
369   ListOfShape& theResuts)
370 {
371
372   ListOfShape aResCombinedShapes;
373   ListOfShape aResFreeShapes;
374
375   GeomShapePtr aResult = theCompound;
376
377   if (!theCompound.get()) {
378     return aResult;
379   }
380
381   if (theType != GeomAPI_Shape::SHELL && theType != GeomAPI_Shape::COMPSOLID) {
382     return aResult;
383   }
384
385   TopAbs_ShapeEnum aTS = TopAbs_EDGE;
386   TopAbs_ShapeEnum aTA = TopAbs_FACE;
387   if (theType == GeomAPI_Shape::COMPSOLID) {
388     aTS = TopAbs_FACE;
389     aTA = TopAbs_SOLID;
390   }
391
392   // map from the resulting shapes to minimal index of the used shape from theCompound list
393   std::map<GeomShapePtr, int> anInputOrder;
394   // map from ancestors-shapes to the index of shapes in theCompound
395   NCollection_DataMap<TopoDS_Shape, int> anAncestorsOrder;
396
397   // Get free shapes.
398   int anOrder = 0;
399   const TopoDS_Shape& aShapesComp = theCompound->impl<TopoDS_Shape>();
400   for (TopoDS_Iterator anIter(aShapesComp); anIter.More(); anIter.Next(), anOrder++) {
401     const TopoDS_Shape& aShape = anIter.Value();
402     if (aShape.ShapeType() > aTA) {
403       std::shared_ptr<GeomAPI_Shape> aGeomShape(new GeomAPI_Shape);
404       aGeomShape->setImpl<TopoDS_Shape>(new TopoDS_Shape(aShape));
405       aResFreeShapes.push_back(aGeomShape);
406       anInputOrder[aGeomShape] = anOrder;
407     } else {
408       for (TopExp_Explorer anExp(aShape, aTA); anExp.More(); anExp.Next()) {
409         anAncestorsOrder.Bind(anExp.Current(), anOrder);
410       }
411     }
412   }
413
414   // Map sub-shapes and shapes.
415   TopTools_IndexedDataMapOfShapeListOfShape aMapSA;
416   TopExp::MapShapesAndAncestors(aShapesComp, aTS, aTA, aMapSA);
417   if (aMapSA.IsEmpty()) {
418     return aResult;
419   }
420   theResuts.clear();
421
422   // Get all shapes with common sub-shapes and free shapes.
423   NCollection_Map<TopoDS_Shape> aFreeShapes;
424   NCollection_Vector<NCollection_Map<TopoDS_Shape>> aShapesWithCommonSubshapes;
425   for (TopTools_IndexedDataMapOfShapeListOfShape::Iterator
426       anIter(aMapSA); anIter.More(); anIter.Next()) {
427     TopTools_ListOfShape& aListOfShape = anIter.ChangeValue();
428     if (aListOfShape.IsEmpty()) {
429       continue;
430     }
431     else if (aListOfShape.Size() == 1) {
432       const TopoDS_Shape& aF = aListOfShape.First();
433       aFreeShapes.Add(aF);
434       aListOfShape.Clear();
435     } else {
436       NCollection_List<TopoDS_Shape> aTempList;
437       NCollection_Map<TopoDS_Shape> aTempMap;
438       for (TopTools_ListOfShape::Iterator aListIt(aListOfShape); aListIt.More(); aListIt.Next()) {
439         aTempList.Append(aListIt.Value());
440         aTempMap.Add(aListIt.Value());
441         aFreeShapes.Remove(aListIt.Value());
442       }
443       aListOfShape.Clear();
444       for (NCollection_List<TopoDS_Shape>::Iterator
445           aTempIter(aTempList); aTempIter.More(); aTempIter.Next()) {
446         const TopoDS_Shape& aTempShape = aTempIter.Value();
447         for (TopTools_IndexedDataMapOfShapeListOfShape::Iterator
448             anIter2(aMapSA); anIter2.More(); anIter2.Next()) {
449           TopTools_ListOfShape& aTempListOfShape = anIter2.ChangeValue();
450           if (aTempListOfShape.IsEmpty()) {
451             continue;
452           } else if (aTempListOfShape.Size() == 1 && aTempListOfShape.First() == aTempShape) {
453             aTempListOfShape.Clear();
454           } else if (aTempListOfShape.Size() > 1) {
455             TopTools_ListOfShape::Iterator anIt1(aTempListOfShape);
456             for (; anIt1.More(); anIt1.Next()) {
457               if (anIt1.Value() == aTempShape) {
458                 TopTools_ListOfShape::Iterator anIt2(aTempListOfShape);
459                 for (; anIt2.More(); anIt2.Next())
460                 {
461                   if (anIt2.Value() != anIt1.Value()) {
462                     if (aTempMap.Add(anIt2.Value())) {
463                       aTempList.Append(anIt2.Value());
464                       aFreeShapes.Remove(anIt2.Value());
465                     }
466                   }
467                 }
468                 aTempListOfShape.Clear();
469                 break;
470               }
471             }
472           }
473         }
474       }
475       aShapesWithCommonSubshapes.Append(aTempMap);
476     }
477   }
478
479   // Combine shapes with common sub-shapes.
480   for (NCollection_Vector<NCollection_Map<TopoDS_Shape>>::Iterator
481       anIter(aShapesWithCommonSubshapes); anIter.More(); anIter.Next()) {
482     TopoDS_Shell aShell;
483     TopoDS_CompSolid aCSolid;
484     TopoDS_Builder aBuilder;
485     anOrder = -1;
486     theType ==
487       GeomAPI_Shape::COMPSOLID ? aBuilder.MakeCompSolid(aCSolid) : aBuilder.MakeShell(aShell);
488     NCollection_Map<TopoDS_Shape>& aShapesMap = anIter.ChangeValue();
489     for (TopExp_Explorer anExp(aShapesComp, aTA); anExp.More(); anExp.Next()) {
490       const TopoDS_Shape& aShape = anExp.Current();
491       if (aShapesMap.Contains(aShape)) {
492         theType ==
493           GeomAPI_Shape::COMPSOLID ? aBuilder.Add(aCSolid, aShape) : aBuilder.Add(aShell, aShape);
494         aShapesMap.Remove(aShape);
495         int aThisOrder = anAncestorsOrder.Find(aShape);
496         if (anOrder == -1 || aThisOrder < anOrder)
497           anOrder = aThisOrder; // take the minimum order position
498       }
499     }
500     std::shared_ptr<GeomAPI_Shape> aGeomShape(new GeomAPI_Shape);
501     TopoDS_Shape* aSh = theType == GeomAPI_Shape::COMPSOLID ? new TopoDS_Shape(aCSolid) :
502                                                               new TopoDS_Shape(aShell);
503     aGeomShape->setImpl<TopoDS_Shape>(aSh);
504     aResCombinedShapes.push_back(aGeomShape);
505     anInputOrder[aGeomShape] = anOrder;
506   }
507
508   // Adding free shapes.
509   for (TopExp_Explorer anExp(aShapesComp, aTA); anExp.More(); anExp.Next()) {
510     const TopoDS_Shape& aShape = anExp.Current();
511     if (aFreeShapes.Contains(aShape)) {
512       std::shared_ptr<GeomAPI_Shape> aGeomShape(new GeomAPI_Shape);
513       aGeomShape->setImpl<TopoDS_Shape>(new TopoDS_Shape(aShape));
514       aResFreeShapes.push_back(aGeomShape);
515       anInputOrder[aGeomShape] = anAncestorsOrder.Find(aShape);
516     }
517   }
518
519   if (aResCombinedShapes.size() == 1 && aResFreeShapes.size() == 0) {
520     aResult = aResCombinedShapes.front();
521     theResuts.push_back(aResult);
522   } else if (aResCombinedShapes.size() == 0 && aResFreeShapes.size() == 1) {
523     aResult = aResFreeShapes.front();
524     theResuts.push_back(aResult);
525   } else {
526     TopoDS_Compound aResultComp;
527     TopoDS_Builder aBuilder;
528     aBuilder.MakeCompound(aResultComp);
529     // put to result compound and result list in accordance to the order numbers
530     std::map<GeomShapePtr, int>::iterator anInputIter = anInputOrder.begin();
531     std::map<int, GeomShapePtr> aNums;
532     for (; anInputIter != anInputOrder.end(); anInputIter++)
533       aNums[anInputIter->second] = anInputIter->first;
534     std::map<int, GeomShapePtr>::iterator aNumsIter = aNums.begin();
535     for (; aNumsIter != aNums.end(); aNumsIter++) {
536       aBuilder.Add(aResultComp, (aNumsIter->second)->impl<TopoDS_Shape>());
537       theResuts.push_back(aNumsIter->second);
538     }
539     aResult->setImpl(new TopoDS_Shape(aResultComp));
540   }
541
542   return aResult;
543 }
544
545 //==================================================================================================
546 static void addSimpleShapeToList(const TopoDS_Shape& theShape,
547                                  NCollection_List<TopoDS_Shape>& theList)
548 {
549   if (theShape.IsNull()) {
550     return;
551   }
552
553   if (theShape.ShapeType() == TopAbs_COMPOUND) {
554     for (TopoDS_Iterator anIt(theShape); anIt.More(); anIt.Next()) {
555       addSimpleShapeToList(anIt.Value(), theList);
556     }
557   } else {
558     theList.Append(theShape);
559   }
560 }
561
562 //==================================================================================================
563 static TopoDS_Compound makeCompound(const NCollection_List<TopoDS_Shape> theShapes)
564 {
565   TopoDS_Compound aCompound;
566
567   BRep_Builder aBuilder;
568   aBuilder.MakeCompound(aCompound);
569
570   for (NCollection_List<TopoDS_Shape>::Iterator anIt(theShapes); anIt.More(); anIt.Next()) {
571     aBuilder.Add(aCompound, anIt.Value());
572   }
573
574   return aCompound;
575 }
576
577 //==================================================================================================
578 std::shared_ptr<GeomAPI_Shape> GeomAlgoAPI_ShapeTools::groupSharedTopology(
579   const std::shared_ptr<GeomAPI_Shape> theCompound)
580 {
581   GeomShapePtr aResult = theCompound;
582
583   if (!theCompound.get()) {
584     return aResult;
585   }
586
587   TopoDS_Shape anInShape = aResult->impl<TopoDS_Shape>();
588   NCollection_List<TopoDS_Shape> anUngroupedShapes, aStillUngroupedShapes;
589   addSimpleShapeToList(anInShape, anUngroupedShapes);
590
591   // Iterate over all shapes and find shapes with shared vertices.
592   TopTools_ListOfShape allVertices;
593   TopTools_DataMapOfShapeListOfShape aVertexShapesMap;
594   for (NCollection_List<TopoDS_Shape>::Iterator aShapesIt(anUngroupedShapes);
595     aShapesIt.More();
596     aShapesIt.Next()) {
597     const TopoDS_Shape& aShape = aShapesIt.Value();
598     for (TopExp_Explorer aShapeExp(aShape, TopAbs_VERTEX);
599       aShapeExp.More();
600       aShapeExp.Next()) {
601       const TopoDS_Shape& aVertex = aShapeExp.Current();
602       if (!aVertexShapesMap.IsBound(aVertex)) {
603         NCollection_List<TopoDS_Shape> aList;
604         aList.Append(aShape);
605         allVertices.Append(aVertex);
606         aVertexShapesMap.Bind(aVertex, aList);
607       }
608       else {
609         if (!aVertexShapesMap.ChangeFind(aVertex).Contains(aShape)) {
610           aVertexShapesMap.ChangeFind(aVertex).Append(aShape);
611         }
612       }
613     }
614   }
615
616   // Iterate over the map and group shapes.
617   NCollection_Vector<TopTools_MapOfShape> aGroups; // groups of shapes connected by vertices
618   while (!allVertices.IsEmpty()) {
619     // Get first group of shapes in map, and then unbind it.
620     const TopoDS_Shape& aKey = allVertices.First();
621     TopTools_ListOfShape aConnectedShapes = aVertexShapesMap.Find(aKey);
622     aVertexShapesMap.UnBind(aKey);
623     allVertices.Remove(aKey);
624     // Iterate over shapes in this group and add to it shapes from groups in map.
625     for (TopTools_ListOfShape::Iterator aConnectedIt(aConnectedShapes);
626       aConnectedIt.More(); aConnectedIt.Next()) {
627       const TopoDS_Shape& aConnected = aConnectedIt.Value();
628       TopTools_ListOfShape aKeysToUnbind;
629       for (TopTools_ListOfShape::Iterator aKeysIt(allVertices);
630         aKeysIt.More();
631         aKeysIt.Next()) {
632         const TopTools_ListOfShape& anOtherConnected = aVertexShapesMap(aKeysIt.Value());
633         if (!anOtherConnected.Contains(aConnected)) {
634           // Other connected group does not contain shape from our connected group
635           continue;
636         }
637         // Other is connected to our, so add them to our connected
638         for (TopTools_ListOfShape::Iterator anOtherIt(anOtherConnected);
639           anOtherIt.More();
640           anOtherIt.Next()) {
641           const TopoDS_Shape& aShape = anOtherIt.Value();
642           if (!aConnectedShapes.Contains(aShape)) {
643             aConnectedShapes.Append(aShape);
644           }
645         }
646         // Save key to unbind from this map.
647         aKeysToUnbind.Append(aKeysIt.Value());
648       }
649       // Unbind groups from map that we added to our group.
650       for (TopTools_ListOfShape::Iterator aKeysIt(aKeysToUnbind);
651         aKeysIt.More();
652         aKeysIt.Next()) {
653         aVertexShapesMap.UnBind(aKeysIt.Value());
654         allVertices.Remove(aKeysIt.Value());
655       }
656     }
657     // Sort shapes from the most complicated to the simplest ones
658     TopTools_MapOfShape aSortedGroup;
659     for (int aST = TopAbs_COMPOUND; aST <= TopAbs_SHAPE; ++aST) {
660       TopTools_ListOfShape::Iterator anIt(aConnectedShapes);
661       while (anIt.More()) {
662         if (anIt.Value().ShapeType() == aST) {
663           aSortedGroup.Add(anIt.Value());
664           aConnectedShapes.Remove(anIt);
665         }
666         else {
667           anIt.Next();
668         }
669       }
670     }
671     aGroups.Append(aSortedGroup);
672   }
673
674   TopoDS_Compound aCompound;
675   BRep_Builder aBuilder;
676   aBuilder.MakeCompound(aCompound);
677   ListOfShape aSolids;
678   for (NCollection_Vector<TopTools_MapOfShape>::Iterator anIt(aGroups); anIt.More(); anIt.Next()) {
679     const TopTools_MapOfShape& aGroup = anIt.ChangeValue();
680     GeomShapePtr aGeomShape(new GeomAPI_Shape());
681     if (aGroup.Size() == 1) {
682       TopTools_MapOfShape::Iterator aOneShapeIter(aGroup);
683       aGeomShape->setImpl(new TopoDS_Shape(aOneShapeIter.Value()));
684     } else {
685       // make sub-shapes in the group have order same as in original shape
686       TopTools_ListOfShape anOrderedGoup;
687       NCollection_List<TopoDS_Shape>::Iterator anUngrouped(anUngroupedShapes);
688       for (; anUngrouped.More(); anUngrouped.Next()) {
689         if (aGroup.Contains(anUngrouped.Value()))
690           anOrderedGoup.Append(anUngrouped.Value());
691       }
692       aGeomShape->setImpl(new TopoDS_Shape(makeCompound(anOrderedGoup)));
693       aGeomShape = GeomAlgoAPI_ShapeTools::combineShapes(aGeomShape,
694                                                          GeomAPI_Shape::COMPSOLID,
695                                                          aSolids);
696     }
697     aBuilder.Add(aCompound, aGeomShape->impl<TopoDS_Shape>());
698   }
699
700   if (!aCompound.IsNull()) {
701     aResult->setImpl(new TopoDS_Shape(aCompound));
702   }
703
704   return aResult;
705 }
706
707 //==================================================================================================
708 bool GeomAlgoAPI_ShapeTools::hasSharedTopology(const ListOfShape& theShapes,
709                                                const GeomAPI_Shape::ShapeType theShapeType)
710 {
711   TopTools_IndexedMapOfShape aSubs;
712   for (ListOfShape::const_iterator anIt = theShapes.begin(); anIt != theShapes.end(); ++anIt) {
713     TopTools_IndexedMapOfShape aCurSubs;
714     TopExp::MapShapes((*anIt)->impl<TopoDS_Shape>(), (TopAbs_ShapeEnum)theShapeType, aCurSubs);
715     for (TopTools_IndexedMapOfShape::Iterator aSubIt(aCurSubs); aSubIt.More(); aSubIt.Next()) {
716       if (aSubs.Contains(aSubIt.Value()))
717         return true;
718       else
719         aSubs.Add(aSubIt.Value());
720     }
721   }
722   return false;
723 }
724
725 //==================================================================================================
726 std::list<std::shared_ptr<GeomAPI_Pnt> >
727   GeomAlgoAPI_ShapeTools::getBoundingBox(const ListOfShape& theShapes, const double theEnlarge)
728 {
729   // Bounding box of all objects.
730   Bnd_Box aBndBox;
731
732   // Getting box.
733   for (ListOfShape::const_iterator
734     anObjectsIt = theShapes.begin(); anObjectsIt != theShapes.end(); anObjectsIt++) {
735     const TopoDS_Shape& aShape = (*anObjectsIt)->impl<TopoDS_Shape>();
736     BRepBndLib::Add(aShape, aBndBox);
737   }
738
739   if (theEnlarge != 0.0) {
740     // We enlarge bounding box just to be sure that plane will be large enough to cut all objects.
741     aBndBox.Enlarge(theEnlarge);
742   }
743
744   Standard_Real aXArr[2] = {aBndBox.CornerMin().X(), aBndBox.CornerMax().X()};
745   Standard_Real aYArr[2] = {aBndBox.CornerMin().Y(), aBndBox.CornerMax().Y()};
746   Standard_Real aZArr[2] = {aBndBox.CornerMin().Z(), aBndBox.CornerMax().Z()};
747   std::list<std::shared_ptr<GeomAPI_Pnt> > aResultPoints;
748   for (int i = 0; i < 2; i++) {
749     for (int j = 0; j < 2; j++) {
750       for (int k = 0; k < 2; k++) {
751         std::shared_ptr<GeomAPI_Pnt> aPnt(new GeomAPI_Pnt(aXArr[i], aYArr[j], aZArr[k]));
752         aResultPoints.push_back(aPnt);
753       }
754     }
755   }
756
757   return aResultPoints;
758 }
759
760 //==================================================================================================
761 std::shared_ptr<GeomAPI_Face> GeomAlgoAPI_ShapeTools::fitPlaneToBox(
762   const std::shared_ptr<GeomAPI_Shape> thePlane,
763   const std::list<std::shared_ptr<GeomAPI_Pnt> >& thePoints)
764 {
765   std::shared_ptr<GeomAPI_Face> aResultFace;
766
767   if (!thePlane.get()) {
768     return aResultFace;
769   }
770
771   const TopoDS_Shape& aShape = thePlane->impl<TopoDS_Shape>();
772   if (aShape.ShapeType() != TopAbs_FACE) {
773     return aResultFace;
774   }
775
776   TopoDS_Face aFace = TopoDS::Face(aShape);
777   Handle(Geom_Surface) aSurf = BRep_Tool::Surface(aFace);
778   if (aSurf.IsNull()) {
779     return aResultFace;
780   }
781
782   GeomLib_IsPlanarSurface isPlanar(aSurf);
783   if (!isPlanar.IsPlanar()) {
784     return aResultFace;
785   }
786
787   if (thePoints.size() != 8) {
788     return aResultFace;
789   }
790
791   const gp_Pln& aFacePln = isPlanar.Plan();
792   Handle(Geom_Plane) aFacePlane = new Geom_Plane(aFacePln);
793   IntAna_Quadric aQuadric(aFacePln);
794   Standard_Real UMin, UMax, VMin, VMax;
795   UMin = UMax = VMin = VMax = 0;
796   for (std::list<std::shared_ptr<GeomAPI_Pnt> >::const_iterator
797        aPointsIt = thePoints.begin(); aPointsIt != thePoints.end(); aPointsIt++) {
798     const gp_Pnt& aPnt = (*aPointsIt)->impl<gp_Pnt>();
799     gp_Lin aLin(aPnt, aFacePln.Axis().Direction());
800     IntAna_IntConicQuad anIntAna(aLin, aQuadric);
801     const gp_Pnt& aPntOnFace = anIntAna.Point(1);
802     Standard_Real aPntU(0), aPntV(0);
803     GeomLib_Tool::Parameters(aFacePlane, aPntOnFace, Precision::Confusion(), aPntU, aPntV);
804     if (aPntU < UMin) UMin = aPntU;
805     if (aPntU > UMax) UMax = aPntU;
806     if (aPntV < VMin) VMin = aPntV;
807     if (aPntV > VMax) VMax = aPntV;
808   }
809   aResultFace.reset(new GeomAPI_Face());
810   aResultFace->setImpl(new TopoDS_Face(BRepLib_MakeFace(aFacePln, UMin, UMax, VMin, VMax).Face()));
811
812   return aResultFace;
813 }
814
815 //==================================================================================================
816 void GeomAlgoAPI_ShapeTools::findBounds(const std::shared_ptr<GeomAPI_Shape> theShape,
817                                         std::shared_ptr<GeomAPI_Vertex>& theV1,
818                                         std::shared_ptr<GeomAPI_Vertex>& theV2)
819 {
820   static GeomVertexPtr aVertex;
821   if (!aVertex) {
822     aVertex = GeomVertexPtr(new GeomAPI_Vertex);
823     aVertex->setImpl(new TopoDS_Vertex());
824   }
825
826   theV1 = aVertex;
827   theV2 = aVertex;
828
829   if (theShape) {
830     const TopoDS_Shape& aShape = theShape->impl<TopoDS_Shape>();
831     TopoDS_Vertex aV1, aV2;
832     ShapeAnalysis::FindBounds(aShape, aV1, aV2);
833
834     std::shared_ptr<GeomAPI_Vertex> aGeomV1(new GeomAPI_Vertex()), aGeomV2(new GeomAPI_Vertex());
835     aGeomV1->setImpl(new TopoDS_Vertex(aV1));
836     aGeomV2->setImpl(new TopoDS_Vertex(aV2));
837     theV1 = aGeomV1;
838     theV2 = aGeomV2;
839   }
840 }
841
842 //==================================================================================================
843 void GeomAlgoAPI_ShapeTools::makeFacesWithHoles(const std::shared_ptr<GeomAPI_Pnt> theOrigin,
844                                                 const std::shared_ptr<GeomAPI_Dir> theDirection,
845                                                 const ListOfShape& theWires,
846                                                 ListOfShape& theFaces)
847 {
848   BRepBuilderAPI_MakeFace aMKFace(gp_Pln(theOrigin->impl<gp_Pnt>(),
849                                           theDirection->impl<gp_Dir>()));
850   TopoDS_Face aFace = aMKFace.Face();
851
852   BRepAlgo_FaceRestrictor aFRestrictor;
853   aFRestrictor.Init(aFace, Standard_False, Standard_True);
854   for (ListOfShape::const_iterator anIt = theWires.cbegin();
855       anIt != theWires.cend();
856       ++anIt) {
857     TopoDS_Wire aWire = TopoDS::Wire((*anIt)->impl<TopoDS_Shape>());
858     aFRestrictor.Add(aWire);
859   }
860
861   aFRestrictor.Perform();
862
863   if (!aFRestrictor.IsDone()) {
864     return;
865   }
866
867   for (; aFRestrictor.More(); aFRestrictor.Next()) {
868     GeomShapePtr aShape(new GeomAPI_Shape());
869     aShape->setImpl(new TopoDS_Shape(aFRestrictor.Current()));
870     theFaces.push_back(aShape);
871   }
872 }
873
874 //==================================================================================================
875 std::shared_ptr<GeomAPI_Pln> GeomAlgoAPI_ShapeTools::findPlane(const ListOfShape& theShapes)
876 {
877   TopoDS_Compound aCompound;
878   BRep_Builder aBuilder;
879   aBuilder.MakeCompound(aCompound);
880
881   for (ListOfShape::const_iterator anIt = theShapes.cbegin(); anIt != theShapes.cend(); ++anIt) {
882     aBuilder.Add(aCompound, (*anIt)->impl<TopoDS_Shape>());
883   }
884   BRepBuilderAPI_FindPlane aFindPlane(aCompound);
885
886   if (aFindPlane.Found() != Standard_True) {
887     return std::shared_ptr<GeomAPI_Pln>();
888   }
889
890   Handle(Geom_Plane) aPlane = aFindPlane.Plane();
891   gp_Pnt aLoc = aPlane->Location();
892   gp_Dir aDir = aPlane->Axis().Direction();
893
894   std::shared_ptr<GeomAPI_Pnt> aGeomPnt(new GeomAPI_Pnt(aLoc.X(), aLoc.Y(), aLoc.Z()));
895   std::shared_ptr<GeomAPI_Dir> aGeomDir(new GeomAPI_Dir(aDir.X(), aDir.Y(), aDir.Z()));
896
897   std::shared_ptr<GeomAPI_Pln> aPln(new GeomAPI_Pln(aGeomPnt, aGeomDir));
898
899   return aPln;
900 }
901
902 //==================================================================================================
903 bool GeomAlgoAPI_ShapeTools::isSubShapeInsideShape(
904   const std::shared_ptr<GeomAPI_Shape> theSubShape,
905   const std::shared_ptr<GeomAPI_Shape> theBaseShape)
906 {
907   if (!theSubShape.get() || !theBaseShape.get()) {
908     return false;
909   }
910
911   const TopoDS_Shape& aSubShape = theSubShape->impl<TopoDS_Shape>();
912   const TopoDS_Shape& aBaseShape = theBaseShape->impl<TopoDS_Shape>();
913
914   if (aSubShape.ShapeType() == TopAbs_VERTEX) {
915     // If sub-shape is a vertex check distance to shape. If it is <= Precision::Confusion() then OK.
916     BRepExtrema_DistShapeShape aDist(aBaseShape, aSubShape);
917     aDist.Perform();
918     if (!aDist.IsDone() || aDist.Value() > Precision::Confusion()) {
919       return false;
920     }
921   } else if (aSubShape.ShapeType() == TopAbs_EDGE) {
922     if (aBaseShape.ShapeType() == TopAbs_FACE) {
923       // Check that edge is on face surface.
924       TopoDS_Face aFace = TopoDS::Face(aBaseShape);
925       TopoDS_Edge anEdge = TopoDS::Edge(aSubShape);
926       BRepLib_CheckCurveOnSurface aCheck(anEdge, aFace);
927       aCheck.Perform();
928       if (!aCheck.IsDone() || aCheck.MaxDistance() > Precision::Confusion()) {
929         return false;
930       }
931
932       // Check intersections.
933       TopoDS_Vertex aV1, aV2;
934       ShapeAnalysis::FindBounds(anEdge, aV1, aV2);
935       gp_Pnt aPnt1 = BRep_Tool::Pnt(aV1);
936       gp_Pnt aPnt2 = BRep_Tool::Pnt(aV2);
937       for (TopExp_Explorer anExp(aBaseShape, TopAbs_EDGE); anExp.More(); anExp.Next()) {
938         const TopoDS_Shape& anEdgeOnFace = anExp.Current();
939         BRepExtrema_DistShapeShape aDist(anEdgeOnFace, anEdge);
940         aDist.Perform();
941         if (aDist.IsDone() && aDist.Value() <= Precision::Confusion()) {
942           // Edge intersect face bound. Check that it is not on edge begin or end.
943           for (Standard_Integer anIndex = 1; anIndex <= aDist.NbSolution(); ++anIndex) {
944             gp_Pnt aPntOnSubShape = aDist.PointOnShape2(anIndex);
945             if (aPntOnSubShape.Distance(aPnt1) > Precision::Confusion()
946                 && aPntOnSubShape.Distance(aPnt2) > Precision::Confusion()) {
947               return false;
948             }
949           }
950         }
951       }
952
953       // No intersections found. Edge is inside or outside face. Check it.
954       BRepAdaptor_Curve aCurveAdaptor(anEdge);
955       gp_Pnt aPointToCheck =
956         aCurveAdaptor.Value((aCurveAdaptor.FirstParameter() +
957                               aCurveAdaptor.LastParameter()) / 2.0);
958       Handle(Geom_Surface) aSurface = BRep_Tool::Surface(aFace);
959       ShapeAnalysis_Surface aSAS(aSurface);
960       gp_Pnt2d aPointOnFace = aSAS.ValueOfUV(aPointToCheck, Precision::Confusion());
961       BRepTopAdaptor_FClass2d aFClass2d(aFace, Precision::Confusion());
962       if (aFClass2d.Perform(aPointOnFace) == TopAbs_OUT) {
963         return false;
964       }
965
966     } else {
967       return false;
968     }
969   } else {
970     return false;
971   }
972
973   return true;
974 }
975
976 //==================================================================================================
977 bool GeomAlgoAPI_ShapeTools::isShapeValid(const std::shared_ptr<GeomAPI_Shape> theShape)
978 {
979   if (!theShape.get()) {
980     return false;
981   }
982
983   BRepCheck_Analyzer aChecker(theShape->impl<TopoDS_Shape>());
984   return (aChecker.IsValid() == Standard_True);
985 }
986
987 //==================================================================================================
988 std::shared_ptr<GeomAPI_Shape>
989   GeomAlgoAPI_ShapeTools::getFaceOuterWire(const std::shared_ptr<GeomAPI_Shape> theFace)
990 {
991   GeomShapePtr anOuterWire;
992
993   if (!theFace.get() || !theFace->isFace()) {
994     return anOuterWire;
995   }
996
997   TopoDS_Face aFace = TopoDS::Face(theFace->impl<TopoDS_Shape>());
998   TopoDS_Wire aWire = BRepTools::OuterWire(aFace);
999
1000   anOuterWire.reset(new GeomAPI_Shape());
1001   anOuterWire->setImpl(new TopoDS_Shape(aWire));
1002
1003   return anOuterWire;
1004 }
1005
1006 //==================================================================================================
1007 static bool boundaryOfEdge(const std::shared_ptr<GeomAPI_Edge> theEdge,
1008                           const std::shared_ptr<GeomAPI_Vertex> theVertex,
1009                           double& theParam)
1010 {
1011   GeomPointPtr aPoint = theVertex->point();
1012   GeomPointPtr aFirstPnt = theEdge->firstPoint();
1013   double aFirstPntTol = theEdge->firstPointTolerance();
1014   GeomPointPtr aLastPnt = theEdge->lastPoint();
1015   double aLastPntTol = theEdge->lastPointTolerance();
1016
1017   double aFirst, aLast;
1018   theEdge->getRange(aFirst, aLast);
1019
1020   bool isFirst = aPoint->distance(aFirstPnt) <= aFirstPntTol;
1021   bool isLast = aPoint->distance(aLastPnt) <= aLastPntTol;
1022   if (isFirst)
1023     theParam = aFirst;
1024   else if (isLast)
1025     theParam = aLast;
1026
1027   return isFirst != isLast;
1028 }
1029
1030 bool GeomAlgoAPI_ShapeTools::isTangent(const std::shared_ptr<GeomAPI_Edge> theEdge1,
1031                                        const std::shared_ptr<GeomAPI_Edge> theEdge2,
1032                                        const std::shared_ptr<GeomAPI_Vertex> theTgPoint)
1033 {
1034   double aParE1 = 0, aParE2 = 0;
1035   if (!boundaryOfEdge(theEdge1, theTgPoint, aParE1) ||
1036       !boundaryOfEdge(theEdge2, theTgPoint, aParE2))
1037     return false;
1038
1039   BRepAdaptor_Curve aC1(theEdge1->impl<TopoDS_Edge>());
1040   BRepAdaptor_Curve aC2(theEdge2->impl<TopoDS_Edge>());
1041   return BRepLProp::Continuity(aC1, aC2, aParE1, aParE2) >= GeomAbs_G1;
1042 }
1043
1044 //==================================================================================================
1045 bool GeomAlgoAPI_ShapeTools::isParallel(const std::shared_ptr<GeomAPI_Edge> theEdge,
1046                                         const std::shared_ptr<GeomAPI_Face> theFace)
1047 {
1048   if (!theEdge.get() || !theFace.get()) {
1049     return false;
1050   }
1051
1052   TopoDS_Edge anEdge = TopoDS::Edge(theEdge->impl<TopoDS_Shape>());
1053   TopoDS_Face aFace  = TopoDS::Face(theFace->impl<TopoDS_Shape>());
1054
1055   BRepExtrema_ExtCF anExt(anEdge, aFace);
1056   return anExt.IsParallel() == Standard_True;
1057 }
1058
1059 //==================================================================================================
1060 std::list<std::shared_ptr<GeomAPI_Vertex> > GeomAlgoAPI_ShapeTools::intersect(
1061   const std::shared_ptr<GeomAPI_Edge> theEdge, const std::shared_ptr<GeomAPI_Face> theFace)
1062 {
1063   std::list<std::shared_ptr<GeomAPI_Vertex> > aResult;
1064   if (!theEdge.get() || !theFace.get()) {
1065     return aResult;
1066   }
1067
1068   TopoDS_Edge anEdge = TopoDS::Edge(theEdge->impl<TopoDS_Shape>());
1069   double aFirstOnCurve, aLastOnCurve;
1070   Handle(Geom_Curve) aCurve = BRep_Tool::Curve(anEdge, aFirstOnCurve, aLastOnCurve);
1071
1072   TopoDS_Face aFace  = TopoDS::Face(theFace->impl<TopoDS_Shape>());
1073   Handle(Geom_Surface) aSurf = BRep_Tool::Surface(aFace);
1074
1075   GeomAPI_IntCS anIntAlgo(aCurve, aSurf);
1076   if (!anIntAlgo.IsDone())
1077     return aResult;
1078   // searching for points-intersection
1079   for (int anIntNum = 1; anIntNum <= anIntAlgo.NbPoints() + anIntAlgo.NbSegments(); anIntNum++) {
1080     gp_Pnt anInt;
1081     if (anIntNum <= anIntAlgo.NbPoints()) {
1082       anInt = anIntAlgo.Point(anIntNum);
1083     } else { // take the middle point on the segment of the intersection
1084       Handle(Geom_Curve) anIntCurve = anIntAlgo.Segment(anIntNum - anIntAlgo.NbPoints());
1085       anIntCurve->D0((anIntCurve->FirstParameter() + anIntCurve->LastParameter()) / 2., anInt);
1086     }
1087     aResult.push_back(std::shared_ptr<GeomAPI_Vertex>(
1088       new GeomAPI_Vertex(anInt.X(), anInt.Y(), anInt.Z())));
1089   }
1090   return aResult;
1091 }
1092
1093 //==================================================================================================
1094 void GeomAlgoAPI_ShapeTools::splitShape(const std::shared_ptr<GeomAPI_Shape>& theBaseShape,
1095                                       const GeomAlgoAPI_ShapeTools::PointToRefsMap& thePointsInfo,
1096                                       std::set<std::shared_ptr<GeomAPI_Shape> >& theShapes)
1097 {
1098   // to split shape at least one point should be presented in the points container
1099   if (thePointsInfo.empty())
1100     return;
1101
1102     // General Fuse to split edge by vertices
1103   BOPAlgo_Builder aBOP;
1104   TopoDS_Edge aBaseEdge = theBaseShape->impl<TopoDS_Edge>();
1105   // Rebuild closed edge to place vertex to one of split points.
1106   // This will prevent edge to be split on same vertex.
1107   if (BRep_Tool::IsClosed(aBaseEdge))
1108   {
1109     Standard_Real aFirst, aLast;
1110     Handle(Geom_Curve) aCurve = BRep_Tool::Curve(aBaseEdge, aFirst, aLast);
1111
1112     PointToRefsMap::const_iterator aPIt = thePointsInfo.begin();
1113     std::shared_ptr<GeomAPI_Pnt> aPnt = aPIt->first;
1114     gp_Pnt aPoint(aPnt->x(), aPnt->y(), aPnt->z());
1115
1116     TopAbs_Orientation anOrientation = aBaseEdge.Orientation();
1117     aBaseEdge = BRepBuilderAPI_MakeEdge(aCurve, aPoint, aPoint).Edge();
1118     aBaseEdge.Orientation(anOrientation);
1119   }
1120   aBOP.AddArgument(aBaseEdge);
1121
1122   PointToRefsMap::const_iterator aPIt = thePointsInfo.begin();
1123   for (; aPIt != thePointsInfo.end(); ++aPIt) {
1124     std::shared_ptr<GeomAPI_Pnt> aPnt = aPIt->first;
1125     TopoDS_Vertex aV = BRepBuilderAPI_MakeVertex(gp_Pnt(aPnt->x(), aPnt->y(), aPnt->z()));
1126     aBOP.AddArgument(aV);
1127   }
1128
1129   aBOP.Perform();
1130   if (aBOP.HasErrors())
1131     return;
1132
1133   // Collect splits
1134   const TopTools_ListOfShape& aSplits = aBOP.Modified(aBaseEdge);
1135   TopTools_ListIteratorOfListOfShape anIt(aSplits);
1136   for (; anIt.More(); anIt.Next()) {
1137     std::shared_ptr<GeomAPI_Shape> anEdge(new GeomAPI_Shape);
1138     anEdge->setImpl(new TopoDS_Shape(anIt.Value()));
1139     theShapes.insert(anEdge);
1140   }
1141 }
1142
1143 //==================================================================================================
1144 void GeomAlgoAPI_ShapeTools::splitShape_p(const std::shared_ptr<GeomAPI_Shape>& theBaseShape,
1145                                           const std::list<std::shared_ptr<GeomAPI_Pnt> >& thePoints,
1146                                           std::set<std::shared_ptr<GeomAPI_Shape> >& theShapes)
1147 {
1148   // General Fuse to split edge by vertices
1149   BOPAlgo_Builder aBOP;
1150   TopoDS_Edge aBaseEdge = theBaseShape->impl<TopoDS_Edge>();
1151   // Rebuild closed edge to place vertex to one of split points.
1152   // This will prevent edge to be split on seam vertex.
1153   if (BRep_Tool::IsClosed(aBaseEdge))
1154   {
1155     Standard_Real aFirst, aLast;
1156     Handle(Geom_Curve) aCurve = BRep_Tool::Curve(aBaseEdge, aFirst, aLast);
1157
1158     std::list<std::shared_ptr<GeomAPI_Pnt> >::const_iterator aPIt = thePoints.begin();
1159     gp_Pnt aPoint((*aPIt)->x(), (*aPIt)->y(), (*aPIt)->z());
1160
1161     TopAbs_Orientation anOrientation = aBaseEdge.Orientation();
1162     aBaseEdge = BRepBuilderAPI_MakeEdge(aCurve, aPoint, aPoint).Edge();
1163     aBaseEdge.Orientation(anOrientation);
1164   }
1165   aBOP.AddArgument(aBaseEdge);
1166
1167   std::list<std::shared_ptr<GeomAPI_Pnt> >::const_iterator aPtIt = thePoints.begin();
1168   for (; aPtIt != thePoints.end(); ++aPtIt) {
1169     std::shared_ptr<GeomAPI_Pnt> aPnt = *aPtIt;
1170     TopoDS_Vertex aV = BRepBuilderAPI_MakeVertex(gp_Pnt(aPnt->x(), aPnt->y(), aPnt->z()));
1171     aBOP.AddArgument(aV);
1172   }
1173
1174   aBOP.Perform();
1175   if (aBOP.HasErrors())
1176     return;
1177
1178   // Collect splits
1179   const TopTools_ListOfShape& aSplits = aBOP.Modified(aBaseEdge);
1180   TopTools_ListIteratorOfListOfShape anIt(aSplits);
1181   for (; anIt.More(); anIt.Next()) {
1182     std::shared_ptr<GeomAPI_Shape> anEdge(new GeomAPI_Shape);
1183     anEdge->setImpl(new TopoDS_Shape(anIt.Value()));
1184     theShapes.insert(anEdge);
1185   }
1186 }
1187
1188 //==================================================================================================
1189 std::shared_ptr<GeomAPI_Shape> GeomAlgoAPI_ShapeTools::findShape(
1190                                   const std::list<std::shared_ptr<GeomAPI_Pnt> >& thePoints,
1191                                   const std::set<std::shared_ptr<GeomAPI_Shape> >& theShapes)
1192 {
1193   std::shared_ptr<GeomAPI_Shape> aResultShape;
1194
1195   if (thePoints.size() == 2) {
1196     std::list<std::shared_ptr<GeomAPI_Pnt> >::const_iterator aPntIt = thePoints.begin();
1197     std::shared_ptr<GeomAPI_Pnt> aFirstPoint = *aPntIt;
1198     aPntIt++;
1199     std::shared_ptr<GeomAPI_Pnt> aLastPoint = *aPntIt;
1200
1201     std::set<std::shared_ptr<GeomAPI_Shape> >::const_iterator anIt = theShapes.begin(),
1202                                                               aLast = theShapes.end();
1203     for (; anIt != aLast; anIt++) {
1204       GeomShapePtr aShape = *anIt;
1205       std::shared_ptr<GeomAPI_Edge> anEdge(new GeomAPI_Edge(aShape));
1206       if (anEdge.get()) {
1207         std::shared_ptr<GeomAPI_Pnt> anEdgeFirstPoint = anEdge->firstPoint();
1208         std::shared_ptr<GeomAPI_Pnt> anEdgeLastPoint = anEdge->lastPoint();
1209         if (anEdgeFirstPoint->isEqual(aFirstPoint) &&
1210             anEdgeLastPoint->isEqual(aLastPoint))
1211             aResultShape = aShape;
1212       }
1213     }
1214   }
1215
1216   return aResultShape;
1217 }
1218
1219 //==================================================================================================
1220 #ifdef FEATURE_MULTIROTATION_TWO_DIRECTIONS
1221 std::shared_ptr<GeomAPI_Dir> GeomAlgoAPI_ShapeTools::buildDirFromAxisAndShape(
1222                                     const std::shared_ptr<GeomAPI_Shape> theBaseShape,
1223                                     const std::shared_ptr<GeomAPI_Ax1> theAxis)
1224 {
1225   gp_Pnt aCentreOfMassPoint =
1226     GeomAlgoAPI_ShapeTools::centreOfMass(theBaseShape)->impl<gp_Pnt>();
1227   Handle(Geom_Line) aLine = new Geom_Line(theAxis->impl<gp_Ax1>());
1228   GeomAPI_ProjectPointOnCurve aPrjTool(aCentreOfMassPoint, aLine);
1229   gp_Pnt aPoint = aPrjTool.NearestPoint();
1230
1231   std::shared_ptr<GeomAPI_Dir> aDir(new GeomAPI_Dir(aCentreOfMassPoint.X()-aPoint.X(),
1232                                                     aCentreOfMassPoint.Y()-aPoint.Y(),
1233                                                     aCentreOfMassPoint.Z()-aPoint.Z()));
1234   return aDir;
1235 }
1236 #endif
1237
1238 //==================================================================================================
1239 static TopoDS_Wire fixParametricGaps(const TopoDS_Wire& theWire)
1240 {
1241   TopoDS_Wire aFixedWire;
1242   Handle(Geom_Curve) aPrevCurve;
1243   double aPrevLastParam = -Precision::Infinite();
1244
1245   BRep_Builder aBuilder;
1246   aBuilder.MakeWire(aFixedWire);
1247
1248   BRepTools_WireExplorer aWExp(theWire);
1249   for (; aWExp.More(); aWExp.Next()) {
1250     TopoDS_Edge anEdge = aWExp.Current();
1251     double aFirst, aLast;
1252     Handle(Geom_Curve) aCurve = BRep_Tool::Curve(anEdge, aFirst, aLast);
1253     if (aCurve == aPrevCurve && Abs(aFirst - aPrevLastParam) > Precision::Confusion()) {
1254       // if parametric gap occurs, create new edge based on the copied curve
1255       aCurve = Handle(Geom_Curve)::DownCast(aCurve->Copy());
1256       TopoDS_Vertex aV1, aV2;
1257       TopExp::Vertices(anEdge, aV1, aV2);
1258       anEdge = TopoDS::Edge(anEdge.EmptyCopied());
1259       aBuilder.UpdateEdge(anEdge, aCurve, BRep_Tool::Tolerance(anEdge));
1260       aBuilder.Add(anEdge, aV1);
1261       aBuilder.Add(anEdge, aV2);
1262     }
1263
1264     aBuilder.Add(aFixedWire, anEdge);
1265
1266     aPrevCurve = aCurve;
1267     aPrevLastParam = aLast;
1268   }
1269
1270   return aFixedWire;
1271 }
1272
1273 //==================================================================================================
1274 std::shared_ptr<GeomAPI_Edge> GeomAlgoAPI_ShapeTools::wireToEdge(
1275       const std::shared_ptr<GeomAPI_Wire>& theWire)
1276 {
1277   GeomEdgePtr anEdge;
1278   if (theWire) {
1279     TopoDS_Wire aWire = theWire->impl<TopoDS_Wire>();
1280     BRepTools_WireExplorer aWExp(aWire);
1281     TopoDS_Edge aNewEdge = aWExp.Current();
1282     aWExp.Next();
1283     if (aWExp.More()) {
1284       // Workaround for the closed wire to avoid jumping of its start point:
1285       // split this wire for two parts, convert them to edges, then compose together
1286       if (BRep_Tool::IsClosed(aWire)) {
1287         aWire = TopoDS::Wire(BRepBuilderAPI_Copy(aWire).Shape());
1288         aWExp.Init(aWire);
1289         aNewEdge = aWExp.Current();
1290
1291         BRep_Builder().Remove(aWire, aNewEdge);
1292         GeomWirePtr aSplitWire(new GeomAPI_Wire);
1293         aSplitWire->setImpl(new TopoDS_Wire(aWire));
1294         GeomEdgePtr aMergedEdge = wireToEdge(aSplitWire);
1295
1296         aWire = BRepBuilderAPI_MakeWire(aNewEdge, aMergedEdge->impl<TopoDS_Edge>());
1297       }
1298
1299       // Workaround: when concatenate a wire consisting of two edges based on the same B-spline
1300       // curve (non-periodic, but having equal start and end points), first of which is placed
1301       // at the end on the curve and second is placed at the start, this workaround copies
1302       // second curve to avoid treating these edges as a single curve by setting trim parameters.
1303       aWire = fixParametricGaps(aWire);
1304       aWire = BRepAlgo::ConcatenateWire(aWire, GeomAbs_G1); // join smooth parts of wire
1305       aNewEdge = BRepAlgo::ConcatenateWireC0(aWire); // join C0 parts of wire
1306
1307       // Reapproximate the result edge to have the parameter equal to curvilinear abscissa.
1308       static const int THE_MAX_DEGREE = 14;
1309       static const int THE_MAX_INTERVALS = 32;
1310       double aFirst, aLast;
1311       Handle(Geom_Curve) aCurve = BRep_Tool::Curve(aNewEdge, aFirst, aLast);
1312 #if OCC_VERSION_LARGE < 0x07070000
1313       Handle(GeomAdaptor_HCurve) aHCurve = new GeomAdaptor_HCurve(aCurve);
1314 #else
1315       Handle(GeomAdaptor_Curve) aHCurve = new GeomAdaptor_Curve(aCurve);
1316 #endif
1317       Approx_CurvilinearParameter anApprox(aHCurve, Precision::Confusion(), aCurve->Continuity(),
1318                                            THE_MAX_DEGREE, THE_MAX_INTERVALS);
1319       if (anApprox.HasResult()) {
1320         Handle(Geom_BSplineCurve) aNewCurve = anApprox.Curve3d();
1321         TColStd_Array1OfReal aKnots = aNewCurve->Knots();
1322         BSplCLib::Reparametrize(aFirst, aLast, aKnots);
1323         aNewCurve->SetKnots(aKnots);
1324         BRep_Builder().UpdateEdge(aNewEdge, aNewCurve, BRep_Tool::Tolerance(aNewEdge));
1325       }
1326     }
1327     anEdge = GeomEdgePtr(new GeomAPI_Edge);
1328     anEdge->setImpl(new TopoDS_Edge(aNewEdge));
1329   }
1330   return anEdge;
1331 }
1332
1333 //==================================================================================================
1334 ListOfShape GeomAlgoAPI_ShapeTools::getLowLevelSubShapes(const GeomShapePtr& theShape)
1335 {
1336   ListOfShape aSubShapes;
1337
1338   if (!theShape->isCompound() && !theShape->isCompSolid() &&
1339       !theShape->isShell() && !theShape->isWire()) {
1340     return aSubShapes;
1341   }
1342
1343   for (GeomAPI_ShapeIterator anIt(theShape); anIt.more(); anIt.next()) {
1344     GeomShapePtr aSubShape = anIt.current();
1345     if (aSubShape->isVertex() || aSubShape->isEdge() ||
1346         aSubShape->isFace() || aSubShape->isSolid()) {
1347       aSubShapes.push_back(aSubShape);
1348     } else {
1349       aSubShapes.splice(aSubShapes.end(), getLowLevelSubShapes(aSubShape));
1350     }
1351   }
1352
1353   return aSubShapes;
1354 }
1355
1356 //==================================================================================================
1357 static void getMinMaxPointsOnLine(const std::list<std::shared_ptr<GeomAPI_Pnt> >& thePoints,
1358                                   const gp_Dir theDir,
1359                                   double& theMin, double& theMax)
1360 {
1361   theMin = RealLast();
1362   theMax = RealFirst();
1363   // Project bounding points on theDir
1364   for (std::list<std::shared_ptr<GeomAPI_Pnt> >::const_iterator
1365          aPointsIt = thePoints.begin(); aPointsIt != thePoints.end(); aPointsIt++) {
1366     const gp_Pnt& aPnt = (*aPointsIt)->impl<gp_Pnt>();
1367     gp_Dir aPntDir (aPnt.XYZ());
1368     Standard_Real proj = (theDir*aPntDir) * aPnt.XYZ().Modulus();
1369     if (proj < theMin) theMin = proj;
1370     if (proj > theMax) theMax = proj;
1371   }
1372 }
1373
1374 //==================================================================================================
1375 void GeomAlgoAPI_ShapeTools::computeThroughAll(const ListOfShape& theObjects,
1376                                                const ListOfShape& theBaseShapes,
1377                                                const std::shared_ptr<GeomAPI_Dir> theDir,
1378                                                double& theToSize, double& theFromSize)
1379 {
1380   // Bounding box of objects
1381   std::list<std::shared_ptr<GeomAPI_Pnt> > aBndObjs =
1382       GeomAlgoAPI_ShapeTools::getBoundingBox(theObjects);
1383   if (aBndObjs.size() != 8) {
1384     return;
1385   }
1386
1387   // the value to enlarge the bounding box of each object to make the extruded shape
1388   // a little bit larger than overall objects to get the correct result of Boolean CUT operation
1389   double anEnlargement = 0.1 * aBndObjs.front()->distance(aBndObjs.back());
1390
1391   // Prism direction
1392   if (theDir.get()) {
1393     // One direction for all prisms
1394     gp_Dir aDir = theDir->impl<gp_Dir>();
1395
1396     // Bounding box of the base
1397     std::list<std::shared_ptr<GeomAPI_Pnt> > aBndBases =
1398         GeomAlgoAPI_ShapeTools::getBoundingBox(theBaseShapes);
1399     if (aBndBases.size() != 8) {
1400       return;
1401     }
1402
1403     // Objects bounds
1404     Standard_Real lowBnd, upperBnd;
1405     getMinMaxPointsOnLine(aBndObjs, aDir, lowBnd, upperBnd);
1406
1407     // Base bounds
1408     Standard_Real lowBase, upperBase;
1409     getMinMaxPointsOnLine(aBndBases, aDir, lowBase, upperBase);
1410
1411     // ----------.-----.---------.--------------.-----------> theDir
1412     //       lowBnd   lowBase   upperBase    upperBnd
1413
1414     theToSize = upperBnd - lowBase;
1415     theFromSize = upperBase - lowBnd;
1416   } else {
1417     // Direction is a normal to each base shape (different normals to bases)
1418     // So we calculate own sizes for each base shape
1419     theToSize = 0.0;
1420     theFromSize = 0.0;
1421
1422     for (ListOfShape::const_iterator anIt = theBaseShapes.begin();
1423          anIt != theBaseShapes.end(); ++anIt) {
1424       const GeomShapePtr& aBaseShape_i = (*anIt);
1425       ListOfShape aBaseShapes_i;
1426       aBaseShapes_i.push_back(aBaseShape_i);
1427
1428       // Bounding box of the base
1429       std::list<std::shared_ptr<GeomAPI_Pnt> > aBndBases =
1430           GeomAlgoAPI_ShapeTools::getBoundingBox(aBaseShapes_i, anEnlargement);
1431       if (aBndBases.size() != 8) {
1432         return;
1433       }
1434
1435       // Direction (normal to aBaseShapes_i)
1436       // Code like in GeomAlgoAPI_Prism
1437       gp_Dir aDir;
1438       const TopoDS_Shape& aBaseShape = aBaseShape_i->impl<TopoDS_Shape>();
1439       BRepBuilderAPI_FindPlane aFindPlane(aBaseShape);
1440       if (aFindPlane.Found() == Standard_True) {
1441         Handle(Geom_Plane) aPlane;
1442         if (aBaseShape.ShapeType() == TopAbs_FACE || aBaseShape.ShapeType() == TopAbs_SHELL) {
1443           TopExp_Explorer anExp(aBaseShape, TopAbs_FACE);
1444           const TopoDS_Shape& aFace = anExp.Current();
1445           Handle(Geom_Surface) aSurface = BRep_Tool::Surface(TopoDS::Face(aFace));
1446           if (aSurface->DynamicType() == STANDARD_TYPE(Geom_RectangularTrimmedSurface)) {
1447             Handle(Geom_RectangularTrimmedSurface) aTrimSurface =
1448               Handle(Geom_RectangularTrimmedSurface)::DownCast(aSurface);
1449             aSurface = aTrimSurface->BasisSurface();
1450           }
1451           if (aSurface->DynamicType() != STANDARD_TYPE(Geom_Plane)) {
1452             return;
1453           }
1454           aPlane = Handle(Geom_Plane)::DownCast(aSurface);
1455         } else {
1456           aPlane = aFindPlane.Plane();
1457         }
1458         aDir = aPlane->Axis().Direction();
1459       } else {
1460         return;
1461       }
1462
1463       // Objects bounds
1464       Standard_Real lowBnd, upperBnd;
1465       getMinMaxPointsOnLine(aBndObjs, aDir, lowBnd, upperBnd);
1466
1467       // Base bounds
1468       Standard_Real lowBase, upperBase;
1469       getMinMaxPointsOnLine(aBndBases, aDir, lowBase, upperBase);
1470
1471       // ----------.-----.---------.--------------.-----------> theDir
1472       //       lowBnd   lowBase   upperBase    upperBnd
1473
1474       double aToSize_i = upperBnd - lowBase;
1475       double aFromSize_i = upperBase - lowBnd;
1476
1477       if (aToSize_i > theToSize) theToSize = aToSize_i;
1478       if (aFromSize_i > theFromSize) theFromSize = aFromSize_i;
1479     }
1480   }
1481 }