Mantid
|
A generalized description of a N-dimensional hyperplane. More...
#include <MDPlane.h>
Public Member Functions | |
bool | doesLineIntersect (const coord_t *pointA, const coord_t *pointB) const |
Given two points defining the start and end point of a line, is there an intersection between the hyperplane and the line defined by the points? More... | |
coord_t | getInequality () const |
const coord_t * | getNormal () const |
size_t | getNumDims () const |
bool | isPointBounded (const coord_t *coords) const |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... >= b)? More... | |
bool | isPointBounded (const Mantid::Kernel::VMD &coords) const |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... >= b)? More... | |
bool | isPointBounded (const std::vector< coord_t > &coords) const |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... >= b)? More... | |
bool | isPointInside (const coord_t *coords) const |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... > b)? False for points that lie on the hyperplane, this is used to detect when two volumes (for example an MDBox and a mask) touch but do not share a finite volume. More... | |
bool | isPointInside (const std::vector< coord_t > &coords) const |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... > b)? False for points that lie on the hyperplane, this is used to detect when two volumes (for example an MDBox and a mask) touch but do not share a finite volume. More... | |
MDPlane (const Mantid::Kernel::VMD &normal, const Mantid::Kernel::VMD &point) | |
Constructor with normal and point. More... | |
MDPlane (const MDPlane &other) | |
Copy constructor. More... | |
MDPlane (const size_t nd, const double *normal, const double *point) | |
Constructor with normal and point. More... | |
MDPlane (const size_t nd, const float *normal, const float *point) | |
Constructor with normal and point. More... | |
MDPlane (const std::vector< coord_t > &normal, const std::vector< coord_t > &point) | |
Constructor with normal and point. More... | |
MDPlane (const std::vector< Mantid::Kernel::VMD > &vectors, const Mantid::Kernel::VMD &origin, const Mantid::Kernel::VMD &insidePoint) | |
Constructor with N-1 vectors on the plane's surface. More... | |
MDPlane & | operator= (const MDPlane &other) |
Assignment operator. More... | |
~MDPlane () | |
Destructor. More... | |
Protected Attributes | |
coord_t | m_inequality |
Right-hand side of the linear equation. aka b in : a1*x1 + a2*x2 + ... < b. More... | |
size_t | m_nd |
Number of dimensions in the MDEventWorkspace. More... | |
coord_t * | m_normal |
Coefficients to multiply each coordinate, sized m_nd. More... | |
Private Member Functions | |
template<typename IterableType > | |
void | construct (IterableType normal, IterableType point) |
Templated method for constructing the MDPlane. More... | |
A generalized description of a N-dimensional hyperplane.
To be used in MDImplicitFunction.
This has be to fully general, with: nd : number of dimensions of space
The general equation for a hyperplane is: a1*x1 + a2*x2 + ... > b
where x1, x2, .. are the n-th coordinate of the point. a1, a2, .. are coefficients (can be 0).
Any plane can be defined with:
Points that are in the direction of the normal of the plane are considered to be bounded by it.
Mantid::Geometry::MDPlane::MDPlane | ( | const Mantid::Kernel::VMD & | normal, |
const Mantid::Kernel::VMD & | point | ||
) |
Constructor with normal and point.
normal | :: normal to the plane. Points that are in the direction of the normal of the plane are considered to be bounded by it. |
point | :: any point that is on the plane |
Definition at line 42 of file MDPlane.cpp.
References construct(), Mantid::Kernel::VMDBase< TYPE >::getNumDims(), and m_nd.
Mantid::Geometry::MDPlane::MDPlane | ( | const std::vector< coord_t > & | normal, |
const std::vector< coord_t > & | point | ||
) |
Constructor with normal and point.
normal | :: normal to the plane. Points that are in the direction of the normal of the plane are considered to be bounded by it. |
point | :: any point that is on the plane |
Definition at line 25 of file MDPlane.cpp.
References construct(), and m_nd.
Mantid::Geometry::MDPlane::MDPlane | ( | const size_t | nd, |
const float * | normal, | ||
const float * | point | ||
) |
Constructor with normal and point.
nd | :: number of dimensions |
normal | :: normal to the plane. Points that are in the direction of the normal of the plane are considered to be bounded by it. |
point | :: any point that is on the plane |
Definition at line 60 of file MDPlane.cpp.
References construct(), and m_nd.
Mantid::Geometry::MDPlane::MDPlane | ( | const size_t | nd, |
const double * | normal, | ||
const double * | point | ||
) |
Constructor with normal and point.
nd | :: number of dimensions |
normal | :: normal to the plane. Points that are in the direction of the normal of the plane are considered to be bounded by it. |
point | :: any point that is on the plane |
Definition at line 74 of file MDPlane.cpp.
References construct(), and m_nd.
Mantid::Geometry::MDPlane::MDPlane | ( | const std::vector< Mantid::Kernel::VMD > & | vectors, |
const Mantid::Kernel::VMD & | origin, | ||
const Mantid::Kernel::VMD & | insidePoint | ||
) |
Constructor with N-1 vectors on the plane's surface.
vectors | :: vector of N-1 vectors with N dimensions. |
origin | :: any point on the plane |
insidePoint | :: coordinate of a point that is known to be inside the plane described |
if | the vectors are collinear |
Definition at line 89 of file MDPlane.cpp.
References construct(), Mantid::Geometry::d, Mantid::Kernel::VMDBase< VMD_t >::getNormalVector(), Mantid::Kernel::VMDBase< TYPE >::getNumDims(), isPointBounded(), m_nd, and m_normal.
Mantid::Geometry::MDPlane::MDPlane | ( | const MDPlane & | other | ) |
Copy constructor.
other | :: MDPlane to copy |
Definition at line 228 of file MDPlane.cpp.
References Mantid::Geometry::d, m_nd, and m_normal.
Mantid::Geometry::MDPlane::~MDPlane | ( | ) |
|
inlineprivate |
Templated method for constructing the MDPlane.
m_nd must have been set before this call!
normal | :: vector giving the normal to the surface |
point | :: vector giving one point on the surface |
Definition at line 170 of file MDPlane.h.
References Mantid::Geometry::d.
Referenced by MDPlane().
|
inline |
Given two points defining the start and end point of a line, is there an intersection between the hyperplane and the line defined by the points?
pointA | :: first point/vertex; nd-sized array of coordinates |
pointB | :: last point/vertex; nd-sized array of coordinates |
|
inline |
|
inline |
Definition at line 57 of file MDPlane.h.
Referenced by Mantid::Geometry::MDPlaneImplicitFunction::toXMLString().
|
inline |
Definition at line 54 of file MDPlane.h.
Referenced by Mantid::Geometry::MDImplicitFunction::addPlane().
|
inline |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... >= b)?
coords | :: nd-sized array of coordinates |
Definition at line 71 of file MDPlane.h.
References Mantid::Geometry::d.
Referenced by MDPlane().
|
inline |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... >= b)?
coords | :: VMD vector giving the point of the right size |
Definition at line 86 of file MDPlane.h.
References Mantid::Geometry::d.
|
inline |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... >= b)?
coords | :: nd-sized vector of coordinates. No size check is made! |
Definition at line 101 of file MDPlane.h.
References Mantid::Geometry::d.
|
inline |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... > b)? False for points that lie on the hyperplane, this is used to detect when two volumes (for example an MDBox and a mask) touch but do not share a finite volume.
coords | :: nd-sized array of coordinates |
Definition at line 119 of file MDPlane.h.
References Mantid::Geometry::d.
|
inline |
Is a point in MDimensions bounded by this hyperplane, that is, is (a1*x1 + a2*x2 + ... > b)? False for points that lie on the hyperplane, this is used to detect when two volumes (for example an MDBox and a mask) touch but do not share a finite volume.
coords | :: nd-sized vector of coordinates. No size check is made! |
Definition at line 137 of file MDPlane.h.
References Mantid::Geometry::d.
Assignment operator.
other | :: MDPlane to copy |
Definition at line 239 of file MDPlane.cpp.
References Mantid::Geometry::d, m_inequality, m_nd, and m_normal.
|
protected |
Right-hand side of the linear equation. aka b in : a1*x1 + a2*x2 + ... < b.
Definition at line 189 of file MDPlane.h.
Referenced by operator=().
|
protected |
Number of dimensions in the MDEventWorkspace.
Definition at line 181 of file MDPlane.h.
Referenced by MDPlane(), and operator=().
|
protected |
Coefficients to multiply each coordinate, sized m_nd.
This is the normal to the plane.
Definition at line 186 of file MDPlane.h.
Referenced by MDPlane(), operator=(), and ~MDPlane().