]> OCCT Git - occt-copy.git/commitdiff
0030599: Visualization - implement StdPrs_SectionLines tool CR30599
authortiv <tiv@opencascade.com>
Fri, 19 Apr 2019 17:27:36 +0000 (20:27 +0300)
committertiv <tiv@opencascade.com>
Fri, 19 Apr 2019 17:27:36 +0000 (20:27 +0300)
StdPrs_SectionLines tool is implemented: it is possible to compute a cross section of the shape with the plane.
# Some features have not been implemented yet:
#   * the method StdPrs_SectionLines::Add() that should add objects inside section to the Prs3d_Presentation;
#   * corresponding Draw Harness command;
#   * interactive object that could represent the cross section;
#   * triangles that lie inside the section plane have not been processed yet.

src/Bnd/Bnd_Box.cxx
src/Bnd/Bnd_Box.hxx
src/StdPrs/FILES
src/StdPrs/StdPrs_SectionLines.cxx [new file with mode: 0644]
src/StdPrs/StdPrs_SectionLines.hxx [new file with mode: 0644]
src/gp/gp_Pln.hxx
src/gp/gp_Pln.lxx

index 3f9ea999ffbd32bc8d797e002328b740c5a35a99..67c711314b67e0264910ef63f3522ae74ce80223 100644 (file)
@@ -210,6 +210,20 @@ gp_Pnt Bnd_Box::CornerMin() const
   return aCornerMin;
 }
 
+//=======================================================================
+// function : CornerMin
+// purpose  :
+//=======================================================================
+
+gp_Pnt Bnd_Box::Center() const
+{
+  if (IsVoid())
+  {
+    throw Standard_ConstructionError("Bnd_Box is void");
+  }
+  return (CornerMin().Coord() + CornerMax().Coord()) / 2.0;
+}
+
 //=======================================================================
 //function : CornerMax
 //purpose  :
index 834338f5bf66aaa4f7ca76376d2b404a484dbaef..56582ab1ac669a2159e4ac5a416f78ca810cb10d 100644 (file)
@@ -142,6 +142,12 @@ public:
   //! Standard_ConstructionError exception will be thrown if the box is void.
   //! if IsVoid()
   Standard_EXPORT gp_Pnt CornerMax() const;
+  
+  //! Returns the center of this bounding box.
+  //! If this bounding box is infinite (i.e. "open"), returned values
+  //! may be equal to +/- Precision::Infinite().
+  //! Standard_ConstructionError exception will be thrown if the box is void.
+  Standard_EXPORT gp_Pnt Center() const;
 
   //! The   Box will be   infinitely   long  in the Xmin
   //! direction.
index 94ffef116e245a94c7b200c8a08fa3dd55c44a70..c97dce05b1593172b926c07d4d855aae968f41ae 100644 (file)
@@ -17,6 +17,8 @@ StdPrs_Plane.hxx
 StdPrs_Point.hxx
 StdPrs_PoleCurve.cxx
 StdPrs_PoleCurve.hxx
+StdPrs_SectionLines.cxx
+StdPrs_SectionLines.hxx
 StdPrs_ShadedShape.cxx
 StdPrs_ShadedShape.hxx
 StdPrs_ShadedSurface.cxx
diff --git a/src/StdPrs/StdPrs_SectionLines.cxx b/src/StdPrs/StdPrs_SectionLines.cxx
new file mode 100644 (file)
index 0000000..653de49
--- /dev/null
@@ -0,0 +1,376 @@
+// Created on: 2019-04-11
+// Created by: Timur Izmaylov
+// Copyright (c) 2019 OPEN CASCADE SAS
+//
+// This file is part of Open CASCADE Technology software library.
+//
+// This library is free software; you can redistribute it and/or modify it under
+// the terms of the GNU Lesser General Public License version 2.1 as published
+// by the Free Software Foundation, with special exception defined in the file
+// OCCT_LGPL_EXCEPTION.txt. Consult the file LICENSE_LGPL_21.txt included in OCCT
+// distribution for complete text of the license and disclaimer of any warranty.
+//
+// Alternatively, this file may be used under the terms of Open CASCADE
+// commercial license or contractual agreement.
+
+#include <StdPrs_SectionLines.hxx>
+
+#include <BRepBndLib.hxx>
+#include <BRepBuilderAPI_MakeEdge.hxx>
+#include <BRep_Tool.hxx>
+#include <Bnd_Box.hxx>
+#include <Graphic3d_ArrayOfSegments.hxx>
+#include <Poly_Triangulation.hxx>
+#include <Precision.hxx>
+#include <TColgp_SequenceOfPnt.hxx>
+#include <TopExp_Explorer.hxx>
+#include <TopoDS.hxx>
+#include <TopoDS_Builder.hxx>
+#include <TopoDS_Edge.hxx>
+#include <TopoDS_Face.hxx>
+#include <gp_Pln.hxx>
+
+#include <vector>
+
+namespace
+{
+
+  //! Pair of points
+  typedef std::pair<gp_Pnt, gp_Pnt> PointPair;
+
+  //! Segment of line
+  typedef PointPair Segment;
+
+  //! Array of segments
+  typedef std::vector<Segment> Segments;
+
+  //! Triple of points
+  typedef gp_Pnt PointTriple[3];
+
+  //! Triangle (a set of 3 points)
+  typedef PointTriple Triangle;
+
+  //! Array of triangles
+  typedef std::vector<Triangle> Triangles;
+
+  //! Adds the segment to the array of segments
+  //! @param thePoint1 the first end of segment
+  //! @param thePoint2 the second end of segment
+  //! @param theSegments the array of segments
+  static void addSegment (const gp_Pnt& thePoint1, const gp_Pnt& thePoint2, Segments& theSegments)
+  {
+    if (thePoint1.Distance (thePoint2) <= Precision::Confusion())
+    {
+      return;
+    }
+    const PointPair anEdge (thePoint1, thePoint2);
+    theSegments.push_back (anEdge);
+  }
+
+  //! Adds the point to the section
+  //! @param theVertex1 the first vertex of the crossed edge
+  //! @param theDistanceToPlane1 the signed distance from the first vertex of the crossed edge to the plane
+  //! @param theVertex2 the second vertex of the crossed edge
+  //! @param theDistanceToPlane2 the signed distance from the second vertex of the crossed edge to the plane
+  //! @param theSectionPoints the container of section points where found point should be added
+  //! @param theNumberOfSectionPoints the number of points which are already in theSectionPoints container
+  static void addSectionPoint (const gp_Pnt&       theVertex1,
+                               const Standard_Real theDistanceToPlane1,
+                               const gp_Pnt&       theVertex2,
+                               const Standard_Real theDistanceToPlane2,
+                               PointTriple&        theSectionPoints,
+                               std::size_t&        theNumberOfSectionPoints)
+  {
+    if (theDistanceToPlane1 * theDistanceToPlane2 > 0.0)
+    {
+      return;
+    }
+    const Standard_Real aDistanceParameter     = theDistanceToPlane1 / (theDistanceToPlane1 - theDistanceToPlane2);
+    theSectionPoints[theNumberOfSectionPoints] = theVertex1.Coord()
+                                                 + aDistanceParameter * (theVertex2.Coord() - theVertex1.Coord());
+    ++theNumberOfSectionPoints;
+  }
+
+  //! Add points found as the result of crossing the triangle with the plane to the section
+  //! @param theTriangle the triangle that is crossed with the plane
+  //! @param areVerticesOnSectionPlane the array of flags that shows whether the triangle vertex
+  //! with corresponding index is on the section plane
+  //! @param theDistancesToPlane the array of signed distances from the plane to the triangle vertex with corresponding
+  //! index
+  //! @param theSectionPoints the container of section points where found point should be added
+  //! @param theNumberOfSectionPoints the number of points which are already in theSectionPoints container
+  static void addSectionPoints (const Triangle& theTriangle,
+                                const bool (&areVerticesOnSectionPlane)[3],
+                                const Standard_Real (&theDistancesToPlane)[3],
+                                PointTriple& theSectionPoints,
+                                std::size_t& theNumberOfSectionPoints)
+  {
+    for (std::size_t aVertexIndex1 = 0; aVertexIndex1 < 2; ++aVertexIndex1)
+    {
+      const bool isVertex1OnSectionPlane = areVerticesOnSectionPlane[aVertexIndex1];
+      if (isVertex1OnSectionPlane)
+      {
+        continue;
+      }
+      for (std::size_t aVertexIndex2 = aVertexIndex1 + 1; aVertexIndex2 < 3; ++aVertexIndex2)
+      {
+        const bool isVertex2OnSectionPlane = areVerticesOnSectionPlane[aVertexIndex2];
+        if (isVertex2OnSectionPlane)
+        {
+          continue;
+        }
+        addSectionPoint (theTriangle[aVertexIndex1],
+                         theDistancesToPlane[aVertexIndex1],
+                         theTriangle[aVertexIndex2],
+                         theDistancesToPlane[aVertexIndex2],
+                         theSectionPoints,
+                         theNumberOfSectionPoints);
+      }
+    }
+  }
+
+  //! Add all segments between points (2 or 3) from point triple to the section
+  //! @param theSectionPoints the points to be added
+  //! @param theNumberOfSectionPoints the number of section points to be added
+  //! @param theSectionSegments the container of section segments where the new segments should be added
+  static void addAllSectionPoints (const PointTriple& theSectionPoints,
+                                   const std::size_t  theNumberOfSectionPoints,
+                                   Segments&          theSectionSegments)
+  {
+    if (theNumberOfSectionPoints < 2)
+    {
+      return;
+    }
+    for (std::size_t aPointIndex1 = 0; aPointIndex1 < theNumberOfSectionPoints - 1; ++aPointIndex1)
+    {
+      for (std::size_t aPointIndex2 = aPointIndex1 + 1; aPointIndex2 < theNumberOfSectionPoints; ++aPointIndex2)
+      {
+        addSegment (theSectionPoints[aPointIndex1], theSectionPoints[aPointIndex2], theSectionSegments);
+      }
+    }
+  }
+
+  //! Add all segments between triangle vertices that are on the section plane
+  //! @param theTriangle the triangle which vertices should be added
+  //! @param areVerticesOnSectionPlane the array of flags that shows whether the triangle vertex
+  //! with corresponding index is on the section plane
+  //! @param theSectionSegments the container of section segments where the new segments should be added
+  static void addAllTriangleVerticesOnSectionPlane (const Triangle& theTriangle,
+                                                    const bool (&areVerticesOnSectionPlane)[3],
+                                                    Segments& theSectionSegments)
+  {
+    for (std::size_t aVertexIndex1 = 0; aVertexIndex1 < 2; ++aVertexIndex1)
+    {
+      const bool isVertex1OnSectionPlane = areVerticesOnSectionPlane[aVertexIndex1];
+      if (!isVertex1OnSectionPlane)
+      {
+        continue;
+      }
+      for (std::size_t aVertexIndex2 = aVertexIndex1 + 1; aVertexIndex2 < 3; ++aVertexIndex2)
+      {
+        const bool isVertex2OnSectionPlane = areVerticesOnSectionPlane[aVertexIndex2];
+        if (!isVertex2OnSectionPlane)
+        {
+          continue;
+        }
+        addSegment (theTriangle[aVertexIndex1], theTriangle[aVertexIndex2], theSectionSegments);
+      }
+    }
+  }
+
+  //! Computes a cross section of the triangle with the plane
+  //! @param theTriangle the triangle to be sectioned
+  //! @param theSectionPlane the section plane
+  //! @param theSectionSegments the container of section segments where the new segments should be added
+  static void crossSection (const Triangle& theTriangle, const gp_Pln& theSectionPlane, Segments& theSectionSegments)
+  {
+    Standard_Real aDistancesToPlane[3];
+    bool          areVerticesOnSectionPlane[3];
+    std::size_t   aNumberOfVerticesOnSectionPlane = 0;
+    for (std::size_t aVertexIndex = 0; aVertexIndex < 3; ++aVertexIndex)
+    {
+      const gp_Pnt   aNode            = theTriangle[aVertexIndex];
+      Standard_Real& aDistanceToPlane = aDistancesToPlane[aVertexIndex];
+      aDistanceToPlane                = theSectionPlane.SignedDistance (aNode);
+      bool& isVertexOnSectionPlane    = areVerticesOnSectionPlane[aVertexIndex];
+      isVertexOnSectionPlane          = (Abs (aDistanceToPlane) <= Precision::Confusion());
+      if (isVertexOnSectionPlane)
+      {
+        ++aNumberOfVerticesOnSectionPlane;
+      }
+    }
+    if (aNumberOfVerticesOnSectionPlane <= 1)
+    {
+      PointTriple aSectionPoints;
+      std::size_t aNumberOfSectionPoints = 0;
+      for (std::size_t aVertexIndex = 0; aVertexIndex < 3; ++aVertexIndex)
+      {
+        if (areVerticesOnSectionPlane[aVertexIndex])
+        {
+          aSectionPoints[aNumberOfSectionPoints] = theTriangle[aVertexIndex];
+          ++aNumberOfSectionPoints;
+        }
+      }
+      addSectionPoints (
+        theTriangle, areVerticesOnSectionPlane, aDistancesToPlane, aSectionPoints, aNumberOfSectionPoints);
+      addAllSectionPoints (aSectionPoints, aNumberOfSectionPoints, theSectionSegments);
+    }
+    else
+    {
+      addAllTriangleVerticesOnSectionPlane (theTriangle, areVerticesOnSectionPlane, theSectionSegments);
+    }
+  }
+
+  //! Gets a triangle of the triangulated face
+  //! @param theNodes the nodes of the face
+  //! @param thePolyTriangle the indices of the triangle nodes
+  //! @param theFaceLocation the location of the face
+  //! @param theTriangle a triangle that should be extracted from the triangulated face
+  static void getTriangleOnFace (const TColgp_Array1OfPnt& theNodes,
+                                 const Poly_Triangle&      thePolyTriangle,
+                                 const TopLoc_Location&    theFaceLocation,
+                                 Triangle&                 theTriangle)
+  {
+    Standard_Integer aVertexNumbers[3];
+    thePolyTriangle.Get (aVertexNumbers[0], aVertexNumbers[1], aVertexNumbers[2]);
+    for (std::size_t aVertexIndex = 0; aVertexIndex < 3; ++aVertexIndex)
+    {
+      const Standard_Integer aNodeNumber      = aVertexNumbers[aVertexIndex];
+      const gp_Pnt           aNode            = theNodes (aNodeNumber);
+      const gp_Pnt           aTransformedNode = aNode.Transformed (theFaceLocation);
+      theTriangle[aVertexIndex]               = aTransformedNode;
+    }
+  }
+
+  //! Computes a cross section of the face triangulation with the plane
+  //! @param theFaceTriangulation the triangulation of the face to be sectioned
+  //! @param theFaceLocation the location of the face
+  //! @param theSectionPlane the section plane
+  //! @param theSectionSegments the container of section segments where the new segments should be added
+  static void crossSection (const Handle (Poly_Triangulation) & theFaceTriangulation,
+                            const TopLoc_Location& theFaceLocation,
+                            const gp_Pln&          theSectionPlane,
+                            Segments&              theSectionSegments)
+  {
+    const Poly_Array1OfTriangle& aTriangles = theFaceTriangulation->Triangles();
+    const TColgp_Array1OfPnt&    aFaceNodes = theFaceTriangulation->Nodes();
+    for (Standard_Integer aTriangleIndex = 1; aTriangleIndex <= aTriangles.Size(); ++aTriangleIndex)
+    {
+      const Poly_Triangle& aPolyTriangle = aTriangles (aTriangleIndex);
+      Triangle             aTriangle;
+      getTriangleOnFace (aFaceNodes, aPolyTriangle, theFaceLocation, aTriangle);
+      crossSection (aTriangle, theSectionPlane, theSectionSegments);
+    }
+  }
+
+  //! Checks whether the shape's bounding box intersects with the plane
+  //! @param theShape the shape which bounding box takes part in checking on intersection with the given plane
+  //! @param thePlane the plane that takes part in checking on intersection with the given shape's bounding box
+  //! @return true if the shape's bounding box intersects with the plane, or false otherwise
+  static bool intersects (const TopoDS_Shape& theShape, const gp_Pln& thePlane)
+  {
+    Bnd_Box aBoundingBox;
+    BRepBndLib::Add (theShape, aBoundingBox);
+    return !aBoundingBox.IsOut (thePlane);
+  }
+
+  //! Computes a cross section of the face triangulation with the plane
+  //! @param theFaceTriangulation the triangulation of the face to be sectioned
+  //! @param theFaceLocation the location of the face
+  //! @param theSectionPlane the section plane
+  //! @param theSectionSegments the container of section segments where the new segments should be added
+  static void crossSection (const TopoDS_Face& theFace, const gp_Pln& theSectionPlane, Segments& theSectionSegments)
+  {
+    if (!intersects (theFace, theSectionPlane))
+    {
+      return;
+    }
+    TopLoc_Location aFaceLocation;
+    const Handle (Poly_Triangulation)& aFaceTriangulation = BRep_Tool::Triangulation (theFace, aFaceLocation);
+    if (aFaceTriangulation.IsNull())
+    {
+      return;
+    }
+    crossSection (aFaceTriangulation, aFaceLocation, theSectionPlane, theSectionSegments);
+  }
+
+  //! Computes a cross section of the shape with the plane
+  //! @param theShape the shape consisting of triangulated faces
+  //! @param theSectionPlane the plane that intersects the shape
+  //! @param theSectionSegments the container of section segments where the section segments should be added
+  static void crossSection (const TopoDS_Shape& theShape, const gp_Pln& theSectionPlane, Segments& theSectionSegments)
+  {
+    if (!intersects (theShape, theSectionPlane))
+    {
+      return;
+    }
+    for (TopExp_Explorer aFaceIterator (theShape, TopAbs_FACE); aFaceIterator.More(); aFaceIterator.Next())
+    {
+      const TopoDS_Face& aFace = TopoDS::Face (aFaceIterator.Current());
+      crossSection (aFace, theSectionPlane, theSectionSegments);
+    }
+  }
+
+} // namespace
+
+//=======================================================================
+// function : CrossSectionAsShape
+// purpose  :
+//=======================================================================
+TopoDS_Compound StdPrs_SectionLines::CrossSectionAsShape (const TopoDS_Shape& theShape, const gp_Pln& theSectionPlane)
+{
+  Segments aSectionSegments;
+  crossSection (theShape, theSectionPlane, aSectionSegments);
+  TopoDS_Builder  aShapeBuilder;
+  TopoDS_Compound aSectionShape;
+  aShapeBuilder.MakeCompound (aSectionShape);
+  for (Segments::const_iterator aSectionSegmentIterator = aSectionSegments.begin();
+       aSectionSegmentIterator != aSectionSegments.end();
+       ++aSectionSegmentIterator)
+  {
+    const Segment&    aSectionSegment = *aSectionSegmentIterator;
+    const TopoDS_Edge aSectionEdge    = BRepBuilderAPI_MakeEdge (aSectionSegment.first, aSectionSegment.second);
+    aShapeBuilder.Add (aSectionShape, aSectionEdge);
+  }
+  return aSectionShape;
+}
+
+//=======================================================================
+// function : CrossSectionAsArraysOfSegments
+// purpose  :
+//=======================================================================
+Handle (Graphic3d_ArrayOfPrimitives) StdPrs_SectionLines::CrossSectionAsArrayOfSegments (const TopoDS_Shape& theShape,
+                                                                                         const gp_Pln& theSectionPlane)
+{
+  Segments aSectionSegments;
+  crossSection (theShape, theSectionPlane, aSectionSegments);
+  const Standard_Integer aNumberOfSectionEdges = static_cast<Standard_Integer> (2 * aSectionSegments.size());
+  const Handle (Graphic3d_ArrayOfSegments)
+    anArrayOfSectionEdges = new Graphic3d_ArrayOfSegments (aNumberOfSectionEdges);
+  for (Segments::const_iterator aSectionSegmentIterator = aSectionSegments.begin();
+       aSectionSegmentIterator != aSectionSegments.end();
+       ++aSectionSegmentIterator)
+  {
+    const Segment& aSectionSegment = *aSectionSegmentIterator;
+    anArrayOfSectionEdges->AddVertex (aSectionSegment.first);
+    anArrayOfSectionEdges->AddVertex (aSectionSegment.second);
+  }
+  return anArrayOfSectionEdges;
+}
+
+//=======================================================================
+// function : CrossSectionAsArraysOfTrianglesAndEdges
+// purpose  :
+//=======================================================================
+void StdPrs_SectionLines::CrossSectionAsArraysOfTrianglesAndSegments (const TopoDS_Shape& theShape,
+                                                                      const gp_Pln&       theSectionPlane,
+                                                                      Handle (Graphic3d_ArrayOfPrimitives)
+                                                                        & theArrayOfSectionTriangles,
+                                                                      Handle (Graphic3d_ArrayOfPrimitives)
+                                                                        & theArrayOfSectionSegments)
+{
+  (void)theShape;
+  (void)theSectionPlane;
+  (void)theArrayOfSectionTriangles;
+  (void)theArrayOfSectionSegments;
+}
diff --git a/src/StdPrs/StdPrs_SectionLines.hxx b/src/StdPrs/StdPrs_SectionLines.hxx
new file mode 100644 (file)
index 0000000..07c50ab
--- /dev/null
@@ -0,0 +1,57 @@
+// Created on: 2019-04-11
+// Created by: Timur Izmaylov
+// Copyright (c) 2019 OPEN CASCADE SAS
+//
+// This file is part of Open CASCADE Technology software library.
+//
+// This library is free software; you can redistribute it and/or modify it under
+// the terms of the GNU Lesser General Public License version 2.1 as published
+// by the Free Software Foundation, with special exception defined in the file
+// OCCT_LGPL_EXCEPTION.txt. Consult the file LICENSE_LGPL_21.txt included in OCCT
+// distribution for complete text of the license and disclaimer of any warranty.
+//
+// Alternatively, this file may be used under the terms of Open CASCADE
+// commercial license or contractual agreement.
+
+#ifndef _StdPrs_SectionLines_HeaderFile
+#define _StdPrs_SectionLines_HeaderFile
+
+#include <Standard_Handle.hxx>
+
+class Graphic3d_ArrayOfPrimitives;
+class TopoDS_Shape;
+class gp_Pln;
+class TopoDS_Compound;
+
+//! Tool for computing a cross section of a TopoDS_Shape with some plane.
+class StdPrs_SectionLines
+{
+public:
+  //! Computes a cross section of the shape with the plane
+  //! @param theShape the shape consisting of triangulated faces
+  //! @param theSectionPlane the plane that intersects the shape
+  //! @return a resulting shape of a cross section
+  Standard_EXPORT static TopoDS_Compound CrossSectionAsShape (const TopoDS_Shape& theShape,
+                                                              const gp_Pln&       theSectionPlane);
+
+  //! Computes a cross section of the shape with the plane
+  //! @param theShape the shape consisting of triangulated faces
+  //! @param theSectionPlane the plane that intersects the shape
+  //! @return a resulting list of segments
+  Standard_EXPORT static Handle (Graphic3d_ArrayOfPrimitives)
+    CrossSectionAsArrayOfSegments (const TopoDS_Shape& theShape, const gp_Pln& theSectionPlane);
+
+  //! Computes a cross section of the shape with the plane
+  //! @param theShape the shape consisting of triangulated faces
+  //! @param theSectionPlane the plane that intersects the shape
+  //! @param theArrayOfSectionTriangles a resulting list of triangles
+  //! @param theArrayOfSectionSegments a resulting list of line segments
+  Standard_EXPORT static void CrossSectionAsArraysOfTrianglesAndSegments (const TopoDS_Shape& theShape,
+                                                                          const gp_Pln&       theSectionPlane,
+                                                                          Handle (Graphic3d_ArrayOfPrimitives)
+                                                                            & theArrayOfSectionTriangles,
+                                                                          Handle (Graphic3d_ArrayOfPrimitives)
+                                                                            & theArrayOfSectionSegments);
+};
+
+#endif // _StdPrs_SectionLines_HeaderFile
index 7306114afac25eedec40d5b4d92e8d8d2b0f4cf5..87b6cc33a1e385597b0866de9c8e21d7740b3fbc 100644 (file)
 #include <Standard_DefineAlloc.hxx>
 #include <Standard_Handle.hxx>
 
-#include <gp_Ax3.hxx>
-#include <Standard_Real.hxx>
 #include <Standard_Boolean.hxx>
+#include <Standard_Real.hxx>
 #include <gp_Ax1.hxx>
+#include <gp_Ax3.hxx>
+
 class Standard_ConstructionError;
 class gp_Ax3;
 class gp_Pnt;
@@ -33,7 +34,6 @@ class gp_Ax2;
 class gp_Trsf;
 class gp_Vec;
 
-
 //! Describes a plane.
 //! A plane is positioned in space with a coordinate system
 //! (a gp_Ax3 object), such that the plane is defined by the
@@ -55,17 +55,14 @@ class gp_Vec;
 //! Geom_Plane which provides additional functions for
 //! constructing planes and works, in particular, with the
 //! parametric equations of planes
-class gp_Pln 
+class gp_Pln
 {
 public:
-
   DEFINE_STANDARD_ALLOC
 
-  
   //! Creates a plane coincident with OXY plane of the
   //! reference coordinate system.
-    gp_Pln();
-  
+  gp_Pln();
 
   //! The coordinate system of the plane is defined with the axis
   //! placement A3.
@@ -73,83 +70,84 @@ public:
   //! The "Location" of A3 defines the location (origin) of the plane.
   //! The "XDirection" and "YDirection" of A3 define the "XAxis" and
   //! the "YAxis" of the plane used to parametrize the plane.
-    gp_Pln(const gp_Ax3& A3);
-  
+  gp_Pln (const gp_Ax3& A3);
 
   //! Creates a plane with the  "Location" point <P>
   //! and the normal direction <V>.
-  Standard_EXPORT gp_Pln(const gp_Pnt& P, const gp_Dir& V);
-  
+  Standard_EXPORT gp_Pln (const gp_Pnt& P, const gp_Dir& V);
 
   //! Creates a plane from its cartesian equation :
   //! A * X + B * Y + C * Z + D = 0.0
   //! Raises ConstructionError if Sqrt (A*A + B*B + C*C) <= Resolution from gp.
-  Standard_EXPORT gp_Pln(const Standard_Real A, const Standard_Real B, const Standard_Real C, const Standard_Real D);
-  
+  Standard_EXPORT gp_Pln (const Standard_Real A, const Standard_Real B, const Standard_Real C, const Standard_Real D);
 
   //! Returns the coefficients of the plane's cartesian equation :
   //! A * X + B * Y + C * Z + D = 0.
-    void Coefficients (Standard_Real& A, Standard_Real& B, Standard_Real& C, Standard_Real& D) const;
-  
+  void Coefficients (Standard_Real& A, Standard_Real& B, Standard_Real& C, Standard_Real& D) const;
+
   //! Modifies this plane, by redefining its local coordinate system so that
   //! -   its origin and "main Direction" become those of the
   //! axis A1 (the "X Direction" and "Y Direction" are then recomputed).
   //! Raises ConstructionError if the A1 is parallel to the "XAxis" of the plane.
-    void SetAxis (const gp_Ax1& A1);
-  
+  void SetAxis (const gp_Ax1& A1);
+
   //! Changes the origin of the plane.
-    void SetLocation (const gp_Pnt& Loc);
-  
+  void SetLocation (const gp_Pnt& Loc);
+
   //! Changes the local coordinate system of the plane.
-    void SetPosition (const gp_Ax3& A3);
-  
+  void SetPosition (const gp_Ax3& A3);
+
   //! Reverses the   U   parametrization of   the  plane
   //! reversing the XAxis.
-    void UReverse();
-  
+  void UReverse();
+
   //! Reverses the   V   parametrization of   the  plane
   //! reversing the YAxis.
-    void VReverse();
-  
+  void VReverse();
+
   //! returns true if the Ax3 is right handed.
-    Standard_Boolean Direct() const;
-  
+  Standard_Boolean Direct() const;
+
   //! Returns the plane's normal Axis.
-    const gp_Ax1& Axis() const;
-  
+  const gp_Ax1& Axis() const;
+
   //! Returns the plane's location (origin).
-    const gp_Pnt& Location() const;
-  
+  const gp_Pnt& Location() const;
+
   //! Returns the local coordinate system of the plane .
-    const gp_Ax3& Position() const;
-  
-  //! Computes the distance between <me> and the point <P>.
-    Standard_Real Distance (const gp_Pnt& P) const;
-  
+  const gp_Ax3& Position() const;
+
+  //! Computes the signed distance between this plane and the given point
+  //! @param thePoint the point the signed distance to which from this plane is computed
+  //! @return the signed distance between the given point and this plane
+  Standard_Real SignedDistance (const gp_Pnt& thePoint) const;
+
+  //! Computes the distance between this plane and the given point
+  //! @param thePoint the point the distance to which from this plane is computed
+  //! @return the distance between the given point and this plane
+  Standard_Real Distance (const gp_Pnt& thePoint) const;
+
   //! Computes the distance between <me> and the line <L>.
-    Standard_Real Distance (const gp_Lin& L) const;
-  
+  Standard_Real Distance (const gp_Lin& L) const;
+
   //! Computes the distance between two planes.
-    Standard_Real Distance (const gp_Pln& Other) const;
-  
+  Standard_Real Distance (const gp_Pln& Other) const;
 
   //! Computes the square distance between <me> and the point <P>.
-    Standard_Real SquareDistance (const gp_Pnt& P) const;
-  
+  Standard_Real SquareDistance (const gp_Pnt& P) const;
 
   //! Computes the square distance between <me> and the line <L>.
-    Standard_Real SquareDistance (const gp_Lin& L) const;
-  
+  Standard_Real SquareDistance (const gp_Lin& L) const;
 
   //! Computes the square distance between two planes.
-    Standard_Real SquareDistance (const gp_Pln& Other) const;
-  
+  Standard_Real SquareDistance (const gp_Pln& Other) const;
+
   //! Returns the X axis of the plane.
-    gp_Ax1 XAxis() const;
-  
+  gp_Ax1 XAxis() const;
+
   //! Returns the Y axis  of the plane.
-    gp_Ax1 YAxis() const;
-  
+  gp_Ax1 YAxis() const;
+
   //! Returns true if this plane contains the point P. This means that
   //! -   the distance between point P and this plane is less
   //! than or equal to LinearTolerance, or
@@ -158,8 +156,8 @@ public:
   //! AngularTolerance, and the distance between the origin
   //! of line L and this plane is less than or equal to
   //! LinearTolerance.
-    Standard_Boolean Contains (const gp_Pnt& P, const Standard_Real LinearTolerance) const;
-  
+  Standard_Boolean Contains (const gp_Pnt& P, const Standard_Real LinearTolerance) const;
+
   //! Returns true if this plane contains the line L. This means that
   //! -   the distance between point P and this plane is less
   //! than or equal to LinearTolerance, or
@@ -168,10 +166,11 @@ public:
   //! AngularTolerance, and the distance between the origin
   //! of line L and this plane is less than or equal to
   //! LinearTolerance.
-    Standard_Boolean Contains (const gp_Lin& L, const Standard_Real LinearTolerance, const Standard_Real AngularTolerance) const;
-  
+  Standard_Boolean Contains (const gp_Lin&       L,
+                             const Standard_Real LinearTolerance,
+                             const Standard_Real AngularTolerance) const;
+
   Standard_EXPORT void Mirror (const gp_Pnt& P);
-  
 
   //! Performs the symmetrical transformation of a plane with respect
   //! to the point <P> which is the center of the symmetry
@@ -179,9 +178,9 @@ public:
   //! The normal direction to the plane is not changed.
   //! The "XAxis" and the "YAxis" are reversed.
   Standard_EXPORT Standard_NODISCARD gp_Pln Mirrored (const gp_Pnt& P) const;
-  
+
   Standard_EXPORT void Mirror (const gp_Ax1& A1);
-  
+
   //! Performs   the symmetrical transformation  of a
   //! plane with respect to an axis placement  which is the axis
   //! of  the symmetry.  The  transformation is performed on the
@@ -191,9 +190,9 @@ public:
   //! if  the  initial plane was right  handed,  else  it is the
   //! opposite.
   Standard_EXPORT Standard_NODISCARD gp_Pln Mirrored (const gp_Ax1& A1) const;
-  
+
   Standard_EXPORT void Mirror (const gp_Ax2& A2);
-  
+
   //! Performs the  symmetrical transformation  of  a
   //! plane    with respect to    an axis  placement.   The axis
   //! placement  <A2> locates the plane  of  the symmetry.   The
@@ -203,66 +202,42 @@ public:
   //! and the "YDirection"  after  transformation if the initial
   //! plane was right handed, else it is the opposite.
   Standard_EXPORT Standard_NODISCARD gp_Pln Mirrored (const gp_Ax2& A2) const;
-  
-    void Rotate (const gp_Ax1& A1, const Standard_Real Ang);
-  
+
+  void Rotate (const gp_Ax1& A1, const Standard_Real Ang);
 
   //! rotates a plane. A1 is the axis of the rotation.
   //! Ang is the angular value of the rotation in radians.
-    Standard_NODISCARD gp_Pln Rotated (const gp_Ax1& A1, const Standard_Real Ang) const;
-  
-    void Scale (const gp_Pnt& P, const Standard_Real S);
-  
+  Standard_NODISCARD gp_Pln Rotated (const gp_Ax1& A1, const Standard_Real Ang) const;
+
+  void Scale (const gp_Pnt& P, const Standard_Real S);
 
   //! Scales a plane. S is the scaling value.
-    Standard_NODISCARD gp_Pln Scaled (const gp_Pnt& P, const Standard_Real S) const;
-  
-    void Transform (const gp_Trsf& T);
-  
+  Standard_NODISCARD gp_Pln Scaled (const gp_Pnt& P, const Standard_Real S) const;
+
+  void Transform (const gp_Trsf& T);
 
   //! Transforms a plane with the transformation T from class Trsf.
   //! The transformation is performed on the "Location"
   //! point, on the "XAxis" and the "YAxis".
   //! The resulting normal direction is the cross product between
   //! the "XDirection" and the "YDirection" after transformation.
-    Standard_NODISCARD gp_Pln Transformed (const gp_Trsf& T) const;
-  
-    void Translate (const gp_Vec& V);
-  
+  Standard_NODISCARD gp_Pln Transformed (const gp_Trsf& T) const;
+
+  void Translate (const gp_Vec& V);
 
   //! Translates a plane in the direction of the vector V.
   //! The magnitude of the translation is the vector's magnitude.
-    Standard_NODISCARD gp_Pln Translated (const gp_Vec& V) const;
-  
-    void Translate (const gp_Pnt& P1, const gp_Pnt& P2);
-  
-
-  //! Translates a plane from the point P1 to the point P2.
-    Standard_NODISCARD gp_Pln Translated (const gp_Pnt& P1, const gp_Pnt& P2) const;
-
-
-
-
-protected:
-
-
+  Standard_NODISCARD gp_Pln Translated (const gp_Vec& V) const;
 
+  void Translate (const gp_Pnt& P1, const gp_Pnt& P2);
 
+  //! Translates a plane from the point P1 to the point P2.
+  Standard_NODISCARD gp_Pln Translated (const gp_Pnt& P1, const gp_Pnt& P2) const;
 
 private:
-
-
-
   gp_Ax3 pos;
-
-
 };
 
-
 #include <gp_Pln.lxx>
 
-
-
-
-
 #endif // _gp_Pln_HeaderFile
index c2b056fac8fb1d5f62477e20e3d24225ddfa6f3f..fdb18f0b2851865211bf1ff62eb3be117a5b2cd3 100644 (file)
@@ -69,15 +69,20 @@ inline const gp_Pnt& gp_Pln::Location() const
 inline   const gp_Ax3& gp_Pln::Position() const
 { return pos; }
 
-inline Standard_Real gp_Pln::Distance(const gp_Pnt& P) const
+inline Standard_Real gp_Pln::SignedDistance (const gp_Pnt& thePoint) const
 {
-  const gp_Pnt& loc = pos.Location ();
-  const gp_Dir& dir = pos.Direction();
-  Standard_Real D = (dir.X() * (P.X() - loc.X()) +
-                    dir.Y() * (P.Y() - loc.Y()) +
-                    dir.Z() * (P.Z() - loc.Z()));
-  if (D < 0) D = - D;
-  return D;
+  const gp_Pnt&       aPlaneLocation      = pos.Location();
+  const gp_Dir&       aPlaneAxisDirection = pos.Direction();
+  const Standard_Real aSignedDistance     = aPlaneAxisDirection.X() * (thePoint.X() - aPlaneLocation.X())
+                                        + aPlaneAxisDirection.Y() * (thePoint.Y() - aPlaneLocation.Y())
+                                        + aPlaneAxisDirection.Z() * (thePoint.Z() - aPlaneLocation.Z());
+  return aSignedDistance;
+}
+
+inline Standard_Real gp_Pln::Distance (const gp_Pnt& thePoint) const
+{
+  const Standard_Real aSignedDistance = SignedDistance (thePoint);
+  return Abs( aSignedDistance);
 }
 
 inline Standard_Real gp_Pln::Distance (const gp_Lin& L)  const