12#include <framework/logging/Logger.h>
15#include <Math/Vector3D.h>
41 template<
typename DataType>
45 static_assert(std::is_floating_point<DataType>::value,
"B2Vector3 only works with floating point types");
62 B2Vector3(
const TVector3& tVec3):
m_coordinates {static_cast<DataType>(tVec3.
X()), static_cast<DataType>(tVec3.
Y()), static_cast<DataType>(tVec3.
Z())} {};
65 B2Vector3(
const TVector3* tVec3):
m_coordinates {static_cast<DataType>(tVec3->
X()), static_cast<DataType>(tVec3->
Y()), static_cast<DataType>(tVec3->
Z())} {};
73 m_coordinates {static_cast<DataType>(b2Vec3.
X()), static_cast<DataType>(b2Vec3.
Y()), static_cast<DataType>(b2Vec3.
Z())} {};
76 m_coordinates {static_cast<DataType>(b2Vec3->
X()), static_cast<DataType>(b2Vec3->
Y()), static_cast<DataType>(b2Vec3->
Z())} {};
79 B2Vector3(
const ROOT::Math::XYZVector& xyzVec):
m_coordinates {static_cast<DataType>(xyzVec.
X()), static_cast<DataType>(xyzVec.
Y()), static_cast<DataType>(xyzVec.
Z())} {};
82 B2Vector3(
const ROOT::Math::XYZVector* xyzVec):
m_coordinates {static_cast<DataType>(xyzVec->
X()), static_cast<DataType>(xyzVec->
Y()), static_cast<DataType>(xyzVec->
Z())} {};
108 bool operator == (
const TVector3& b)
const {
return X() == b.X() &&
Y() == b.Y() &&
Z() == b.Z(); }
110 bool operator == (
const ROOT::Math::XYZVector& b)
const {
return X() == b.X() &&
Y() == b.Y() &&
Z() == b.Z(); }
114 bool operator != (
const TVector3& b)
const {
return !(*
this == b); }
116 bool operator != (
const ROOT::Math::XYZVector& b)
const {
return !(*
this == b); }
151 DataType
Phi()
const {
return X() == 0 &&
Y() == 0 ? 0 : atan2(
Y(),
X()); }
153 DataType
Theta()
const {
return X() == 0 &&
Y() == 0 &&
Z() == 0 ? 0 : atan2(
Perp(),
Z()); }
155 DataType
CosTheta()
const {
const double pTot =
Mag();
return pTot == 0 ? 1 :
Z() / pTot; }
157 DataType
Mag2()
const {
return X() *
X() +
Y() *
Y() +
Z() *
Z(); }
159 DataType
Mag()
const {
return std::hypot((
double)
Perp(), (
double)
Z()); }
164 const double perp =
Perp();
165 SetX(perp * cos((
double)phi));
166 SetY(perp * sin((
double)phi));
172 const double ma =
Mag();
173 const double ph =
Phi();
174 const double ctheta = std::cos((
double) theta);
175 const double stheta = std::sin((
double) theta);
176 SetX(ma * stheta * std::cos(ph));
177 SetY(ma * stheta * std::cos(ph));
184 double factor =
Mag();
186 B2WARNING(
name() <<
"::SetMag: zero vector can't be stretched");
188 factor = mag / factor;
196 DataType
Perp2()
const {
return X() *
X() +
Y() *
Y(); }
200 DataType
Perp()
const {
return std::hypot((
double)
X(), (
double)
Y()); }
205 const double p =
Perp();
215 const double tot = axis.Mag2();
216 const double ss =
Dot(axis);
218 if (tot > 0.0) per -= ss * ss / tot;
219 if (per < 0) per = 0;
234 if (std::isnan(angle)) {
235 B2ERROR(
name() <<
"::Mpi_pi: function called with NaN");
238 angle = std::remainder(angle, 2 * M_PI);
240 if (angle == M_PI) angle = -M_PI;
247 const double deta =
Eta() - v.
Eta();
249 return std::hypot(deta, dphi);
261 const double amag = std::abs(mag);
262 const double sinTheta = std::sin((
double)theta);
271 const double tot =
Mag2();
273 return tot > 0.0 ? p *= (1.0 / std::sqrt(tot)) : p;
279 const double xVal = std::abs((
double)
X());
280 const double yVal = std::abs((
double)
Y());
281 const double zVal = std::abs((
double)
Z());
292 return X() * p.X() +
Y() * p.Y() +
Z() * p.Z();
304 const double ptot2 =
Mag2() * q.
Mag2();
308 double arg =
Dot(q) / std::sqrt(ptot2);
309 if (arg > 1.0) arg = 1.0;
310 if (arg < -1.0) arg = -1.0;
311 return std::acos(arg);
322 if (std::abs(cosTheta) < 1)
return -0.5 * std::log((1.0 - cosTheta) / (1.0 + cosTheta));
323 if (
Z() == 0)
return 0;
325 if (
Z() > 0)
return 10e10;
338 const double s = std::sin((
double)angle);
339 const double c = std::cos((
double)angle);
340 const double yOld =
Y();
350 const double s = std::sin((
double)angle);
351 const double c = std::cos((
double)angle);
352 const double zOld =
Z();
362 const double s = std::sin((
double)angle);
363 const double c = std::cos((
double)angle);
364 const double xOld =
X();
374 const double u1 = NewUzVector.
X();
375 const double u2 = NewUzVector.
Y();
376 const double u3 = NewUzVector.
Z();
377 double up = u1 * u1 + u2 * u2;
381 DataType px =
X(), py =
Y(), pz =
Z();
382 m_coordinates[0] = (u1 * u3 * px - u2 * py + u1 * up * pz) / up;
383 m_coordinates[1] = (u2 * u3 * px + u1 * py + u2 * up * pz) / up;
385 }
else if (u3 < 0.) {
402 *
this = (n * (n.Dot(*
this)) + cos(alpha) * ((n.Cross(*
this)).
Cross(n)) + sin(alpha) * (n.Cross(*
this)));
423 DataType
at(
unsigned i)
const;
431 DataType
X()
const {
return x(); }
433 DataType
Y()
const {
return y(); }
435 DataType
Z()
const {
return z(); }
437 DataType
Px()
const {
return x(); }
439 DataType
Py()
const {
return y(); }
441 DataType
Pz()
const {
return z(); }
450 void GetXYZ(ROOT::Math::XYZVector* xyzVec)
const;
473 void SetXYZ(
const ROOT::Math::XYZVector& xyzVec);
475 void SetXYZ(
const ROOT::Math::XYZVector* xyzVec);
489 std::ostringstream output;
490 output <<
"(x,y,z)=("
491 << std::fixed << std::setprecision(precision)
492 <<
X() <<
"," <<
Y() <<
"," <<
Z() <<
")";
499 std::ostringstream output;
500 output <<
"(rho, theta, phi)=("
501 << std::fixed << std::setprecision(precision)
502 <<
Mag() <<
"," <<
Theta() * 180. / M_PI <<
"," <<
Phi() * 180. / M_PI <<
")";
522 template <
typename DataType>
525 return (a.X() == b.X() && a.Y() == b.Y() && a.Z() == b.Z());
529 template <
typename DataType>
536 template <
typename DataType>
543 template <
typename DataType>
550 template <
typename DataType>
557 template <
typename DataType>
564 template <
typename DataType>
571 template <
typename DataType>
578 template <
typename DataType>
585 template <
typename DataType>
592 template <
typename DataType>
600 template<
typename DataType >
603 m_coordinates[0] = b.X();
604 m_coordinates[1] = b.Y();
605 m_coordinates[2] = b.Z();
610 template<
typename DataType >
613 m_coordinates[0] = b.X();
614 m_coordinates[1] = b.Y();
615 m_coordinates[2] = b.Z();
620 template<
typename DataType >
623 m_coordinates[0] = b.X();
624 m_coordinates[1] = b.Y();
625 m_coordinates[2] = b.Z();
630 template<
typename DataType >
633 m_coordinates[0] += b.X();
634 m_coordinates[1] += b.Y();
635 m_coordinates[2] += b.Z();
641 template<
typename DataType >
644 m_coordinates[0] -= b.X();
645 m_coordinates[1] -= b.Y();
646 m_coordinates[2] -= b.Z();
651 template<
typename DataType >
654 m_coordinates[0] *= a;
655 m_coordinates[1] *= a;
656 m_coordinates[2] *= a;
661 template<
typename DataType >
664 m_coordinates[0] =
static_cast<Double_t
>(tVec.X());
665 m_coordinates[1] =
static_cast<Double_t
>(tVec.Y());
666 m_coordinates[2] =
static_cast<Double_t
>(tVec.Z());
670 template<
typename DataType >
673 m_coordinates[0] =
static_cast<Double_t
>(tVec->X());
674 m_coordinates[1] =
static_cast<Double_t
>(tVec->Y());
675 m_coordinates[2] =
static_cast<Double_t
>(tVec->Z());
679 template<
typename DataType >
682 m_coordinates[0] =
static_cast<Double_t
>(xyzVec.X());
683 m_coordinates[1] =
static_cast<Double_t
>(xyzVec.Y());
684 m_coordinates[2] =
static_cast<Double_t
>(xyzVec.Z());
688 template<
typename DataType >
691 m_coordinates[0] =
static_cast<Double_t
>(xyzVec->X());
692 m_coordinates[1] =
static_cast<Double_t
>(xyzVec->Y());
693 m_coordinates[2] =
static_cast<Double_t
>(xyzVec->Z());
696 template<
typename DataType >
705 template<
typename DataType >
708 tVec->SetXYZ(
static_cast<Double_t
>(X()),
709 static_cast<Double_t
>(Y()),
710 static_cast<Double_t
>(Z()));
714 template<
typename DataType >
717 xyzVec->SetXYZ(
static_cast<Double_t
>(X()),
718 static_cast<Double_t
>(Y()),
719 static_cast<Double_t
>(Z()));
724 template<
typename DataType >
729 static_cast<Double_t
>(X()),
730 static_cast<Double_t
>(Y()),
731 static_cast<Double_t
>(Z())
737 template<
typename DataType >
741 ROOT::Math::XYZVector(
742 static_cast<Double_t
>(X()),
743 static_cast<Double_t
>(Y()),
744 static_cast<Double_t
>(Z())
750 template <
typename DataType>
761 B2FATAL(this->name() <<
"::access operator: given index (i=" << i <<
") is out of bounds!");
766 template <
typename DataType>
769 return std::string(
"B2Vector3<") +
typeid(DataType).name() + std::string(
">");
A fast and root compatible alternative to TVector3.
DataType Phi() const
The azimuth angle.
void SetMag(DataType mag)
Set magnitude keeping theta and phi constant.
B2Vector3< DataType > operator-() const
unary minus
DataType Pz() const
access variable Z (= .at(2) without boundary check)
DataType m_coordinates[3]
Make sure that we only have floating point vectors.
B2Vector3(const B2Vector3< DataType > &b2Vec3)
Constructor expecting a B2Vector3 of same type.
DataType operator()(unsigned i) const
member access without boundary check
DataType Z() const
access variable Z (= .at(2) without boundary check)
DataType Perp2() const
The transverse component squared (R^2 in cylindrical coordinate system).
void SetPerp(DataType r)
Set the transverse component keeping phi and z constant.
DataType y() const
access variable Y (= .at(1) without boundary check)
void SetX(DataType x)
set X/1st-coordinate
DataType Theta() const
The polar angle.
DataType CosTheta() const
Cosine of the polar angle.
B2Vector3< DataType > Cross(const B2Vector3< DataType > &p) const
Cross product.
void GetXYZ(Float_t *carray) const
directly copies coordinates to an array of float
DataType z() const
access variable Z (= .at(2) without boundary check)
B2Vector3(const B2Vector3< DataType > *b2Vec3)
Constructor expecting a pointer to a B2Vector3.
void SetMagThetaPhi(DataType mag, DataType theta, DataType phi)
setter with mag, theta, phi
DataType X() const
access variable X (= .at(0) without boundary check)
std::string PrintStringXYZ(unsigned precision=4) const
create a string containing vector in cartesian coordinates
B2Vector3< DataType > Orthogonal() const
Vector orthogonal to this one.
DataType Eta() const
Returns the pseudo-rapidity.
B2Vector3< DataType > operator*(DataType a) const
Scaling of 3-vectors with a real number.
DataType DeltaPhi(const B2Vector3< DataType > &v) const
returns phi in the interval [-PI,PI)
void RotateY(DataType angle)
Rotates the B2Vector3 around the y-axis.
DataType Y() const
access variable Y (= .at(1) without boundary check)
DataType operator[](unsigned i) const
member access without boundary check
std::string PrintString(unsigned precision=4) const
create a string containing vector in cartesian and spherical coordinates
DataType Mag() const
The magnitude (rho in spherical coordinate system).
DataType value_type
storage type of the vector
B2Vector3(const ROOT::Math::XYZVector &xyzVec)
Constructor expecting a XYZVector.
bool operator==(const B2Vector3< DataType > &b) const
Comparison for equality with a B2Vector3.
DataType x() const
access variable X (= .at(0) without boundary check)
B2Vector3(const DataType(*coords)[3])
Constructor using a pointer.
DataType DeltaR(const B2Vector3< DataType > &v) const
return deltaR with respect to input-vector
DataType Perp(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
DataType & operator()(unsigned i)
member access without boundary check
void RotateX(DataType angle)
Rotates the B2Vector3 around the x-axis.
void Abs()
calculates the absolute value of the coordinates element-wise
B2Vector3(const DataType(&coords)[3])
Constructor using a reference.
B2Vector3(const ROOT::Math::XYZVector *xyzVec)
Constructor expecting a pointer to a XYZVector.
DataType Pt(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
DataType Mag2() const
The magnitude squared (rho^2 in spherical coordinate system).
void SetTheta(DataType theta)
Set theta keeping mag and phi constant.
B2Vector3< DataType > operator+(const B2Vector3< DataType > &b) const
Addition of 3-vectors.
B2Vector3(const TVector3 &tVec3)
Constructor expecting a TVector3.
void RotateZ(DataType angle)
Rotates the B2Vector3 around the z-axis.
void Print()
just for backward compatibility, should not be used with new code
void Sqrt()
calculates the square root of the absolute values of the coordinates element-wise
void SetZ(DataType z)
set Z/3rd-coordinate
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
void SetY(DataType y)
set Y/2nd-coordinate
B2Vector3(const B2Vector3< OtherType > *b2Vec3)
Constructor expecting a pointer to a B2Vector3 of different type.
DataType PseudoRapidity() const
Returns the pseudo-rapidity, i.e.
DataType DrEtaPhi(const B2Vector3< DataType > &v) const
return DrEtaPhi with respect to input-vector
bool operator!=(const B2Vector3< DataType > &b) const
Comparison != with a B2Vector3.
DataType Perp() const
The transverse component (R in cylindrical coordinate system).
B2Vector3(void)
empty Constructor sets everything to 0
DataType Perp2(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
B2Vector3(const B2Vector3< OtherType > &b2Vec3)
Constructor expecting a B2Vector3 of different type.
void GetXYZ(Double_t *carray) const
directly copies coordinates to an array of double
void RotateUz(const B2Vector3< DataType > &NewUzVector)
Rotates reference frame from Uz to newUz (unit vector).
DataType Py() const
access variable Y (= .at(1) without boundary check)
static DataType Mpi_pi(DataType angle)
returns given angle in the interval [-PI,PI)
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
DataType Pt() const
The transverse component (R in cylindrical coordinate system).
B2Vector3(const TVector3 *tVec3)
Constructor expecting a pointer to a TVector3.
DataType Px() const
access variable X (= .at(0) without boundary check)
void SetXYZ(DataType x, DataType y, DataType z)
set all coordinates using data type
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
DataType & operator[](unsigned i)
member access without boundary check
std::string PrintStringCyl(unsigned precision=4) const
create a string containing vector in spherical coordinates
B2Vector3< DataType > operator/(DataType a) const
Scaling of 3-vectors with a real number.
void Rotate(DataType alpha, const B2Vector3< DataType > &v)
Rotation around an arbitrary axis v with angle alpha.
void SetPhi(DataType phi)
Set phi keeping mag and theta constant.
B2Vector3(const DataType xVal, const DataType yVal, const DataType zVal)
Constructor expecting 3 coordinates.
bool operator==(const DecayNode &node1, const DecayNode &node2)
Compare two Decay Nodes: They are equal if All daughter decay nodes are equal or one of the daughter ...
bool operator!=(const DecayNode &node1, const DecayNode &node2)
Not equal: See operator==.
void SetXYZ(const ROOT::Math::XYZVector *xyzVec)
set all coordinates using a pointer to XYZVector
void GetXYZ(ROOT::Math::XYZVector *xyzVec) const
directly copies coordinates to a XYZVector
void SetXYZ(const TVector3 &tVec)
set all coordinates using a reference to TVector3
B2Vector3< DataType > & operator-=(const B2Vector3< DataType > &b)
subtraction
TVector3 GetTVector3() const
returns a TVector3 containing the same coordinates
B2Vector3< float > B2Vector3F
typedef for common usage with float
void GetXYZ(TVector3 *tVec) const
directly copies coordinates to a TVector3
void SetXYZ(const TVector3 *tVec)
set all coordinates using a pointer to TVector3
ROOT::Math::XYZVector GetXYZVector() const
returns a XYZVector containing the same coordinates
B2Vector3< DataType > & operator+=(const B2Vector3< DataType > &b)
addition
B2Vector3< DataType > operator*(DataType a, const B2Vector3< DataType > &p)
non-memberfunction Scaling of 3-vectors with a real number
B2Vector3< DataType > operator-(const TVector3 &a, const B2Vector3< DataType > &b)
non-memberfunction for substracting a TVector3 from a B2Vector3
B2Vector3< DataType > & operator*=(DataType a)
scaling with real numbers
B2Vector3< DataType > operator+(const TVector3 &a, const B2Vector3< DataType > &b)
non-memberfunction for adding a TVector3 to a B2Vector3
void SetXYZ(const ROOT::Math::XYZVector &xyzVec)
set all coordinates using a reference to XYZVector
B2Vector3< DataType > & operator=(const B2Vector3< DataType > &b)
Assignment via B2Vector3.
B2Vector3< double > B2Vector3D
typedef for common usage with double
DataType at(unsigned i) const
safe member access (with boundary check!)
static std::string name()
Returns the name of the B2Vector.
Abstract base class for different kinds of events.