1
0
mirror of https://git.dev.opencascade.org/repos/occt.git synced 2025-08-09 13:22:24 +03:00

0032338: Visualization - provide straightforward interface for ray-picking

This commit is contained in:
osa
2021-05-18 10:07:59 +03:00
committed by bugmaster
parent e1eb39d29d
commit 0461e7fd03
12 changed files with 1011 additions and 1 deletions

View File

@@ -4,6 +4,8 @@ SelectMgr_AndFilter.cxx
SelectMgr_AndFilter.hxx
SelectMgr_AndOrFilter.cxx
SelectMgr_AndOrFilter.hxx
SelectMgr_AxisIntersector.cxx
SelectMgr_AxisIntersector.hxx
SelectMgr_BaseIntersector.cxx
SelectMgr_BaseIntersector.hxx
SelectMgr_BaseFrustum.cxx

View File

@@ -0,0 +1,541 @@
// Copyright (c) 2021 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 <SelectMgr_AxisIntersector.hxx>
#include <Bnd_Range.hxx>
#include <BVH_Tools.hxx>
#include <Precision.hxx>
#include <SelectBasics_PickResult.hxx>
#include <SelectMgr_ViewClipRange.hxx>
// =======================================================================
// function : Constructor
// purpose :
// =======================================================================
SelectMgr_AxisIntersector::SelectMgr_AxisIntersector()
{
}
// =======================================================================
// function : Init
// purpose :
// =======================================================================
void SelectMgr_AxisIntersector::Init (const gp_Ax1& theAxis)
{
mySelectionType = SelectMgr_SelectionType_Point;
myAxis = theAxis;
}
// =======================================================================
// function : Build
// purpose :
// =======================================================================
void SelectMgr_AxisIntersector::Build()
{
}
// =======================================================================
// function : ScaleAndTransform
// purpose :
// =======================================================================
Handle(SelectMgr_BaseIntersector) SelectMgr_AxisIntersector::ScaleAndTransform (const Standard_Integer theScaleFactor,
const gp_GTrsf& theTrsf,
const Handle(SelectMgr_FrustumBuilder)& theBuilder) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::ScaleAndTransform() should be called after selection axis initialization");
(void )theScaleFactor;
(void )theBuilder;
if (theTrsf.Form() == gp_Identity)
{
return new SelectMgr_AxisIntersector();
}
gp_Pnt aTransformedLoc = myAxis.Location();
theTrsf.Transforms (aTransformedLoc.ChangeCoord());
gp_XYZ aTransformedDir = myAxis.Direction().XYZ();
gp_GTrsf aTrsf = theTrsf;
aTrsf.SetTranslationPart (gp_XYZ(0., 0., 0.));
aTrsf.Transforms (aTransformedDir);
Handle(SelectMgr_AxisIntersector) aRes = new SelectMgr_AxisIntersector();
aRes->myAxis = gp_Ax1(aTransformedLoc, gp_Dir(aTransformedDir));
aRes->mySelectionType = mySelectionType;
return aRes;
}
// =======================================================================
// function : hasIntersection
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::hasIntersection (const SelectMgr_Vec3& theBoxMin,
const SelectMgr_Vec3& theBoxMax,
Standard_Real& theTimeEnter,
Standard_Real& theTimeLeave) const
{
const gp_Pnt& anAxisLoc = myAxis.Location();
const gp_Dir& anAxisDir = myAxis.Direction();
BVH_Ray<Standard_Real, 3> aRay(SelectMgr_Vec3(anAxisLoc.X(), anAxisLoc.Y(), anAxisLoc.Z()),
SelectMgr_Vec3(anAxisDir.X(), anAxisDir.Y(), anAxisDir.Z()));
if (!BVH_Tools<Standard_Real, 3>::RayBoxIntersection (aRay, theBoxMin, theBoxMax, theTimeEnter, theTimeLeave))
{
return Standard_False;
}
return Standard_True;
}
// =======================================================================
// function : hasIntersection
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::hasIntersection (const gp_Pnt& thePnt,
Standard_Real& theDepth) const
{
const gp_Pnt& anAxisLoc = myAxis.Location();
const gp_Dir& anAxisDir = myAxis.Direction();
// Check that vectors are co-directed (thePnt lies on this axis)
gp_Dir aDirToPnt(thePnt.XYZ() - anAxisLoc.XYZ());
if (!anAxisDir.IsEqual (aDirToPnt, Precision::Angular()))
{
return Standard_False;
}
theDepth = anAxisLoc.Distance (thePnt);
return Standard_True;
}
// =======================================================================
// function : raySegmentDistance
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::raySegmentDistance (const gp_Pnt& theSegPnt1,
const gp_Pnt& theSegPnt2,
SelectBasics_PickResult& thePickResult) const
{
gp_XYZ anU = theSegPnt2.XYZ() - theSegPnt1.XYZ();
gp_XYZ aV = myAxis.Direction().XYZ();
gp_XYZ aW = theSegPnt1.XYZ() - myAxis.Location().XYZ();
gp_XYZ anUVNormVec = aV.Crossed (anU);
gp_XYZ anUWNormVec = aW.Crossed (anU);
if (anUVNormVec.Modulus() <= Precision::Confusion() ||
anUWNormVec.Modulus() <= Precision::Confusion())
{
// Lines have no intersection
thePickResult.Invalidate();
return false;
}
Standard_Real aParam = anUWNormVec.Dot (anUVNormVec) / anUVNormVec.SquareModulus();
if (aParam < 0.0)
{
// Intersection is out of axis start point
thePickResult.Invalidate();
return false;
}
gp_XYZ anIntersectPnt = myAxis.Location().XYZ() + aV * aParam;
if ((anIntersectPnt - theSegPnt1.XYZ()).SquareModulus() +
(anIntersectPnt - theSegPnt2.XYZ()).SquareModulus() >
anU.SquareModulus() + Precision::Confusion())
{
// Intersection point doesn't lie on the segment
thePickResult.Invalidate();
return false;
}
thePickResult.SetDepth (myAxis.Location().Distance (anIntersectPnt));
thePickResult.SetPickedPoint (anIntersectPnt);
return true;
}
// =======================================================================
// function : rayPlaneIntersection
// purpose :
// =======================================================================
bool SelectMgr_AxisIntersector::rayPlaneIntersection (const gp_Vec& thePlane,
const gp_Pnt& thePntOnPlane,
SelectBasics_PickResult& thePickResult) const
{
gp_XYZ anU = myAxis.Direction().XYZ();
gp_XYZ aW = myAxis.Location().XYZ() - thePntOnPlane.XYZ();
Standard_Real aD = thePlane.Dot (anU);
Standard_Real aN = -thePlane.Dot (aW);
if (Abs (aD) < Precision::Confusion())
{
thePickResult.Invalidate();
return false;
}
Standard_Real aParam = aN / aD;
if (aParam < 0.0)
{
thePickResult.Invalidate();
return false;
}
gp_Pnt aClosestPnt = myAxis.Location().XYZ() + anU * aParam;
thePickResult.SetDepth (myAxis.Location().Distance (aClosestPnt));
thePickResult.SetPickedPoint (aClosestPnt);
return true;
}
// =======================================================================
// function : Overlaps
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::Overlaps (const SelectMgr_Vec3& theBoxMin,
const SelectMgr_Vec3& theBoxMax,
Standard_Boolean* theInside) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::Overlaps() should be called after selection axis initialization");
(void )theInside;
Standard_Real aTimeEnter, aTimeLeave;
if (!hasIntersection (theBoxMin, theBoxMax, aTimeEnter, aTimeLeave))
{
return Standard_False;
}
if (theInside != NULL)
{
*theInside &= (aTimeEnter >= 0.0);
}
return Standard_True;
}
// =======================================================================
// function : Overlaps
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::Overlaps (const SelectMgr_Vec3& theBoxMin,
const SelectMgr_Vec3& theBoxMax,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::Overlaps() should be called after selection axis initialization");
Standard_Real aTimeEnter, aTimeLeave;
if (!hasIntersection (theBoxMin, theBoxMax, aTimeEnter, aTimeLeave))
{
return Standard_False;
}
Standard_Real aDepth = 0.0;
Bnd_Range aRange(Max (aTimeEnter, 0.0), aTimeLeave);
aRange.GetMin (aDepth);
if (!theClipRange.GetNearestDepth (aRange, aDepth))
{
return Standard_False;
}
thePickResult.SetDepth (aDepth);
return Standard_True;
}
// =======================================================================
// function : Overlaps
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::Overlaps (const gp_Pnt& thePnt,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::Overlaps() should be called after selection axis initialization");
Standard_Real aDepth = 0.0;
if (!hasIntersection (thePnt, aDepth))
{
return Standard_False;
}
thePickResult.SetDepth (aDepth);
thePickResult.SetPickedPoint (thePnt);
return !theClipRange.IsClipped (thePickResult.Depth());
}
// =======================================================================
// function : Overlaps
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::Overlaps (const gp_Pnt& thePnt) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::Overlaps() should be called after selection axis initialization");
Standard_Real aDepth = 0.0;
return hasIntersection (thePnt, aDepth);
}
// =======================================================================
// function : Overlaps
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::Overlaps (const gp_Pnt& thePnt1,
const gp_Pnt& thePnt2,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::Overlaps() should be called after selection axis initialization");
if (!raySegmentDistance (thePnt1, thePnt2, thePickResult))
{
return Standard_False;
}
return !theClipRange.IsClipped (thePickResult.Depth());
}
// =======================================================================
// function : Overlaps
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::Overlaps (const TColgp_Array1OfPnt& theArrayOfPnts,
Select3D_TypeOfSensitivity theSensType,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::Overlaps() should be called after selection axis initialization");
if (theSensType == Select3D_TOS_BOUNDARY)
{
Standard_Integer aMatchingSegmentsNb = -1;
SelectBasics_PickResult aPickResult;
thePickResult.Invalidate();
const Standard_Integer aLower = theArrayOfPnts.Lower();
const Standard_Integer anUpper = theArrayOfPnts.Upper();
for (Standard_Integer aPntIter = aLower; aPntIter <= anUpper; ++aPntIter)
{
const gp_Pnt& aStartPnt = theArrayOfPnts.Value (aPntIter);
const gp_Pnt& aEndPnt = theArrayOfPnts.Value (aPntIter == anUpper ? aLower : (aPntIter + 1));
if (raySegmentDistance (aStartPnt, aEndPnt, aPickResult))
{
aMatchingSegmentsNb++;
thePickResult = SelectBasics_PickResult::Min (thePickResult, aPickResult);
}
}
if (aMatchingSegmentsNb == -1)
{
return Standard_False;
}
}
else if (theSensType == Select3D_TOS_INTERIOR)
{
Standard_Integer aStartIdx = theArrayOfPnts.Lower();
const gp_XYZ& aPnt1 = theArrayOfPnts.Value (aStartIdx).XYZ();
const gp_XYZ& aPnt2 = theArrayOfPnts.Value (aStartIdx + 1).XYZ();
const gp_XYZ& aPnt3 = theArrayOfPnts.Value (aStartIdx + 2).XYZ();
const gp_XYZ aVec1 = aPnt1 - aPnt2;
const gp_XYZ aVec2 = aPnt3 - aPnt2;
gp_Vec aPolyNorm = aVec2.Crossed (aVec1);
if (aPolyNorm.Magnitude() <= Precision::Confusion())
{
// treat degenerated polygon as point
return Overlaps (theArrayOfPnts.First(), theClipRange, thePickResult);
}
else if (!rayPlaneIntersection (aPolyNorm, theArrayOfPnts.First(), thePickResult))
{
return Standard_False;
}
}
return !theClipRange.IsClipped (thePickResult.Depth());
}
// =======================================================================
// function : Overlaps
// purpose :
// =======================================================================
Standard_Boolean SelectMgr_AxisIntersector::Overlaps (const gp_Pnt& thePnt1,
const gp_Pnt& thePnt2,
const gp_Pnt& thePnt3,
Select3D_TypeOfSensitivity theSensType,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::Overlaps() should be called after selection axis initialization");
if (theSensType == Select3D_TOS_BOUNDARY)
{
const gp_Pnt aPntsArrayBuf[4] = { thePnt1, thePnt2, thePnt3, thePnt1 };
const TColgp_Array1OfPnt aPntsArray (aPntsArrayBuf[0], 1, 4);
return Overlaps (aPntsArray, Select3D_TOS_BOUNDARY, theClipRange, thePickResult);
}
else if (theSensType == Select3D_TOS_INTERIOR)
{
gp_Vec aTriangleNormal (gp_XYZ (RealLast(), RealLast(), RealLast()));
const gp_XYZ aTrEdges[3] = { thePnt2.XYZ() - thePnt1.XYZ(),
thePnt3.XYZ() - thePnt2.XYZ(),
thePnt1.XYZ() - thePnt3.XYZ() };
aTriangleNormal = aTrEdges[2].Crossed (aTrEdges[0]);
if (aTriangleNormal.SquareMagnitude() < gp::Resolution())
{
// consider degenerated triangle as point or segment
return aTrEdges[0].SquareModulus() > gp::Resolution()
? Overlaps (thePnt1, thePnt2, theClipRange, thePickResult)
: (aTrEdges[1].SquareModulus() > gp::Resolution()
? Overlaps (thePnt2, thePnt3, theClipRange, thePickResult)
: Overlaps (thePnt1, theClipRange, thePickResult));
}
const gp_Pnt aPnts[3] = {thePnt1, thePnt2, thePnt3};
const Standard_Real anAlpha = aTriangleNormal.XYZ().Dot (myAxis.Direction().XYZ());
if (Abs (anAlpha) < gp::Resolution())
{
// handle the case when triangle normal and selecting frustum direction are orthogonal
SelectBasics_PickResult aPickResult;
thePickResult.Invalidate();
for (Standard_Integer anEdgeIter = 0; anEdgeIter < 3; ++anEdgeIter)
{
const gp_Pnt& aStartPnt = aPnts[anEdgeIter];
const gp_Pnt& anEndPnt = aPnts[anEdgeIter < 2 ? anEdgeIter + 1 : 0];
if (raySegmentDistance (aStartPnt, anEndPnt, aPickResult))
{
thePickResult = SelectBasics_PickResult::Min (thePickResult, aPickResult);
}
}
thePickResult.SetSurfaceNormal (aTriangleNormal);
return !theClipRange.IsClipped (thePickResult.Depth());
}
// check if intersection point belongs to triangle's interior part
const gp_XYZ anEdge = (thePnt1.XYZ() - myAxis.Location().XYZ()) * (1.0 / anAlpha);
const Standard_Real aTime = aTriangleNormal.Dot (anEdge);
const gp_XYZ aVec = myAxis.Direction().XYZ().Crossed (anEdge);
const Standard_Real anU = aVec.Dot (aTrEdges[2]);
const Standard_Real aV = aVec.Dot (aTrEdges[0]);
const Standard_Boolean isInterior = (aTime >= 0.0) && (anU >= 0.0) && (aV >= 0.0) && (anU + aV <= 1.0);
const gp_Pnt aPtOnPlane = myAxis.Location().XYZ() + myAxis.Direction().XYZ() * aTime;
if (isInterior)
{
thePickResult.SetDepth (myAxis.Location().Distance (aPtOnPlane));
thePickResult.SetPickedPoint (aPtOnPlane);
thePickResult.SetSurfaceNormal (aTriangleNormal);
return !theClipRange.IsClipped (thePickResult.Depth());
}
Standard_Real aMinDist = RealLast();
Standard_Integer aNearestEdgeIdx1 = -1;
for (Standard_Integer anEdgeIdx = 0; anEdgeIdx < 3; ++anEdgeIdx)
{
gp_XYZ aW = aPtOnPlane.XYZ() - aPnts[anEdgeIdx].XYZ();
Standard_Real aCoef = aTrEdges[anEdgeIdx].Dot (aW) / aTrEdges[anEdgeIdx].Dot (aTrEdges[anEdgeIdx]);
Standard_Real aDist = aPtOnPlane.Distance (aPnts[anEdgeIdx].XYZ() + aCoef * aTrEdges[anEdgeIdx]);
if (aDist < aMinDist)
{
aMinDist = aDist;
aNearestEdgeIdx1 = anEdgeIdx;
}
}
Standard_Integer aNearestEdgeIdx2 = (aNearestEdgeIdx1 + 1) % 3;
const gp_Vec aVec12 (aPnts[aNearestEdgeIdx1], aPnts[aNearestEdgeIdx2]);
if (aVec12.SquareMagnitude() > gp::Resolution()
&& myAxis.Direction().IsParallel (aVec12, Precision::Angular()))
{
aNearestEdgeIdx2 = aNearestEdgeIdx1 == 0 ? 2 : aNearestEdgeIdx1 - 1;
}
if (raySegmentDistance (aPnts[aNearestEdgeIdx1], aPnts[aNearestEdgeIdx2], thePickResult))
{
thePickResult.SetSurfaceNormal (aTriangleNormal);
}
}
return !theClipRange.IsClipped (thePickResult.Depth());
}
//=======================================================================
// function : GetNearPnt
// purpose :
//=======================================================================
const gp_Pnt& SelectMgr_AxisIntersector::GetNearPnt() const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::GetNearPnt() should be called after selection axis initialization");
return myAxis.Location();
}
//=======================================================================
// function : GetFarPnt
// purpose :
//=======================================================================
const gp_Pnt& SelectMgr_AxisIntersector::GetFarPnt() const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::GetFarPnt() should be called after selection axis initialization");
static gp_Pnt anInfPnt(RealLast(), RealLast(), RealLast());
return anInfPnt;
}
//=======================================================================
// function : GetViewRayDirection
// purpose :
//=======================================================================
const gp_Dir& SelectMgr_AxisIntersector::GetViewRayDirection() const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::GetViewRayDirection() should be called after selection axis initialization");
return myAxis.Direction();
}
// =======================================================================
// function : DistToGeometryCenter
// purpose :
// =======================================================================
Standard_Real SelectMgr_AxisIntersector::DistToGeometryCenter (const gp_Pnt& theCOG) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::DistToGeometryCenter() should be called after selection axis initialization");
return theCOG.Distance (myAxis.Location());
}
// =======================================================================
// function : DetectedPoint
// purpose :
// =======================================================================
gp_Pnt SelectMgr_AxisIntersector::DetectedPoint (const Standard_Real theDepth) const
{
Standard_ASSERT_RAISE(mySelectionType == SelectMgr_SelectionType_Point,
"Error! SelectMgr_AxisIntersector::DetectedPoint() should be called after selection axis initialization");
return myAxis.Location().XYZ() + myAxis.Direction().XYZ() * theDepth;
}
//=======================================================================
//function : DumpJson
//purpose :
//=======================================================================
void SelectMgr_AxisIntersector::DumpJson (Standard_OStream& theOStream, Standard_Integer theDepth) const
{
OCCT_DUMP_CLASS_BEGIN (theOStream, SelectMgr_AxisIntersector)
OCCT_DUMP_BASE_CLASS (theOStream, theDepth, SelectMgr_BaseIntersector)
OCCT_DUMP_FIELD_VALUES_DUMPED (theOStream, theDepth, &myAxis)
}

View File

@@ -0,0 +1,147 @@
// Copyright (c) 2021 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 _SelectMgr_AxisIntersector_HeaderFile
#define _SelectMgr_AxisIntersector_HeaderFile
#include <SelectMgr_BaseIntersector.hxx>
#include <gp_Ax1.hxx>
//! This class contains representation of selecting axis, created in case of point selection
//! and algorithms for overlap detection between this axis and sensitive entities.
class SelectMgr_AxisIntersector : public SelectMgr_BaseIntersector
{
public:
//! Empty constructor
Standard_EXPORT SelectMgr_AxisIntersector();
//! Destructor
virtual ~SelectMgr_AxisIntersector() {}
//! Initializes selecting axis according to the input one
Standard_EXPORT void Init (const gp_Ax1& theAxis);
//! Builds axis according to internal parameters.
//! NOTE: it should be called after Init() method
Standard_EXPORT virtual void Build() Standard_OVERRIDE;
//! IMPORTANT: Scaling doesn't make sense for this intersector.
//! Returns a copy of the intersector transformed using the matrix given.
//! Builder is an optional argument that represents corresponding settings for re-constructing transformed
//! frustum from scratch. Can be null if reconstruction is not expected furthermore.
Standard_EXPORT virtual Handle(SelectMgr_BaseIntersector)
ScaleAndTransform (const Standard_Integer theScaleFactor,
const gp_GTrsf& theTrsf,
const Handle(SelectMgr_FrustumBuilder)& theBuilder = Handle(SelectMgr_FrustumBuilder)()) const Standard_OVERRIDE;
public:
//! Intersection test between defined axis and given axis-aligned box
Standard_EXPORT virtual Standard_Boolean Overlaps (const SelectMgr_Vec3& theBoxMin,
const SelectMgr_Vec3& theBoxMax,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const Standard_OVERRIDE;
//! Returns true if selecting axis intersects axis-aligned bounding box
//! with minimum corner at point theMinPt and maximum at point theMaxPt
Standard_EXPORT virtual Standard_Boolean Overlaps (const SelectMgr_Vec3& theBoxMin,
const SelectMgr_Vec3& theBoxMax,
Standard_Boolean* theInside) const Standard_OVERRIDE;
//! Intersection test between defined axis and given point
Standard_EXPORT virtual Standard_Boolean Overlaps (const gp_Pnt& thePnt,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const Standard_OVERRIDE;
//! Intersection test between defined axis and given point
Standard_EXPORT virtual Standard_Boolean Overlaps (const gp_Pnt& thePnt) const Standard_OVERRIDE;
//! Intersection test between defined axis and given ordered set of points,
//! representing line segments. The test may be considered of interior part or
//! boundary line defined by segments depending on given sensitivity type
Standard_EXPORT virtual Standard_Boolean Overlaps (const TColgp_Array1OfPnt& theArrayOfPnts,
Select3D_TypeOfSensitivity theSensType,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const Standard_OVERRIDE;
//! Checks if selecting axis intersects line segment
Standard_EXPORT virtual Standard_Boolean Overlaps (const gp_Pnt& thePnt1,
const gp_Pnt& thePnt2,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const Standard_OVERRIDE;
//! Intersection test between defined axis and given triangle. The test may
//! be considered of interior part or boundary line defined by triangle vertices
//! depending on given sensitivity type
Standard_EXPORT virtual Standard_Boolean Overlaps (const gp_Pnt& thePnt1,
const gp_Pnt& thePnt2,
const gp_Pnt& thePnt3,
Select3D_TypeOfSensitivity theSensType,
const SelectMgr_ViewClipRange& theClipRange,
SelectBasics_PickResult& thePickResult) const Standard_OVERRIDE;
public:
//! Measures distance between start axis point and given point theCOG.
Standard_EXPORT virtual Standard_Real DistToGeometryCenter (const gp_Pnt& theCOG) const Standard_OVERRIDE;
//! Calculates the point on a axis ray that was detected during the run of selection algo by given depth
Standard_EXPORT virtual gp_Pnt DetectedPoint (const Standard_Real theDepth) const Standard_OVERRIDE;
//! Returns near point along axis.
Standard_EXPORT virtual const gp_Pnt& GetNearPnt() const Standard_OVERRIDE;
//! Returns far point along axis (infinite).
Standard_EXPORT virtual const gp_Pnt& GetFarPnt() const Standard_OVERRIDE;
//! Returns axis direction.
Standard_EXPORT virtual const gp_Dir& GetViewRayDirection() const Standard_OVERRIDE;
//! Dumps the content of me into the stream
Standard_EXPORT virtual void DumpJson (Standard_OStream& theOStream, Standard_Integer theDepth = -1) const Standard_OVERRIDE;
protected:
//! Returns true if selecting axis intersects axis-aligned bounding box
//! with minimum corner at point theBoxMin and maximum at point theBoxMax.
//! Also returns enter and leave time of axis-box intersection.
Standard_EXPORT Standard_Boolean hasIntersection (const SelectMgr_Vec3& theBoxMin,
const SelectMgr_Vec3& theBoxMax,
Standard_Real& theTimeEnter,
Standard_Real& theTimeLeave) const;
//! Returns true if selecting axis intersects point.
//! Also returns time of axis-point intersection.
Standard_EXPORT Standard_Boolean hasIntersection (const gp_Pnt& thePnt,
Standard_Real& theDepth) const;
//! Returns true if selecting axis intersects segment.
//! Also saves time of axis-segment intersection and intersection point as pick result.
Standard_EXPORT Standard_Boolean raySegmentDistance (const gp_Pnt& theSegPnt1,
const gp_Pnt& theSegPnt2,
SelectBasics_PickResult& thePickResult) const;
//! Returns true if selecting axis intersects plane.
//! Also saves time of axis-plane intersection and intersection point as pick result.
Standard_EXPORT Standard_Boolean rayPlaneIntersection (const gp_Vec& thePlane,
const gp_Pnt& thePntOnPlane,
SelectBasics_PickResult& thePickResult) const;
private:
gp_Ax1 myAxis;
};
#endif // _SelectMgr_AxisIntersector_HeaderFile

View File

@@ -16,6 +16,7 @@
#include <SelectMgr_SelectingVolumeManager.hxx>
#include <Graphic3d_SequenceOfHClipPlane.hxx>
#include <SelectMgr_AxisIntersector.hxx>
#include <SelectMgr_RectangularFrustum.hxx>
#include <SelectMgr_TriangularFrustumSet.hxx>
@@ -261,6 +262,21 @@ void SelectMgr_SelectingVolumeManager::InitPolylineSelectingVolume (const TColgp
aPolylineVolume->SetAllowOverlapDetection (IsOverlapAllowed());
}
//=======================================================================
// function : InitAxisSelectingVolume
// purpose :
//=======================================================================
void SelectMgr_SelectingVolumeManager::InitAxisSelectingVolume (const gp_Ax1& theAxis)
{
Handle(SelectMgr_AxisIntersector) anAxisVolume = Handle(SelectMgr_AxisIntersector)::DownCast(myActiveSelectingVolume);
if (anAxisVolume.IsNull())
{
anAxisVolume = new SelectMgr_AxisIntersector();
}
anAxisVolume->Init (theAxis);
myActiveSelectingVolume = anAxisVolume;
}
//=======================================================================
// function : InitSelectingVolume
// purpose :

View File

@@ -51,6 +51,9 @@ public:
//! Creates, initializes and activates set of triangular selecting frustums for polyline selection
Standard_EXPORT void InitPolylineSelectingVolume (const TColgp_Array1OfPnt2d& thePoints);
//! Creates and activates axis selector for point selection
Standard_EXPORT void InitAxisSelectingVolume (const gp_Ax1& theAxis);
//! Sets as active the custom selecting volume
Standard_EXPORT void InitSelectingVolume (const Handle(SelectMgr_BaseIntersector)& theVolume);

View File

@@ -119,6 +119,22 @@ void SelectMgr_ViewerSelector3d::Pick (const TColgp_Array1OfPnt2d& thePolyline,
TraverseSensitives();
}
//=======================================================================
// Function: Pick
// Purpose :
//=======================================================================
void SelectMgr_ViewerSelector3d::Pick (const gp_Ax1& theAxis,
const Handle(V3d_View)& theView)
{
updateZLayers (theView);
mySelectingVolumeMgr.InitAxisSelectingVolume (theAxis);
mySelectingVolumeMgr.BuildSelectingVolume();
mySelectingVolumeMgr.SetViewClipping (theView->ClipPlanes(), Handle(Graphic3d_SequenceOfHClipPlane)(), NULL);
TraverseSensitives();
}
//=======================================================================
// Function: DisplaySensitive.
// Purpose : Display active primitives.

View File

@@ -54,6 +54,12 @@ public:
Standard_EXPORT void Pick (const TColgp_Array1OfPnt2d& thePolyline,
const Handle(V3d_View)& theView);
//! Picks the sensitive entity according to the input axis.
//! This is geometric intersection 3D objects by axis
//! (camera parameters are ignored and objects with transform persistance are skipped).
Standard_EXPORT void Pick (const gp_Ax1& theAxis,
const Handle(V3d_View)& theView);
//! Dump of detection results into image.
//! This method performs axis picking for each pixel in the image
//! and generates a color depending on picking results and selection image type.