Skip to content
Snippets Groups Projects
Commit 3a71dab2 authored by Léon Victor's avatar Léon Victor Committed by Nicolas Vuaille
Browse files

Add Frustum implicit function, representation and widget

parent 5fd6ec84
No related branches found
No related tags found
No related merge requests found
......@@ -84,6 +84,7 @@ set(classes
vtkExtractStructuredGridHelper
vtkFieldData
vtkFindCellStrategy
vtkFrustum
vtkGenericAdaptorCell
vtkGenericAttribute
vtkGenericAttributeCollection
......
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#include "vtkFrustum.h"
#include "vtkImplicitBoolean.h"
#include "vtkMath.h"
#include "vtkObjectFactory.h"
#include "vtkPlane.h"
#include "vtkVector.h"
#include <cmath>
VTK_ABI_NAMESPACE_BEGIN
vtkStandardNewMacro(vtkFrustum);
//------------------------------------------------------------------------------
vtkFrustum::vtkFrustum()
{
this->NearPlane->SetNormal(0, 1, 0);
this->NearPlane->SetOrigin(0, this->NearPlaneDistance, 0);
this->CalculateHorizontalPlanesNormal();
this->CalculateVerticalPlanesNormal();
this->BooleanOp->AddFunction(this->NearPlane);
this->BooleanOp->AddFunction(this->BottomPlane);
this->BooleanOp->AddFunction(this->TopPlane);
this->BooleanOp->AddFunction(this->RightPlane);
this->BooleanOp->AddFunction(this->LeftPlane);
this->BooleanOp->SetOperationTypeToUnion();
}
//------------------------------------------------------------------------------
vtkFrustum::~vtkFrustum() = default;
//------------------------------------------------------------------------------
double vtkFrustum::EvaluateFunction(double x[3])
{
return this->BooleanOp->EvaluateFunction(x);
}
//------------------------------------------------------------------------------
void vtkFrustum::EvaluateGradient(double x[3], double g[3])
{
return this->BooleanOp->EvaluateGradient(x, g);
}
//------------------------------------------------------------------------------
void vtkFrustum::SetHorizontalAngle(double angleInDegrees)
{
angleInDegrees = vtkMath::ClampValue(angleInDegrees, 1., 89.);
if (this->HorizontalAngle == angleInDegrees)
{
return;
}
this->HorizontalAngle = angleInDegrees;
this->CalculateHorizontalPlanesNormal();
this->Modified();
}
//------------------------------------------------------------------------------
void vtkFrustum::SetVerticalAngle(double angleInDegrees)
{
angleInDegrees = vtkMath::ClampValue(angleInDegrees, 1., 89.);
if (this->VerticalAngle == angleInDegrees)
{
return;
}
this->VerticalAngle = angleInDegrees;
this->CalculateVerticalPlanesNormal();
this->Modified();
}
//------------------------------------------------------------------------------
void vtkFrustum::SetNearPlaneDistance(double distance)
{
distance = std::max(distance, 0.0);
if (this->NearPlaneDistance == distance)
{
return;
}
this->NearPlaneDistance = distance;
this->NearPlane->SetOrigin(0, distance, 0);
this->Modified();
}
//------------------------------------------------------------------------------
void vtkFrustum::PrintSelf(ostream& os, vtkIndent indent)
{
this->Superclass::PrintSelf(os, indent);
os << indent << "Near Plane Distance: " << this->NearPlaneDistance << std::endl;
os << indent << "Horizontal Angle: " << this->HorizontalAngle << std::endl;
os << indent << "Vertical Angle: " << this->VerticalAngle << std::endl;
}
//------------------------------------------------------------------------------
void vtkFrustum::CalculateHorizontalPlanesNormal()
{
double angleRadians = vtkMath::RadiansFromDegrees(this->HorizontalAngle);
double cosAngle = std::cos(angleRadians);
double sinAngle = std::sin(angleRadians);
vtkVector3d leftPlaneNormal(cosAngle, sinAngle, 0);
vtkVector3d rightPlaneNormal(-cosAngle, sinAngle, 0);
this->RightPlane->SetNormal(rightPlaneNormal.GetData());
this->LeftPlane->SetNormal(leftPlaneNormal.GetData());
}
//------------------------------------------------------------------------------
void vtkFrustum::CalculateVerticalPlanesNormal()
{
double angleRadians = vtkMath::RadiansFromDegrees(this->VerticalAngle);
double cosAngle = std::cos(angleRadians);
double sinAngle = std::sin(angleRadians);
vtkVector3d topPlaneNormal(0, sinAngle, -cosAngle);
vtkVector3d bottomPlaneNormal(0, sinAngle, cosAngle);
this->TopPlane->SetNormal(topPlaneNormal.GetData());
this->BottomPlane->SetNormal(bottomPlaneNormal.GetData());
}
VTK_ABI_NAMESPACE_END
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
/**
* @class vtkFrustum
* @brief implicit function for a frustum
*
* vtkFrustum represents a 4-sided frustum, with a near plane but infinite on the far side. It is
* defined by the two angles between its forward axis and its horizontal and vertical planes, and
* the distance between its origin and near plane. vtkFrustum is a concrete implementation of
* vtkImplicitFunction. The frustum is oriented toward the Y Axis; its top face facing
* toward the Z Axis and its "right" face facing the X Axis.
*
* @warning
* The frustum is infinite in extent towards its far plane. To truncate the frustum in modeling
* operations use the vtkImplicitBoolean in combination with clipping planes.
*
*/
#ifndef vtkFrustum_h
#define vtkFrustum_h
#include "vtkCommonDataModelModule.h" // For export macro
#include "vtkImplicitFunction.h"
#include "vtkNew.h" // For vtkNew
VTK_ABI_NAMESPACE_BEGIN
class vtkImplicitBoolean;
class vtkPlane;
class VTKCOMMONDATAMODEL_EXPORT vtkFrustum : public vtkImplicitFunction
{
public:
static vtkFrustum* New();
vtkTypeMacro(vtkFrustum, vtkImplicitFunction);
void PrintSelf(ostream& os, vtkIndent indent) override;
using vtkImplicitFunction::EvaluateFunction;
double EvaluateFunction(double x[3]) override;
using vtkImplicitFunction::EvaluateGradient;
void EvaluateGradient(double x[3], double g[3]) override;
///@{
/**
* Get/Set the near plane distance of the frustum, i.e. the distance between its origin and near
* plane along the forward axis. Values below 0 will be clamped. Defaults to 0.5.
*/
vtkGetMacro(NearPlaneDistance, double)
void SetNearPlaneDistance(double distance);
///@}
///@{
/**
* Get/Set the horizontal angle of the frustum in degrees. It represents the angle between its
* forward axis and its right and left planes. Clamped between 1 and 89 degrees. Defaults to 30.
*/
vtkGetMacro(HorizontalAngle, double)
void SetHorizontalAngle(double angleInDegrees);
///@}
///@{
/**
* Get/Set the vertical angle of the frustum in degrees. It represents the angle between its
* forward axis and its top and bottom planes. Clamped between 1 and 89 degrees. Defaults to 30.
*/
vtkGetMacro(VerticalAngle, double)
void SetVerticalAngle(double angleInDegrees);
///@}
///@{
/**
* Get individual planes that make up the frustum.
* @note: Do not attempt to modify ! Use the vertical/horizontal angles and near plane distance to
* parameterize the frustum instead.
*/
VTK_FUTURE_CONST vtkPlane* GetTopPlane() { return this->TopPlane; }
VTK_FUTURE_CONST vtkPlane* GetBottomPlane() { return this->BottomPlane; }
VTK_FUTURE_CONST vtkPlane* GetRightPlane() { return this->RightPlane; }
VTK_FUTURE_CONST vtkPlane* GetLeftPlane() { return this->LeftPlane; }
VTK_FUTURE_CONST vtkPlane* GetNearPlane() { return this->NearPlane; }
///@}
protected:
vtkFrustum();
~vtkFrustum() override;
private:
vtkFrustum(const vtkFrustum&) = delete;
void operator=(const vtkFrustum&) = delete;
/// Compute and set the horizontal or vertical planes' normals according to the defined angle
void CalculateHorizontalPlanesNormal();
void CalculateVerticalPlanesNormal();
double NearPlaneDistance = 0.5;
double VerticalAngle = 30;
double HorizontalAngle = 30;
vtkNew<vtkPlane> NearPlane;
vtkNew<vtkPlane> BottomPlane;
vtkNew<vtkPlane> TopPlane;
vtkNew<vtkPlane> RightPlane;
vtkNew<vtkPlane> LeftPlane;
vtkNew<vtkImplicitBoolean> BooleanOp;
};
VTK_ABI_NAMESPACE_END
#endif
......@@ -86,6 +86,8 @@ set(classes
vtkImplicitConeWidget
vtkImplicitCylinderRepresentation
vtkImplicitCylinderWidget
vtkImplicitFrustumRepresentation
vtkImplicitFrustumWidget
vtkImplicitImageRepresentation
vtkImplicitPlaneRepresentation
vtkImplicitPlaneWidget
......
......@@ -38,6 +38,7 @@ vtk_add_test_cxx(vtkInteractionWidgetsCxxTests tests
TestImplicitConeWidget.cxx
TestImplicitCylinderWidget.cxx
TestImplicitCylinderWidget2.cxx
TestImplicitFrustumWidget.cxx
TestImplicitPlaneWidget.cxx
TestImplicitPlaneWidget2.cxx
TestImplicitPlaneWidget2LockNormalToCamera.cxx
......
This diff is collapsed.
This diff is collapsed.
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#ifndef vtkImplicitFrustumRepresentation_h
#define vtkImplicitFrustumRepresentation_h
#include "vtkEllipseArcSource.h"
#include "vtkInteractionWidgetsModule.h" // For export macro
#include "vtkSetGet.h"
#include "vtkVector.h" // For vtkVector3d
#include "vtkWidgetRepresentation.h"
#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
VTK_ABI_NAMESPACE_BEGIN
class vtkActor;
class vtkBox;
class vtkConeSource;
class vtkCellPicker;
class vtkEllipseArcSource;
class vtkFrustum;
class vtkLineSource;
class vtkProperty;
class vtkPolyData;
class vtkPolyDataMapper;
class vtkSphereSource;
class vtkTransform;
class vtkTubeFilter;
/**
* @class vtkImplicitFrustumRepresentation
* @brief The representation for a vtkImplicitFrustumWidget
*
* This class is a concrete representation for the vtkImplicitFrustumWidget. It represents an
* infinite frustum defined by its origin, its orientation, the two angles between its forward axis
* and its horizontal and vertical planes, and the distance between its origin and near plane. This
* frustum representation can be manipulated by using the vtkImplicitFrustumWidget.
*
* @sa vtkImplicitFrustumWidget vtkFrustum
*/
class VTKINTERACTIONWIDGETS_EXPORT VTK_MARSHALAUTO vtkImplicitFrustumRepresentation
: public vtkWidgetRepresentation
{
public:
// Manage the state of the widget
enum InteractionStateType
{
Outside = 0,
Moving, // Generic state set by the widget
MovingOrigin,
AdjustingHorizontalAngle,
AdjustingVerticalAngle,
AdjustingNearPlaneDistance,
AdjustingYaw,
AdjustingPitch,
AdjustingRoll,
TranslatingOriginOnAxis
};
static vtkImplicitFrustumRepresentation* New();
vtkTypeMacro(vtkImplicitFrustumRepresentation, vtkWidgetRepresentation);
void PrintSelf(ostream& os, vtkIndent indent) override;
///@{
/**
* Get the origin of the frustum representation. The origin is located along the frustum axis.
* Default is (0, 0, 0)
*/
void SetOrigin(double x, double y, double z);
void SetOrigin(double x[3]);
void SetOrigin(const vtkVector3d& xyz);
double* GetOrigin() VTK_SIZEHINT(3);
void GetOrigin(double xyz[3]) const;
///@}
///@{
/**
* Get/Set the orientation of the frustum. Defaults to (0,0,0)
*/
void SetOrientation(double x, double y, double z);
void SetOrientation(const double xyz[3]);
void SetOrientation(const vtkVector3d& xyz);
double* GetOrientation() VTK_SIZEHINT(3);
void GetOrientation(double& x, double& y, double& z);
void GetOrientation(double xyz[3]);
///@}
///@{
/**
* Get/Set the horizontal angle of the frustum in degrees. It represents the angle between its
* forward axis and its right and left planes. Defaults to 30.
*/
double GetHorizontalAngle() const;
void SetHorizontalAngle(double angle);
///@}
///@{
/**
* Get/Set the vertical angle of the frustum in degrees. It represents the angle between its
* forward axis and its top and bottom planes. Defaults to 30.
*/
double GetVerticalAngle() const;
void SetVerticalAngle(double angle);
///@}
///@{
/**
* Get/Set the near plane distance of the frustum, i.e. the distance between its origin and near
* plane along the forward axis. Defaults to 0.5.
*/
double GetNearPlaneDistance() const;
void SetNearPlaneDistance(double angle);
///@}
///@{
/**
* Force the frustum widget to be aligned with one of the x-y-z axes.
* If one axis is set on, the other two will be set off.
* Remember that when the state changes, a ModifiedEvent is invoked.
* This can be used to snap the frustum to the axes if it is originally not aligned.
* Default to false.
*/
void SetAlongXAxis(bool);
vtkGetMacro(AlongXAxis, bool);
vtkBooleanMacro(AlongXAxis, bool);
void SetAlongYAxis(bool);
vtkGetMacro(AlongYAxis, bool);
vtkBooleanMacro(AlongYAxis, bool);
void SetAlongZAxis(bool);
vtkGetMacro(AlongZAxis, bool);
vtkBooleanMacro(AlongZAxis, bool);
///@}
///@{
/**
* Enable/disable the drawing of the frustum. In some cases the frustum interferes with the object
* that it is operating on (e.g., the frustum interferes with the cut surface it produces
* resulting in z-buffer artifacts.) By default it is on.
*/
void SetDrawFrustum(bool draw);
vtkGetMacro(DrawFrustum, bool);
vtkBooleanMacro(DrawFrustum, bool);
///@}
///@{
/**
* Set/Get the bounds of the widget representation. PlaceWidget can also be used to set the bounds
* of the widget but it may also have other effects on the internal state of the representation.
* Use this function when only the widget bounds need to be modified.
*/
vtkSetVector6Macro(WidgetBounds, double);
double* GetWidgetBounds();
///@}
/**
* Grab the polydata that defines the frustum. The polydata contains polygons that are clipped by
* the bounding box.
*/
void GetPolyData(vtkPolyData* pd);
/**
* Satisfies the superclass API. This will change the state of the widget to match changes that
* have been made to the underlying PolyDataSource.
*/
void UpdatePlacement();
/**
* Get properties of the frustum actor.
*/
vtkGetObjectMacro(FrustumProperty, vtkProperty);
///@{
/**
* Get the properties of the edge handles actors. i.e. the properties of the near plane and angles
* handles when selected or not.
*/
vtkGetObjectMacro(EdgeHandleProperty, vtkProperty);
vtkGetObjectMacro(SelectedEdgeHandleProperty, vtkProperty);
///@}
///@{
/**
* Set the color of all the widgets handles (origin, orientations, near plane and angles) and
* their color during interaction. Foreground color applies to the frustum itself
*/
void SetInteractionColor(double, double, double);
void SetInteractionColor(double c[3]) { this->SetInteractionColor(c[0], c[1], c[2]); }
void SetHandleColor(double, double, double);
void SetHandleColor(double c[3]) { this->SetHandleColor(c[0], c[1], c[2]); }
void SetForegroundColor(double, double, double);
void SetForegroundColor(double c[3]) { this->SetForegroundColor(c[0], c[1], c[2]); }
///@}
///@{
/**
* Methods to interface with the vtkImplicitFrustumWidget.
*/
int ComputeInteractionState(int X, int Y, int modify = 0) override;
void PlaceWidget(double bounds[6]) override;
void BuildRepresentation() override;
void StartWidgetInteraction(double eventPos[2]) override;
void WidgetInteraction(double newEventPos[2]) override;
void EndWidgetInteraction(double newEventPos[2]) override;
///@}
///@{
/**
* Methods supporting the rendering process.
*/
double* GetBounds() override;
void GetActors(vtkPropCollection* pc) override;
void ReleaseGraphicsResources(vtkWindow*) override;
int RenderOpaqueGeometry(vtkViewport*) override;
int RenderTranslucentPolygonalGeometry(vtkViewport*) override;
vtkTypeBool HasTranslucentPolygonalGeometry() override;
///@}
/**
* The interaction state may be set from a widget (e.g., vtkImplicitFrustumWidget) or other
* object. This controls how the interaction with the widget proceeds. Normally this method is
* used as part of a handshaking process with the widget: First ComputeInteractionState() is
* invoked and returns a state based on geometric considerations (i.e., cursor near a widget
* feature), then based on events, the widget may modify this further.
*/
vtkSetClampMacro(InteractionState, InteractionStateType, InteractionStateType::Outside,
InteractionStateType::TranslatingOriginOnAxis);
///@{
/**
* Sets the visual appearance of the representation based on the state it is in. This state is
* usually the same as InteractionState.
*/
virtual void SetRepresentationState(InteractionStateType);
vtkGetMacro(RepresentationState, InteractionStateType);
///@}
/*
* Register internal Pickers within PickingManager
*/
void RegisterPickers() override;
///@{
/**
* Gets/Sets the constraint axis for translations.
* Defaults to Axis::NONE
**/
vtkGetMacro(TranslationAxis, int);
vtkSetClampMacro(TranslationAxis, int, Axis::NONE, Axis::ZAxis);
///@}
///@{
/**
* Toggles constraint translation axis on/off.
*/
void SetXTranslationAxisOn() { this->TranslationAxis = Axis::XAxis; }
void SetYTranslationAxisOn() { this->TranslationAxis = Axis::YAxis; }
void SetZTranslationAxisOn() { this->TranslationAxis = Axis::ZAxis; }
void SetTranslationAxisOff() { this->TranslationAxis = Axis::NONE; }
///@}
/**
* Returns true if translation is constrained to an axis
**/
bool IsTranslationConstrained() { return this->TranslationAxis != Axis::NONE; }
/**
* Get the concrete represented frustum
**/
void GetFrustum(vtkFrustum* frustum) const;
protected:
vtkImplicitFrustumRepresentation();
~vtkImplicitFrustumRepresentation() override;
private:
enum class FrustumFace
{
None = -1,
Right = 0,
Left,
Top,
Bottom,
Near,
};
struct SphereHandle
{
vtkNew<vtkSphereSource> Source;
vtkNew<vtkPolyDataMapper> Mapper;
vtkNew<vtkActor> Actor;
SphereHandle();
};
struct EdgeHandle
{
vtkNew<vtkPolyData> PolyData;
vtkNew<vtkTubeFilter> Tuber;
vtkNew<vtkPolyDataMapper> Mapper;
vtkNew<vtkActor> Actor;
EdgeHandle();
};
struct EllipseHandle
{
vtkNew<vtkEllipseArcSource> Source;
vtkNew<vtkTubeFilter> Tuber;
vtkNew<vtkPolyDataMapper> Mapper;
vtkNew<vtkActor> Actor;
EllipseHandle();
};
vtkImplicitFrustumRepresentation(const vtkImplicitFrustumRepresentation&) = delete;
void operator=(const vtkImplicitFrustumRepresentation&) = delete;
// Helpers to get the frustum basis
vtkVector3d GetForwardAxis();
vtkVector3d GetUpAxis();
vtkVector3d GetRightAxis();
void HighlightOriginHandle(bool highlight);
void HighlightFarPlaneVerticalHandle(bool highlight);
void HighlightFarPlaneHorizontalHandle(bool highlight);
void HighlightNearPlaneHandle(bool highlight);
void HighlightRollHandle(bool highlight);
void HighlightYawHandle(bool highlight);
void HighlightPitchHandle(bool highlight);
// Methods to manipulate the frustum
void TranslateOrigin(const vtkVector3d& p1, const vtkVector3d& p2);
void TranslateOriginOnAxis(const vtkVector3d& p1, const vtkVector3d& p2);
void AdjustHorizontalAngle(const vtkVector3d& p1, const vtkVector3d& p2);
void AdjustVerticalAngle(const vtkVector3d& p1, const vtkVector3d& p2);
void AdjustNearPlaneDistance(
const vtkVector2d& eventPosition, const vtkVector3d& p1, const vtkVector3d& p2);
void Rotate(
const vtkVector3d& prevPickPoint, const vtkVector3d& pickPoint, const vtkVector3d& axis);
// Set the frustum transform according to the representation's orientation and position
void UpdateFrustumTransform();
// Re-compute the widget handles' sizes
void SizeHandles();
// Generate the frustum polydata, cropped by the bounding box
void BuildFrustum();
// The actual frustum we're manipulating
vtkNew<vtkFrustum> Frustum;
InteractionStateType RepresentationState = InteractionStateType::Outside;
int TranslationAxis = Axis::NONE;
// Keep track of event positions
vtkVector3d LastEventPosition = { 0., 0., 0. };
bool AlongXAxis = false;
bool AlongYAxis = false;
bool AlongZAxis = false;
vtkVector<double, 6> WidgetBounds;
double Length = 1;
vtkVector3d Origin = { 0, 0, 0 };
vtkNew<vtkTransform> OrientationTransform;
vtkNew<vtkPolyData> FrustumPD;
vtkNew<vtkPolyDataMapper> FrustumMapper;
vtkNew<vtkActor> FrustumActor;
bool DrawFrustum = true;
std::array<EdgeHandle, 4> FarPlaneHandles;
EdgeHandle NearPlaneEdgesHandle;
SphereHandle NearPlaneCenterHandle;
EllipseHandle RollHandle;
EllipseHandle YawHandle;
EllipseHandle PitchHandle;
SphereHandle OriginHandle;
FrustumFace ActiveEdgeHandle = FrustumFace::None;
vtkNew<vtkCellPicker> Picker;
vtkNew<vtkCellPicker> FrustumPicker;
// Properties used to control the appearance of selected objects and the manipulator in general.
vtkNew<vtkProperty> FrustumProperty;
vtkNew<vtkProperty> EdgeHandleProperty;
vtkNew<vtkProperty> SelectedEdgeHandleProperty;
vtkNew<vtkProperty> OriginHandleProperty;
vtkNew<vtkProperty> SelectedOriginHandleProperty;
vtkNew<vtkBox> BoundingBox;
};
VTK_ABI_NAMESPACE_END
#endif
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
#include "vtkImplicitFrustumWidget.h"
#include "vtkCallbackCommand.h"
#include "vtkCommand.h"
#include "vtkEvent.h"
#include "vtkImplicitFrustumRepresentation.h"
#include "vtkRenderWindow.h"
#include "vtkRenderWindowInteractor.h"
#include "vtkWidgetCallbackMapper.h"
#include "vtkWidgetEvent.h"
#include "vtkWidgetEventTranslator.h"
VTK_ABI_NAMESPACE_BEGIN
vtkStandardNewMacro(vtkImplicitFrustumWidget);
//------------------------------------------------------------------------------
vtkImplicitFrustumWidget::vtkImplicitFrustumWidget()
{
// Define widget events
// Left mouse button: widget handles selection
this->CallbackMapper->SetCallbackMethod(vtkCommand::LeftButtonPressEvent, vtkWidgetEvent::Select,
this, vtkImplicitFrustumWidget::SelectAction);
this->CallbackMapper->SetCallbackMethod(vtkCommand::LeftButtonReleaseEvent,
vtkWidgetEvent::EndSelect, this, vtkImplicitFrustumWidget::EndSelectAction);
// Middle mouse button: Translate
this->CallbackMapper->SetCallbackMethod(vtkCommand::MiddleButtonPressEvent,
vtkWidgetEvent::Translate, this, vtkImplicitFrustumWidget::TranslateAction);
this->CallbackMapper->SetCallbackMethod(vtkCommand::MiddleButtonReleaseEvent,
vtkWidgetEvent::EndTranslate, this, vtkImplicitFrustumWidget::EndSelectAction);
// Mouse motion
this->CallbackMapper->SetCallbackMethod(
vtkCommand::MouseMoveEvent, vtkWidgetEvent::Move, this, vtkImplicitFrustumWidget::MoveAction);
// X Key: Lock/Unlock X Axis
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 'x', 1,
"x", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 24, 1,
"x", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 'X', 1,
"X", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 'x',
1, "x", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 24, 1,
"x", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 'X',
1, "X", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
// Y Key: Lock/Unlock Y Axis
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 'y', 1,
"y", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 25, 1,
"y", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 'Y', 1,
"Y", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 'y',
1, "y", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 25, 1,
"y", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 'Y',
1, "Y", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
// Z Key: Lock/Unlock Z Axis
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 'z', 1,
"z", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 26, 1,
"z", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyPressEvent, vtkEvent::AnyModifier, 'Z', 1,
"Z", vtkWidgetEvent::ModifyEvent, this, vtkImplicitFrustumWidget::TranslationAxisLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 'z',
1, "z", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 26, 1,
"z", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
this->CallbackMapper->SetCallbackMethod(vtkCommand::KeyReleaseEvent, vtkEvent::AnyModifier, 'Z',
1, "Z", vtkWidgetEvent::Reset, this, vtkImplicitFrustumWidget::TranslationAxisUnLock);
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::SelectAction(vtkAbstractWidget* w)
{
vtkImplicitFrustumWidget* self = reinterpret_cast<vtkImplicitFrustumWidget*>(w);
vtkImplicitFrustumRepresentation* repr = self->GetFrustumRepresentation();
// Get the event position
int X = self->Interactor->GetEventPosition()[0];
int Y = self->Interactor->GetEventPosition()[1];
// We want to update the angle, axis and origin as appropriate
repr->SetInteractionState(vtkImplicitFrustumRepresentation::InteractionStateType::Moving);
int interactionState = repr->ComputeInteractionState(X, Y);
self->UpdateCursorShape(interactionState);
if (repr->GetInteractionState() ==
vtkImplicitFrustumRepresentation::InteractionStateType::Outside)
{
return;
}
if (self->Interactor->GetControlKey() &&
interactionState == vtkImplicitFrustumRepresentation::InteractionStateType::MovingOrigin)
{
repr->SetInteractionState(
vtkImplicitFrustumRepresentation::InteractionStateType::TranslatingOriginOnAxis);
}
// We are definitely selected
self->GrabFocus(self->EventCallbackCommand);
double eventPos[2] = { static_cast<double>(X), static_cast<double>(Y) };
self->State = WidgetStateType::Active;
repr->StartWidgetInteraction(eventPos);
self->EventCallbackCommand->SetAbortFlag(true);
self->StartInteraction();
self->InvokeEvent(vtkCommand::StartInteractionEvent, nullptr);
self->Render();
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::TranslateAction(vtkAbstractWidget* w)
{
vtkImplicitFrustumWidget* self = reinterpret_cast<vtkImplicitFrustumWidget*>(w);
vtkImplicitFrustumRepresentation* repr = self->GetFrustumRepresentation();
// Get the event position
int X = self->Interactor->GetEventPosition()[0];
int Y = self->Interactor->GetEventPosition()[1];
// We want to compute an orthogonal vector to the plane that has been selected
repr->SetInteractionState(vtkImplicitFrustumRepresentation::InteractionStateType::MovingOrigin);
int interactionState = repr->ComputeInteractionState(X, Y);
self->UpdateCursorShape(interactionState);
if (repr->GetInteractionState() ==
vtkImplicitFrustumRepresentation::InteractionStateType::Outside)
{
return;
}
// We are definitely selected
self->GrabFocus(self->EventCallbackCommand);
self->State = WidgetStateType::Active;
double eventPos[2] = { static_cast<double>(X), static_cast<double>(Y) };
repr->StartWidgetInteraction(eventPos);
self->EventCallbackCommand->SetAbortFlag(true);
self->StartInteraction();
self->InvokeEvent(vtkCommand::StartInteractionEvent, nullptr);
self->Render();
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::MoveAction(vtkAbstractWidget* w)
{
vtkImplicitFrustumWidget* self = reinterpret_cast<vtkImplicitFrustumWidget*>(w);
vtkImplicitFrustumRepresentation* representation = self->GetFrustumRepresentation();
// Change the cursor shape when the mouse is hovering
// the widget. Unfortunately, this results in a few extra picks
// due to the cell picker. However given that its picking simple geometry
// like the handles/arrows, this should be very quick
int X = self->Interactor->GetEventPosition()[0];
int Y = self->Interactor->GetEventPosition()[1];
bool changed = false;
if (self->ManagesCursor && self->State != WidgetStateType::Active)
{
auto oldState = static_cast<vtkImplicitFrustumRepresentation::InteractionStateType>(
representation->GetInteractionState());
representation->SetInteractionState(
vtkImplicitFrustumRepresentation::InteractionStateType::Moving);
auto newState = static_cast<vtkImplicitFrustumRepresentation::InteractionStateType>(
representation->ComputeInteractionState(X, Y));
changed = self->UpdateCursorShape(newState);
representation->SetInteractionState(oldState);
changed |= newState != oldState;
}
// See whether we're active
if (self->State == WidgetStateType::Idle)
{
if (changed && self->ManagesCursor)
{
self->Render();
}
return;
}
// Okay, adjust the representation
double eventPos[2] = { static_cast<double>(X), static_cast<double>(Y) };
representation->WidgetInteraction(eventPos);
// moving something
self->EventCallbackCommand->SetAbortFlag(true);
self->InvokeEvent(vtkCommand::InteractionEvent, nullptr);
self->Render();
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::EndSelectAction(vtkAbstractWidget* w)
{
vtkImplicitFrustumWidget* self = reinterpret_cast<vtkImplicitFrustumWidget*>(w);
vtkImplicitFrustumRepresentation* repr = self->GetFrustumRepresentation();
if (self->State != WidgetStateType::Active ||
repr->GetInteractionState() == vtkImplicitFrustumRepresentation::InteractionStateType::Outside)
{
return;
}
// Return state to not selected
double eventPos[2];
repr->EndWidgetInteraction(eventPos);
self->State = WidgetStateType::Idle;
self->ReleaseFocus();
// Update cursor if managed
self->UpdateCursorShape(repr->GetRepresentationState());
self->EventCallbackCommand->SetAbortFlag(true);
self->EndInteraction();
self->InvokeEvent(vtkCommand::EndInteractionEvent, nullptr);
self->Render();
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::CreateDefaultRepresentation()
{
if (!this->WidgetRep)
{
this->WidgetRep = vtkImplicitFrustumRepresentation::New();
}
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::SetRepresentation(vtkImplicitFrustumRepresentation* rep)
{
this->Superclass::SetWidgetRepresentation(reinterpret_cast<vtkWidgetRepresentation*>(rep));
}
//------------------------------------------------------------------------------
bool vtkImplicitFrustumWidget::UpdateCursorShape(int state)
{
// So as to change the cursor shape when the mouse is poised over
// the widget.
if (this->ManagesCursor)
{
if (state == vtkImplicitFrustumRepresentation::InteractionStateType::Outside)
{
return this->RequestCursorShape(VTK_CURSOR_DEFAULT);
}
else
{
return this->RequestCursorShape(VTK_CURSOR_HAND);
}
}
return false;
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::TranslationAxisLock(vtkAbstractWidget* widget)
{
vtkImplicitFrustumWidget* self = reinterpret_cast<vtkImplicitFrustumWidget*>(widget);
vtkImplicitFrustumRepresentation* repr = self->GetFrustumRepresentation();
char* cKeySym = self->Interactor->GetKeySym();
std::string keySym = cKeySym != nullptr ? cKeySym : "";
std::transform(keySym.begin(), keySym.end(), keySym.begin(), ::toupper);
if (keySym == "X")
{
repr->SetXTranslationAxisOn();
}
else if (keySym == "Y")
{
repr->SetYTranslationAxisOn();
}
else if (keySym == "Z")
{
repr->SetZTranslationAxisOn();
}
}
//------------------------------------------------------------------------------
void vtkImplicitFrustumWidget::TranslationAxisUnLock(vtkAbstractWidget* widget)
{
vtkImplicitFrustumWidget* self = reinterpret_cast<vtkImplicitFrustumWidget*>(widget);
self->GetFrustumRepresentation()->SetTranslationAxisOff();
}
VTK_ABI_NAMESPACE_END
// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
// SPDX-License-Identifier: BSD-3-Clause
/**
* @class vtkImplicitFrustumWidget
* @brief 3D widget for manipulating an infinite frustum
*
* This 3D widget defines an infinite frustum that can be
* interactively placed in a scene.
*
* To use this widget, you generally pair it with a vtkImplicitFrustumRepresentation
* (or a subclass). Various options are available for controlling how the
* representation appears, and how the widget functions.
*
* @par Event Bindings:
* By default, the widget responds to the following VTK events (i.e., it
* watches the vtkRenderWindowInteractor for these events):
*
* - LeftButtonPressEvent - select a widget handle
* - LeftButtonReleaseEvent - release the currently held widget handle
* - MouseMoveEvent - dependant on the current manipulation mode:
* - Origin handle: Translate the frustum origin (constrained to the the x, y or z axis one of the
* corresponding key is held)
* - Near plane edges handle: Adjust the near plane distance
* - Far plane edges handle: Adjust the horizontal/vertical frustum angles
*
* In all the cases, independent of what is picked, the widget responds to the
* following VTK events:
*
* - MiddleButtonPressEvent - grab the frustum
* - MiddleButtonReleaseEvent - release the frustum
* - MouseMoveEvent - move the widget (if middle button is pressed)
*
* @sa
* vtk3DWidget vtkImplicitFrustumRepresentation vtkFrustum
*/
#ifndef vtkImplicitFrustumWidget_h
#define vtkImplicitFrustumWidget_h
#include "vtkAbstractWidget.h"
#include "vtkInteractionWidgetsModule.h" // For export macro
#include "vtkWrappingHints.h" // For VTK_MARSHALAUTO
VTK_ABI_NAMESPACE_BEGIN
class vtkImplicitFrustumRepresentation;
class VTKINTERACTIONWIDGETS_EXPORT VTK_MARSHALAUTO vtkImplicitFrustumWidget
: public vtkAbstractWidget
{
public:
static vtkImplicitFrustumWidget* New();
vtkTypeMacro(vtkImplicitFrustumWidget, vtkAbstractWidget);
/**
* Specify an instance of vtkWidgetRepresentation used to represent this
* widget in the scene. Note that the representation is a subclass of vtkProp
* so it can be added to the renderer independent of the widget.
*/
void SetRepresentation(vtkImplicitFrustumRepresentation* rep);
/**
* Return the representation as a vtkImplicitFrustumRepresentation.
*/
vtkImplicitFrustumRepresentation* GetFrustumRepresentation()
{
return reinterpret_cast<vtkImplicitFrustumRepresentation*>(this->WidgetRep);
}
/**
* Create the default widget representation if one is not set.
*/
void CreateDefaultRepresentation() override;
private:
enum class WidgetStateType
{
Idle,
Active
};
vtkImplicitFrustumWidget();
~vtkImplicitFrustumWidget() override = default;
vtkImplicitFrustumWidget(const vtkImplicitFrustumWidget&) = delete;
void operator=(const vtkImplicitFrustumWidget&) = delete;
///@{
/**
* Callbacks for widget events
*/
static void SelectAction(vtkAbstractWidget* widget);
static void EndSelectAction(vtkAbstractWidget* widget);
static void TranslateAction(vtkAbstractWidget* widget);
static void MoveAction(vtkAbstractWidget* widget);
static void TranslationAxisLock(vtkAbstractWidget* widget);
static void TranslationAxisUnLock(vtkAbstractWidget* widget);
/// @}
/**
* Update the cursor shape based on the interaction state. Returns true
* if the cursor shape requested is different from the existing one.
*/
bool UpdateCursorShape(int interactionState);
WidgetStateType State = WidgetStateType::Idle;
};
VTK_ABI_NAMESPACE_END
#endif
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment