Model#
Overview#
The model is the heart of the Polyhedral Gravity Model since it contains the the two major classes:
The class polyhedralGravity::Polyhedron
stores the mesh and density
information and governs the compliance with Tsoulis et al.’s gravity model’s preconditions
It ensures that all plane unit normals of a constructed Polyhedron are consistently
pointing polyhedralGravity::NormalOrientation::OUTWARDS
or
polyhedralGravity::NormalOrientation::INWARDS
.
Depending on the value of polyhedralGravity::PolyhedronIntegrity
, these
constraints are either enforced by modifying input mesh data or by throwing
an std::invalid_argument
exception.
For this purpose, it uses the Möller–Trumbore intersection algorithm.
The class polyhedralGravity::GravityEvaluable
is provides as a way to
perform the evaluation of the polyhedral gravity model repeatedly
without the need to re-initialize the polyhedron and the gravity model as
caching is performed.
It takes a polyhedralGravity::Polyhedron
and provides an
polyhedralGravity::GravityEvaluable::operator()()
to evaluate the
model at computation point(s) \(P\) optionally using parallelization.
A polyhedralGravity::GravityModel::evaluate()
summarizes the
functionality polyhedralGravity::GravityEvaluable
, but does not
provide any caching throughout multiple calls.
Polyhedron#
-
class Polyhedron#
Data structure containing the model data of one polyhedron. This includes nodes, edges (faces) and elements. The index always starts with zero!
Public Functions
-
Polyhedron(const std::vector<Array3> &vertices, const std::vector<IndexArray3> &faces, double density, const NormalOrientation &orientation = NormalOrientation::OUTWARDS, const PolyhedronIntegrity &integrity = PolyhedronIntegrity::AUTOMATIC)#
Generates a polyhedron from nodes and faces.
Note
ASSERTS PRE-CONDITION that the in the indexing in the faces vector starts with zero!
- Parameters:
vertices – a vector of nodes
faces – a vector of faces containing the formation of faces off vertices
density – the density of the polyhedron (it must match the unit of the mesh, e.g., mesh in \([m]\) requires density in \([kg/m^3]\))
orientation – specify if the plane unit normals point outwards or inwards (defaults: to OUTWARDS)
integrity – specify if the mesh input is checked/ healed to fulfill the constraints of Tsoulis’ algorithm (see PolyhedronIntegrity)
- Throws:
std::invalid_argument – if no face contains the node zero indicating mathematical index
std::invalid_argument – dpending on the integrity flag
-
Polyhedron(const PolyhedralSource &polyhedralSource, double density, const NormalOrientation &orientation = NormalOrientation::OUTWARDS, const PolyhedronIntegrity &integrity = PolyhedronIntegrity::AUTOMATIC)#
Generates a polyhedron from nodes and faces.
Note
ASSERTS PRE-CONDITION that the in the indexing in the faces vector starts with zero!
- Parameters:
polyhedralSource – a tuple of vector containing the nodes and trianglular faces.
density – the density of the polyhedron (it must match the unit of the mesh, e.g., mesh in \([m]\) requires density in \([kg/m^3]\))
orientation – specify if the plane unit normals point outwards or inwards (defaults: to OUTWARDS)
integrity – specify if the mesh input is checked/ healed to fulfill the constraints of Tsoulis’ algorithm (see PolyhedronIntegrity)
- Throws:
std::invalid_argument – if no face contains the node zero indicating mathematical index
std::invalid_argument – dpending on the integrity flag
-
Polyhedron(const PolyhedralFiles &polyhedralFiles, double density, const NormalOrientation &orientation = NormalOrientation::OUTWARDS, const PolyhedronIntegrity &integrity = PolyhedronIntegrity::AUTOMATIC)#
Generates a polyhedron from nodes and faces.
Note
ASSERTS PRE-CONDITION that the in the indexing in the faces vector starts with zero!
- Parameters:
polyhedralFiles – a list of files (see TetgenAdapter
density – the density of the polyhedron (it must match the unit of the mesh, e.g., mesh in \([m]\) requires density in \([kg/m^3]\))
orientation – specify if the plane unit normals point outwards or inwards (defaults: to OUTWARDS)
integrity – specify if the mesh input is checked/ healed to fulfill the constraints of Tsoulis’ algorithm (see PolyhedronIntegrity)
- Throws:
std::invalid_argument – if no face contains the node zero indicating mathematical index
std::invalid_argument – dpending on the integrity flag
-
Polyhedron(const std::variant<PolyhedralSource, PolyhedralFiles> &polyhedralSource, double density, const NormalOrientation &orientation = NormalOrientation::OUTWARDS, const PolyhedronIntegrity &integrity = PolyhedronIntegrity::AUTOMATIC)#
Generates a polyhedron from nodes and faces. This constructor using a variant is maninly utilized from the Python Interface.
Note
ASSERTS PRE-CONDITION that the in the indexing in the faces vector starts with zero!
- Parameters:
polyhedralSource – a list of files (see TetgenAdapter or a tuple of vector containing the nodes and trianglular faces.
density – the density of the polyhedron (it must match the unit of the mesh, e.g., mesh in \([m]\) requires density in \([kg/m^3]\))
orientation – specify if the plane unit normals point outwards or inwards (defaults: to OUTWARDS)
integrity – specify if the mesh input is checked/ healed to fulfill the constraints of Tsoulis’ algorithm (see PolyhedronIntegrity)
- Throws:
std::invalid_argument – if no face contains the node zero indicating mathematical index
std::invalid_argument – dpending on the integrity flag
-
~Polyhedron() = default#
Default destructor
-
const std::vector<Array3> &getVertices() const#
Returns the vertices of this polyhedron
- Returns:
vector of cartesian coordinates
-
const Array3 &getVertex(size_t index) const#
Returns the vertex at a specific index
- Parameters:
index – size_t
- Returns:
cartesian coordinates of the vertex at index
-
size_t countVertices() const#
The number of points (nodes) that make up the polyhedron.
- Returns:
a size_t
-
const std::vector<IndexArray3> &getFaces() const#
Returns the triangular faces of this polyhedron
- Returns:
vector of triangular faces, where each element size_t references a vertex in the vertices vector
-
const IndexArray3 &getFace(size_t index) const#
Returns the indices of the vertices making up the face at the given index.
- Parameters:
index – size_t
- Returns:
triplet of the vertic’es indices forming the face
-
Array3Triplet getResolvedFace(size_t index) const#
Returns the resolved face with its concrete cartesian coordinates at the given index.
- Parameters:
index – size_t
- Returns:
triplet of vertices’ cartesian coordinates
-
size_t countFaces() const#
Returns the number of faces (triangles) that make up the polyhedral.
- Returns:
a size_t
-
double getDensity() const#
Returns the constant density of this polyhedron.
- Returns:
the constant density a double
-
void setDensity(double density)#
Sets the density to a new value. The density’s unit must match to the scaling of the mesh.
- Parameters:
density – the new constant density of the polyhedron
-
NormalOrientation getOrientation() const#
Returns the orientation of the plane unit normals of this polyhedron
- Returns:
OUTWARDS or INWARDS
-
double getOrientationFactor() const#
Retruns the plane unit normal orientation factor. If the unit normals are outwards pointing, it is 1.0 as Tsoulis inteneded. If the unit normals are inwards pointing, it is -1.0 (reversed).
- Returns:
1.0 or -1.0 depending on plane unit orientation
-
std::string toString() const#
Returns a string representation of the Polyhedron. Mainly used for the representation method in the Python interface.
- Returns:
string representation of the Polyhedron
-
std::tuple<std::vector<Array3>, std::vector<IndexArray3>, double, NormalOrientation> getState() const#
Returns the internal data strcuture of Python pickle support.
- Returns:
tuple of vertices, faces, density and normal orientation
-
inline auto transformIterator(const Array3 &offset = {0.0, 0.0, 0.0}) const#
An iterator transforming the polyhedron’s coordinates on demand by a given offset. This function returns a pair of transform iterators (first = begin(), second = end()).
- Parameters:
offset – the offset to apply
- Returns:
pair of transform iterators
-
std::pair<NormalOrientation, std::set<size_t>> checkPlaneUnitNormalOrientation() const#
This method determines the majority vertex ordering of a polyhedron and the set of faces which violate the majority constraint and need to be adpated. Hence, if the set is empty, all faces obey to the returned ordering/ plane unit normal orientation.
- Returns:
a pair consisting of majority ordering (OUTWARDS or INWARDS pointing normals) and a set of face indices which violate the constraint
-
Polyhedron(const std::vector<Array3> &vertices, const std::vector<IndexArray3> &faces, double density, const NormalOrientation &orientation = NormalOrientation::OUTWARDS, const PolyhedronIntegrity &integrity = PolyhedronIntegrity::AUTOMATIC)#
-
enum class polyhedralGravity::NormalOrientation : char#
The orientation of the plane unit normals of the polyhedron. We use this property as the concret definition of the vertices ordering depends on the utilized cooridnate system. However, the normal alignement is independent. Tsoulis et al. equations require the normals to point outwards of the polyhedron. If the opposite hold, the result is negated.
Values:
-
enumerator OUTWARDS#
Outwards pointing plane unit normals
-
enumerator INWARDS#
Inwards pointing plane unit normals
-
enumerator OUTWARDS#
-
enum class polyhedralGravity::PolyhedronIntegrity : char#
The three mode the poylhedron class takes in the constructor in order to determine what initilaization checks to conduct. This enum is exclusivly utilized in the constrcutor of a Polyhedron and its private method runIntegrityMeasures
Values:
-
enumerator DISABLE#
All activities regarding MeshChecking are disabled. No runtime overhead!
-
enumerator VERIFY#
Only verification of the NormalOrientation. A misalignment (e.g. specified OUTWARDS, but is not) leads to a runtime_error. Runtime Cost O(n^2)
-
enumerator AUTOMATIC#
Like VERIFY, but also informs the user about the option in any case on the runtime costs. This is the implicit default option. Runtime Cost: O(n^2) and output to stdout in every case!
-
enumerator HEAL#
Verification and Autmatioc Healing of the NormalOrientation. A misalignemt does not lead to a runtime_error, but to an internal correction. Runtime Cost: O(n^2) and a modification of the mesh input!
-
enumerator DISABLE#
GravityModel#
-
class GravityEvaluable#
Class for evaluating the polyhedrale gravity model for a given constant density polyhedron. Caches the polyhedron and data which is independent of the computation point P. Provides an operator() for evaluating the polyhedrale gravity model for a given constant density polyhedron at computation point P and choosing between parallel and serial evaluation.
Public Functions
-
inline explicit GravityEvaluable(const Polyhedron &polyhedron)#
Instantiates a GravityEvaluable with a given constant density polyhedron. In contrast to the GravityModel::evaluate, this evaluate method on the GravityEvaluable caches intermediate results and input data and subsequent evaluations will be faster.
- Parameters:
polyhedron – the constant density polyhedron
-
inline GravityEvaluable(const Polyhedron &polyhedron, const std::vector<Array3Triplet> &segmentVectors, const std::vector<Array3> &planeUnitNormals, const std::vector<Array3Triplet> &segmentUnitNormals)#
Instantiates a GravityEvaluable with a given constant density polyhedron and caches. This is for restoring a GravityEvaluable from a previous state.
- Parameters:
polyhedron – the polyhedron
segmentVectors – the segment vectors
planeUnitNormals – the plane unit normals
segmentUnitNormals – the segment unit normals
-
inline std::variant<GravityModelResult, std::vector<GravityModelResult>> operator()(const std::variant<Array3, std::vector<Array3>> &computationPoints, bool parallelization = true) const#
Evaluates the polyhedrale gravity model for a given constant density polyhedron at computation point P. Wrapper for evaluate<parallelization>.
- Parameters:
computationPoints – the computation point P or multiple computation points in a vector
parallelization – if true, the calculation is parallelized
- Returns:
the GravityModelResult containing the potential, acceleration and second derivative
-
std::string toString() const#
Returns a string representation of the GravityEvaluable.
- Returns:
string representation of the GravityEvaluable
-
std::tuple<Polyhedron, std::vector<Array3Triplet>, std::vector<Array3>, std::vector<Array3Triplet>> getState() const#
Returns the polyhedron, the density and the internal caches.
- Returns:
tuple of polyhedron, density, segmentVectors, planeUnitNormals and segmentUnitNormals
-
inline explicit GravityEvaluable(const Polyhedron &polyhedron)#
-
namespace GravityModel#
Namespace containing the methods used to evaluate the polyhedrale Gravity Model
Note
Naming scheme corresponds to the following: evaluate() –> main Method for evaluating the gravity model *() –> Methods calculating one property for the evaluation
Functions
-
GravityModelResult evaluate(const Polyhedron &polyhedron, const Array3 &computationPoint, bool parallel)#
Evaluates the polyhedrale gravity model for a given constant density polyhedron at computation point P.
- Parameters:
polyhedron – the polyhedron consisting of vertices and triangular faces
computationPoint – the computation Point P
parallel – whether to evaluate in parallel or serial
- Returns:
the GravityModelResult containing the potential, the acceleration and the change of acceleration at computation Point P
-
std::vector<GravityModelResult> evaluate(const Polyhedron &polyhedron, const std::vector<Array3> &computationPoints, bool parallel)#
Evaluates the polyhedral gravity model for a given constant density polyhedron at multiple computation points.
- Parameters:
polyhedron – the polyhedron consisting of vertices and triangular faces
computationPoints – vector of computation points
parallel – whether to evaluate in parallel or serial
- Returns:
the GravityModelResult containing the potential, the acceleration and the change of acceleration foreach computation Point P
-
namespace detail#
Functions
-
Array3Triplet buildVectorsOfSegments(const Array3 &vertex0, const Array3 &vertex1, const Array3 &vertex2)#
Computes the segment vectors G_ij for one plane of the polyhedron according to Tsoulis (18). The segment vectors G_ij represent the vector from one vertex of the face to the neighboring vertex and depict every line segment of the triangular face (A-B-C)
- Parameters:
vertex0 – the first vertex A
vertex1 – the second vertex B
vertex2 – the third vertex C
- Returns:
the segment vectors for a plane
-
Array3 buildUnitNormalOfPlane(const Array3 &segmentVector1, const Array3 &segmentVector2)#
Computes the plane unit normal N_p for one plane p of the polyhedron according to Tsoulis (19). The plane unit normal is the outward pointing normal of the face from the polyhedron.
- Parameters:
segmentVector1 – first edge
segmentVector2 – second edge
- Returns:
plane unit normal
-
Array3Triplet buildUnitNormalOfSegments(const Array3Triplet &segmentVectors, const Array3 &planeUnitNormal)#
Computes the segment unit normals n_pq for one plane p of the polyhedron according to Tsoulis (20). The segment unit normal n_pq represent the normal of one line segment of a polyhedrale face.
- Parameters:
segmentVectors – the segment vectors of the face G_p(0-2)
planeUnitNormal – the plane unit normal N_p
- Returns:
segment unit normals n_pq for plane p with q = {0, 1, 2}
-
double computeUnitNormalOfPlaneDirection(const Array3 &planeUnitNormal, const Array3 &vertex0)#
Computes the plane unit normal orientation/ direction sigma_p for one plane p of the polyhedron according to Tsoulis (21). The plane unit normal orientation values represents the relative position of computation point P with respect to the pointing direction of N_p. E. g. if N_p points to the half-space containing P, the inner product of N_p and -G_i1 will be positive, leading to a negative sigma_p. If sigma_p is zero than P and P’ lie geometrically in the same plane –> P == P’.
- Parameters:
planeUnitNormal – the plane unit normal N_p
vertex0 – the first vertex of the plane
- Returns:
plane normal orientation
-
HessianPlane computeHessianPlane(const Array3 &p, const Array3 &q, const Array3 &r)#
-
double distanceBetweenOriginAndPlane(const HessianPlane &hessianPlane)#
Calculates the (plane) distances h_p of computation point P to the plane S_p given in Hessian Form according to the following equation: h_p = D / sqrt(A^2+B^2+C^2)
- Parameters:
hessianPlane – Hessian Plane Form of S_p
- Returns:
plane distance h_p
-
Array3 projectPointOrthogonallyOntoPlane(const Array3 &planeUnitNormal, double planeDistance, const HessianPlane &hessianPlane)#
Computes P’ for a given plane p according to equation (22) of Tsoulis paper. P’ is the orthogonal projection of the computation point P onto the plane S_p.
- Parameters:
planeUnitNormal – the plane unit normal N_p
planeDistance – the distance from P to the plane h_p
hessianPlane – the Hessian Plane Form
- Returns:
P’ for this plane
-
Array3 computeUnitNormalOfSegmentsDirections(const Array3Triplet &vertices, const Array3 &projectionPointOnPlane, const Array3Triplet &segmentUnitNormalsForPlane)#
Computes the segment normal orientations/ directions sigma_pq for a given plane p. If sigma_pq is negative, this denotes that n_pq points to the half-plane containing P’. Nn case sigma_pq is positive, P’ resides in the other half-plane and if sigma_pq is zero, then P’ lies directly on the segment pq.
- Parameters:
vertices – the vertices of this plane
projectionPointOnPlane – the projection point P’ for this plane
segmentUnitNormalsForPlane – the segment unit normals sigma_pq for this plane
- Returns:
the segment normal orientations for the plane p
-
Array3Triplet projectPointOrthogonallyOntoSegments(const Array3 &projectionPointOnPlane, const Array3 &segmentNormalOrientations, const Array3Triplet &face)#
Computes the orthogonal projection Points P’’ foreach segment q of a given plane p.
- Parameters:
projectionPointOnPlane – the projection Point P’
segmentNormalOrientations – the segment normal orientations sigma_pq for this plane p
face – the vertices of the plane p
- Returns:
the orthogonal projection points of P on the segment P’’ foreach segment q of p
-
Array3 projectPointOrthogonallyOntoSegment(const Array3 &vertex1, const Array3 &vertex2, const Array3 &orthogonalProjectionPointOnPlane)#
Calculates the point P’’ for a given Segment consisting of vertices v1 and v2 and the orthogonal projection point P’ for the plane consisting of those vertices. Solves the three equations given in (24), (25) and (26).
Note
If sigma_pq is zero then P’’ == P’, this is not checked by this method, but has to be assured first
- Parameters:
vertex1 – first endpoint of segment
vertex2 – second endpoint of segment
orthogonalProjectionPointOnPlane – the orthogonal projection P’ of P on this plane
- Returns:
P’’ for this segment
-
Array3 distancesBetweenProjectionPoints(const Array3 &orthogonalProjectionPointOnPlane, const Array3Triplet &orthogonalProjectionPointOnSegments)#
Computes the (segment) distances h_pq between P’ for a given plane p and P’’ for a given segment q of plane p.
- Parameters:
orthogonalProjectionPointOnPlane – the orthogonal projection point P’ for p
orthogonalProjectionPointOnSegments – the orthogonal projection points P’’ for each segment q of p
- Returns:
distances h_pq for plane p
-
std::array<Distance, 3> distancesToSegmentEndpoints(const Array3Triplet &segmentVectorsForPlane, const Array3Triplet &orthogonalProjectionPointsOnSegmentForPlane, const Array3Triplet &face)#
Computes the 3D distances l1_pq and l2_pq between the computation point P and the line segment endpoints of each polyhedral segment for one plane. Computes the 1D distances s1_pq and s2_pq between orthogonal projection of P on the line segment P’’_pq and the line segment endpoints for each polyhedral segment for one plane
- Parameters:
segmentVectorsForPlane – the segment vectors G_pq for plane p
orthogonalProjectionPointsOnSegmentForPlane – the orthogonal projection Points P’’ for plane p
face – the vertices of plane p
- Returns:
distances l1_pq and l2_pq and s1_pq and s2_pq foreach segment q of plane p
-
std::array<TranscendentalExpression, 3> computeTranscendentalExpressions(const std::array<Distance, 3> &distancesForPlane, double planeDistance, const Array3 &segmentDistancesForPlane, const Array3 &segmentNormalOrientationsForPlane, const Array3 &projectionPointVertexNorms)#
Calculates the Transcendental Expressions LN_pq and AN_pq for every line segment of the polyhedron for a given plane p. LN_pq is calculated according to (14) using the natural logarithm and AN_pq is calculated according to (15) using the arctan.
- Parameters:
distancesForPlane – the distances l1, l2, s1, s2 foreach segment q of plane p
planeDistance – the plane distance h_p for plane p
segmentDistancesForPlane – the segment distance h_pq for segment q of plane p
segmentNormalOrientationsForPlane – the segment normal orientations n_pq for a plane p
orthogonalProjectionPointOnPlane – the orthogonal projection point P’ for plane p
face – the vertices of plane p
- Returns:
LN_pq and AN_pq foreach segment q of plane p
-
std::pair<double, Array3> computeSingularityTerms(const Array3Triplet &segmentVectorsForPlane, const Array3 &segmentNormalOrientationForPlane, const Array3 &projectionPointVertexNorms, const Array3 &planeUnitNormal, double planeDistance, double planeNormalOrientation)#
Calculates the singularities (correction) terms according to the Flow text for a given plane p.
- Parameters:
segmentVectorsForPlane – the segment vectors for a given plane
segmentNormalOrientationForPlane – the segment orientation sigma_pq
projectionPointVertexNorms – the projection point P’
planeUnitNormal – the plane unit normal N_p
planeDistance – the plane distance h_p
planeNormalOrientation – the plane normal orientation sigma_p
face – the vertices of plane p
- Returns:
the singularities for a plane p
-
Array3 computeNormsOfProjectionPointAndVertices(const Array3 &orthogonalProjectionPointOnPlane, const Array3Triplet &face)#
Computes the L2 norms of the orthogonal projection point P’ on a plane p with each vertex of that plane p. The values are later used to determine if P’ is situated at a vertex.
- Parameters:
orthogonalProjectionPointOnPlane – the orthogonal projection point P’
face – the vertices of plane p
- Returns:
the norms of p and each vertex
-
Array3Triplet buildVectorsOfSegments(const Array3 &vertex0, const Array3 &vertex1, const Array3 &vertex2)#
-
GravityModelResult evaluate(const Polyhedron &polyhedron, const Array3 &computationPoint, bool parallel)#
Named Tuple#
-
struct Distance#
Contains the 3D distances l1_pq and l2_pq between P and the endpoints of segment pq and the 1D distances s1_pq and s2_pq between P’’ and the segment endpoints.
Note
This struct is basically a named tuple
Public Functions
Public Members
-
double l1#
the 3D distance between computation point P and the first endpoint of line segment pq
-
double l2#
the 3D distance between computation point P and the second endpoint of line segment pq
-
double s1#
the 1D distance between projection of the computation point on line segment pq and the first endpoint of line segment pq
-
double s2#
the 1D distance between projection of the computation point on line segment pq and the second endpoint of line segment pq
-
double l1#
-
struct TranscendentalExpression#
Contains the Transcendental Expressions LN_pq and AN_pq for a given line segment pq of the polyhedron.
Note
This struct is basically a named tuple
Public Functions
-
inline bool operator==(const TranscendentalExpression &rhs) const#
Checks two TranscendentalExpressions for equality with another one by ensuring that the members are almost equal.
Note
Just used for testing purpose
- Parameters:
rhs – the other TranscendentalExpressions
- Returns:
true if equal
-
inline bool operator!=(const TranscendentalExpression &rhs) const#
Checks two TranscendentalExpressions for inequality with another one by ensuring that the members are not almost equal.
Note
Just used for testing purpose
- Parameters:
rhs – the other TranscendentalExpressions
- Returns:
false if unequal
Public Members
-
double ln#
The LN values for plane p and segment q of this plane is calculated in the following way: LN_pq = ln ((s_2_pq + l_2_pq) / (s_1_pq + l_1_pq))
Note
see Tsoulis Paper Equation (14)
-
double an#
The AN values for plane p and segment q of this plane is calculated in the following way: AN_pq = arctan ((h_p * s_2_pq) / (h_pq * l_2_pq)) - arctan ((h_pq * s_1_pq) / (h_pq * l_1_pq))
Note
see Tsoulis Paper Equation (15)
Friends
-
inline friend std::ostream &operator<<(std::ostream &os, const TranscendentalExpression &expression)#
Pretty output of this struct on the given ostream.
- Parameters:
os – the ostream
expression – a TranscendentalExpression
- Returns:
os
-
inline bool operator==(const TranscendentalExpression &rhs) const#
-
struct HessianPlane#
A struct describing a plane in Hessian Normal Form: ax + by + cz + d = 0 where a,b,c are the plane’s normal and d as the signed distance to the plane from the origin along the normal.
Public Functions
-
inline bool operator==(const HessianPlane &rhs) const#
Checking the equality of two this Hessian Plane with another one by ensuring that the members are almost equal.
Note
Just used for testing purpose
- Parameters:
rhs – other HessianPlane
- Returns:
true if equal
-
inline bool operator!=(const HessianPlane &rhs) const#
Checking the inequality of two this Hessian Plane with another one by ensuring that the members are not almost equal.
Note
Just used for testing purpose
- Parameters:
rhs – other HessianPlane
- Returns:
true if unequal
Public Members
-
double a#
part of the planes normal [a, b, c]
-
double b#
part of the panes normal [a, b, c]
-
double c#
part of the planes normal [a, b, c]
-
double d#
the signed distance to the plane from the origin along the normal
Friends
-
inline friend std::ostream &operator<<(std::ostream &os, const HessianPlane &hessianPlane)#
Pretty output of this struct on the given ostream.
- Parameters:
os – the ostream
hessianPlane – a HessianPlane
- Returns:
os
-
inline bool operator==(const HessianPlane &rhs) const#