Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
Plane3D.h
Go to the documentation of this file.
1// @(#)root/mathcore:$Id$
2// Authors: L. Moneta 12/2005
3
4/**********************************************************************
5 * *
6 * Copyright (c) 2005 , LCG ROOT MathLib Team *
7 * *
8 * *
9 **********************************************************************/
10
11// Header file for class LorentzVector
12//
13// Created by: moneta at Fri Dec 02 2005
14//
15// Last update: $Id$
16//
17#ifndef ROOT_Math_GenVector_Plane3D
18#define ROOT_Math_GenVector_Plane3D 1
19
20#include <type_traits>
21
24
25
26
27namespace ROOT {
28
29namespace Math {
30
31namespace Impl {
32
33//_______________________________________________________________________________
34/**
35 Class describing a geometrical plane in 3 dimensions.
36 A Plane3D is a 2 dimensional surface spanned by two linearly independent vectors.
37 The plane is described by the equation
38 \f$ a*x + b*y + c*z + d = 0 \f$ where (a,b,c) are the components of the
39 normal vector to the plane \f$ n = (a,b,c) \f$ and \f$ d = - n \dot x \f$, where x is any point
40 belonging to plane.
41 More information on the mathematics describing a plane in 3D is available on
42 <A HREF=http://mathworld.wolfram.com/Plane.html>MathWord</A>.
43 The Plane3D class contains the 4 scalar values in T which represent the
44 four coefficients, fA, fB, fC, fD. fA, fB, fC are the normal components normalized to 1,
45 i.e. fA**2 + fB**2 + fC**2 = 1
46
47 @ingroup GenVector
48
49 @sa Overview of the @ref GenVector "physics vector library"
50*/
51
52template <typename T = double>
53class Plane3D {
54
55public:
56 // ------ ctors ------
57
58 typedef T Scalar;
59
62
63 /**
64 default constructor create plane z = 0
65 */
66 Plane3D() : fA(0), fB(0), fC(1), fD(0) {}
67
68 /**
69 generic constructors from the four scalar values describing the plane
70 according to the equation ax + by + cz + d = 0
71 \param a scalar value
72 \param b scalar value
73 \param c scalar value
74 \param d sxcalar value
75 */
76 Plane3D(const Scalar &a, const Scalar &b, const Scalar &c, const Scalar &d) : fA(a), fB(b), fC(c), fD(d)
77 {
78 // renormalize a,b,c to unit
79 Normalize();
80 }
81
82 /**
83 constructor a Plane3D from a normal vector and a point coplanar to the plane
84 \param n normal expressed as a ROOT::Math::DisplacementVector3D<Cartesian3D<T> >
85 \param p point expressed as a ROOT::Math::PositionVector3D<Cartesian3D<T> >
86 */
87 Plane3D(const Vector &n, const Point &p) { BuildFromVecAndPoint(n, p); }
88
89 /**
90 Construct from a generic DisplacementVector3D (normal vector) and PositionVector3D (point coplanar to
91 the plane)
92 \param n normal expressed as a generic ROOT::Math::DisplacementVector3D
93 \param p point expressed as a generic ROOT::Math::PositionVector3D
94 */
95 template <class T1, class T2, class U>
97 {
99 }
100
101 /**
102 constructor from three Cartesian point belonging to the plane
103 \param p1 point1 expressed as a generic ROOT::Math::PositionVector3D
104 \param p2 point2 expressed as a generic ROOT::Math::PositionVector3D
105 \param p3 point3 expressed as a generic ROOT::Math::PositionVector3D
106 */
107 Plane3D(const Point &p1, const Point &p2, const Point &p3) { BuildFrom3Points(p1, p2, p3); }
108
109 /**
110 constructor from three generic point belonging to the plane
111 \param p1 point1 expressed as ROOT::Math::DisplacementVector3D<Cartesian3D<T> >
112 \param p2 point2 expressed as ROOT::Math::DisplacementVector3D<Cartesian3D<T> >
113 \param p3 point3 expressed as ROOT::Math::DisplacementVector3D<Cartesian3D<T> >
114 */
115 template <class T1, class T2, class T3, class U>
117 {
118 BuildFrom3Points(Point(p1.X(), p1.Y(), p1.Z()), Point(p2.X(), p2.Y(), p2.Z()), Point(p3.X(), p3.Y(), p3.Z()));
119 }
120
121 // compiler-generated copy ctor and dtor are fine.
122 Plane3D(const Plane3D &) = default;
123
124 // ------ assignment ------
125
126 /**
127 Assignment operator from other Plane3D class
128 */
129 Plane3D &operator=(const Plane3D &) = default;
130
131 /**
132 Return the a coefficient of the plane equation \f$ a*x + b*y + c*z + d = 0 \f$. It is also the
133 x-component of the vector perpendicular to the plane.
134 */
135 Scalar A() const { return fA; }
136
137 /**
138 Return the b coefficient of the plane equation \f$ a*x + b*y + c*z + d = 0 \f$. It is also the
139 y-component of the vector perpendicular to the plane
140 */
141 Scalar B() const { return fB; }
142
143 /**
144 Return the c coefficient of the plane equation \f$ a*x + b*y + c*z + d = 0 \f$. It is also the
145 z-component of the vector perpendicular to the plane
146 */
147 Scalar C() const { return fC; }
148
149 /**
150 Return the d coefficient of the plane equation \f$ a*x + b*y + c*z + d = 0 \f$. It is also
151 the distance from the origin (HesseDistance)
152 */
153 Scalar D() const { return fD; }
154
155 /**
156 Return normal vector to the plane as Cartesian DisplacementVector
157 */
158 Vector Normal() const { return Vector(fA, fB, fC); }
159
160 /**
161 Return the Hesse Distance (distance from the origin) of the plane or
162 the d coefficient expressed in normalize form
163 */
164 Scalar HesseDistance() const { return fD; }
165
166 /**
167 Return the signed distance to a Point.
168 The distance is signed positive if the Point is in the same side of the
169 normal vector to the plane.
170 \param p Point expressed in Cartesian Coordinates
171 */
172 Scalar Distance(const Point &p) const { return fA * p.X() + fB * p.Y() + fC * p.Z() + fD; }
173
174 /**
175 Return the distance to a Point described with generic coordinates
176 \param p Point expressed as generic ROOT::Math::PositionVector3D
177 */
178 template <class T1, class U>
180 {
181 return Distance(Point(p.X(), p.Y(), p.Z()));
182 }
183
184 /**
185 Return the projection of a Cartesian point to a plane
186 \param p Point expressed as PositionVector3D<Cartesian3D<T> >
187 */
189 {
190 const Scalar d = Distance(p);
191 return XYZPoint(p.X() - fA * d, p.Y() - fB * d, p.Z() - fC * d);
192 }
193
194 /**
195 Return the projection of a point to a plane
196 \param p Point expressed as generic ROOT::Math::PositionVector3D
197 */
198 template <class T1, class U>
200 {
201 const Point pxyz = ProjectOntoPlane(Point(p.X(), p.Y(), p.Z()));
202 return PositionVector3D<T, U>(pxyz.X(), pxyz.Y(), pxyz.Z());
203 }
204
205 // ------------------- Equality -----------------
206
207 /**
208 Exact equality
209 */
210 bool operator==(const Plane3D &rhs) const { return (fA == rhs.fA && fB == rhs.fB && fC == rhs.fC && fD == rhs.fD); }
211 bool operator!=(const Plane3D &rhs) const { return !(operator==(rhs)); }
212
213protected:
214 /**
215 Normalize the normal (a,b,c) plane components
216 */
217 template <typename SCALAR = T, typename std::enable_if<std::is_arithmetic<SCALAR>::value>::type * = nullptr>
219 {
220 // normalize the plane
221 using std::sqrt;
222 const SCALAR s = sqrt(fA * fA + fB * fB + fC * fC);
223 // what to do if s = 0 ?
224 if (s == SCALAR(0)) {
225 fD = SCALAR(0);
226 } else {
227 const SCALAR w = Scalar(1) / s;
228 fA *= w;
229 fB *= w;
230 fC *= w;
231 fD *= w;
232 }
233 }
234
235 /**
236 Normalize the normal (a,b,c) plane components
237 */
238 template <typename SCALAR = T, typename std::enable_if<!std::is_arithmetic<SCALAR>::value>::type * = nullptr>
240 {
241 // normalize the plane
242 using std::sqrt;
243 SCALAR s = sqrt(fA * fA + fB * fB + fC * fC);
244 // what to do if s = 0 ?
245 const auto m = (s == SCALAR(0));
246 // set zero entries to 1 in the vector to avoid /0 later on
247 s(m) = SCALAR(1);
248 fD(m) = SCALAR(0);
249 const SCALAR w = SCALAR(1) / s;
250 fA *= w;
251 fB *= w;
252 fC *= w;
253 fD *= w;
254 }
255
256private:
257 // internal method to construct class from a vector and a point
258 void BuildFromVecAndPoint(const Vector &n, const Point &p)
259 {
260 // build from a normal vector and a point
261 fA = n.X();
262 fB = n.Y();
263 fC = n.Z();
264 fD = -n.Dot(p);
265 Normalize();
266 }
267
268 // internal method to construct class from 3 points
269 void BuildFrom3Points(const Point &p1, const Point &p2, const Point &p3)
270 {
271 // plane from three points
272 // normal is (x3-x1) cross (x2 -x1)
273 const Vector n = (p2 - p1).Cross(p3 - p1);
274 fA = n.X();
275 fB = n.Y();
276 fC = n.Z();
277 fD = -n.Dot(p1);
278 Normalize();
279 }
280
281 // plane data members the four scalar which satisfies fA*x + fB*y + fC*z + fD = 0
282 // for every point (x,y,z) belonging to the plane.
283 // fA**2 + fB**2 + fC** =1 plane is stored in normalized form
288
289 }; // Plane3D<>
290
291 /**
292 Stream Output and Input
293 */
294 // TODO - I/O should be put in the manipulator form
295 template <typename T>
296 std::ostream &operator<<(std::ostream &os, const Plane3D<T> &p)
297 {
298 os << "\n"
299 << p.Normal().X() << " " << p.Normal().Y() << " " << p.Normal().Z() << " " << p.HesseDistance() << "\n";
300 return os;
301 }
302
303 } // end namespace Impl
304
305 // typedefs for double and float versions
308
309} // end namespace Math
310
311} // end namespace ROOT
312
313
314#endif
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
TBuffer & operator<<(TBuffer &buf, const Tmpl *obj)
Definition TBuffer.h:399
winID h TVirtualViewer3D TVirtualGLPainter p
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
DefaultCoordinateSystemTag Default tag for identifying any coordinate system.
Class describing a generic displacement vector in 3 dimensions.
Scalar X() const
Cartesian X, converting if necessary from internal coordinate system.
Class describing a geometrical plane in 3 dimensions.
Definition Plane3D.h:53
Plane3D(const Scalar &a, const Scalar &b, const Scalar &c, const Scalar &d)
generic constructors from the four scalar values describing the plane according to the equation ax + ...
Definition Plane3D.h:76
Vector Normal() const
Return normal vector to the plane as Cartesian DisplacementVector.
Definition Plane3D.h:158
Plane3D & operator=(const Plane3D &)=default
Assignment operator from other Plane3D class.
Plane3D()
default constructor create plane z = 0
Definition Plane3D.h:66
bool operator!=(const Plane3D &rhs) const
Definition Plane3D.h:211
DisplacementVector3D< Cartesian3D< T >, DefaultCoordinateSystemTag > Vector
Definition Plane3D.h:60
Point ProjectOntoPlane(const Point &p) const
Return the projection of a Cartesian point to a plane.
Definition Plane3D.h:188
PositionVector3D< T1, U > ProjectOntoPlane(const PositionVector3D< T1, U > &p) const
Return the projection of a point to a plane.
Definition Plane3D.h:199
Scalar HesseDistance() const
Return the Hesse Distance (distance from the origin) of the plane or the d coefficient expressed in n...
Definition Plane3D.h:164
Scalar C() const
Return the c coefficient of the plane equation .
Definition Plane3D.h:147
Plane3D(const Vector &n, const Point &p)
constructor a Plane3D from a normal vector and a point coplanar to the plane
Definition Plane3D.h:87
Plane3D(const PositionVector3D< T1, U > &p1, const PositionVector3D< T2, U > &p2, const PositionVector3D< T3, U > &p3)
constructor from three generic point belonging to the plane
Definition Plane3D.h:116
Scalar A() const
Return the a coefficient of the plane equation .
Definition Plane3D.h:135
Scalar D() const
Return the d coefficient of the plane equation .
Definition Plane3D.h:153
Scalar Distance(const Point &p) const
Return the signed distance to a Point.
Definition Plane3D.h:172
Scalar B() const
Return the b coefficient of the plane equation .
Definition Plane3D.h:141
void Normalize()
Normalize the normal (a,b,c) plane components.
Definition Plane3D.h:218
PositionVector3D< Cartesian3D< T >, DefaultCoordinateSystemTag > Point
Definition Plane3D.h:61
bool operator==(const Plane3D &rhs) const
Exact equality.
Definition Plane3D.h:210
void BuildFromVecAndPoint(const Vector &n, const Point &p)
Definition Plane3D.h:258
Plane3D(const Plane3D &)=default
Plane3D(const DisplacementVector3D< T1, U > &n, const PositionVector3D< T2, U > &p)
Construct from a generic DisplacementVector3D (normal vector) and PositionVector3D (point coplanar to...
Definition Plane3D.h:96
void BuildFrom3Points(const Point &p1, const Point &p2, const Point &p3)
Definition Plane3D.h:269
Plane3D(const Point &p1, const Point &p2, const Point &p3)
constructor from three Cartesian point belonging to the plane
Definition Plane3D.h:107
Scalar Distance(const PositionVector3D< T1, U > &p) const
Return the distance to a Point described with generic coordinates.
Definition Plane3D.h:179
Class describing a generic position vector (point) in 3 dimensions.
Scalar Y() const
Cartesian Y, converting if necessary from internal coordinate system.
Scalar Z() const
Cartesian Z, converting if necessary from internal coordinate system.
Scalar X() const
Cartesian X, converting if necessary from internal coordinate system.
SVector< T, 3 > Cross(const SVector< T, 3 > &lhs, const SVector< T, 3 > &rhs)
Vector Cross Product (only for 3-dim vectors) .
Definition Functions.h:323
const Int_t n
Definition legend1.C:16
Namespace for new Math classes and functions.
Impl::Plane3D< double > Plane3D
Definition Plane3D.h:306
PositionVector3D< Cartesian3D< double >, DefaultCoordinateSystemTag > XYZPoint
3D Point based on the cartesian coordinates x,y,z in double precision
Definition Point3Dfwd.h:38
VecExpr< UnaryOp< Sqrt< T >, VecExpr< A, T, D >, T >, T, D > sqrt(const VecExpr< A, T, D > &rhs)
Impl::Plane3D< float > Plane3DF
Definition Plane3D.h:307
This file contains a specialised ROOT message handler to test for diagnostic in unit tests.
TMarker m
Definition textangle.C:8