Belle II Software  release-08-01-10
B2Vector3.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 
12 #include <framework/logging/Logger.h>
13 
14 #include <TVector3.h>
15 #include <Math/Vector3D.h>
16 #include <string>
17 #include <iostream> // std::cout, std::fixed
18 #include <iomanip> // std::setprecision
19 #include <typeinfo>
20 #include <cmath>
21 
22 
23 namespace Belle2 {
41  template<typename DataType>
42  class B2Vector3 {
43  protected:
45  static_assert(std::is_floating_point<DataType>::value, "B2Vector3 only works with floating point types");
47  DataType m_coordinates[3];
48  public:
50  typedef DataType value_type;
51 
53  B2Vector3(void) : m_coordinates {static_cast<DataType>(0), static_cast<DataType>(0), static_cast<DataType>(0)} {};
55  B2Vector3(const DataType xVal, const DataType yVal, const DataType zVal): m_coordinates {xVal, yVal, zVal} {};
57  explicit B2Vector3(const DataType(& coords)[3]): m_coordinates {coords[0], coords[1], coords[2]} {};
59  explicit B2Vector3(const DataType(* coords)[3]): m_coordinates {(*coords)[0], (*coords)[1], (*coords)[2]} {};
61  // cppcheck-suppress noExplicitConstructor
62  B2Vector3(const TVector3& tVec3): m_coordinates {static_cast<DataType>(tVec3.X()), static_cast<DataType>(tVec3.Y()), static_cast<DataType>(tVec3.Z())} {};
64  // cppcheck-suppress noExplicitConstructor
65  B2Vector3(const TVector3* tVec3): m_coordinates {static_cast<DataType>(tVec3->X()), static_cast<DataType>(tVec3->Y()), static_cast<DataType>(tVec3->Z())} {};
67  explicit B2Vector3(const B2Vector3<DataType>& b2Vec3): m_coordinates {b2Vec3.X(), b2Vec3.Y(), b2Vec3.Z()} {};
69  explicit B2Vector3(const B2Vector3<DataType>* b2Vec3): m_coordinates {b2Vec3->X(), b2Vec3->Y(), b2Vec3->Z()} {};
71  // cppcheck-suppress noExplicitConstructor
72  template <typename OtherType> B2Vector3(const B2Vector3<OtherType>& b2Vec3):
73  m_coordinates {static_cast<DataType>(b2Vec3.X()), static_cast<DataType>(b2Vec3.Y()), static_cast<DataType>(b2Vec3.Z())} {};
75  template <typename OtherType> explicit B2Vector3(const B2Vector3<OtherType>* b2Vec3):
76  m_coordinates {static_cast<DataType>(b2Vec3->X()), static_cast<DataType>(b2Vec3->Y()), static_cast<DataType>(b2Vec3->Z())} {};
78  // cppcheck-suppress noExplicitConstructor
79  B2Vector3(const ROOT::Math::XYZVector& xyzVec): m_coordinates {static_cast<DataType>(xyzVec.X()), static_cast<DataType>(xyzVec.Y()), static_cast<DataType>(xyzVec.Z())} {};
81  // cppcheck-suppress noExplicitConstructor
82  B2Vector3(const ROOT::Math::XYZVector* xyzVec): m_coordinates {static_cast<DataType>(xyzVec->X()), static_cast<DataType>(xyzVec->Y()), static_cast<DataType>(xyzVec->Z())} {};
83 
85  DataType operator()(unsigned i) const { return m_coordinates[i]; }
87  DataType operator[](unsigned i) const { return m_coordinates[i]; }
89  DataType& operator()(unsigned i) { return m_coordinates[i]; }
91  DataType& operator[](unsigned i) { return m_coordinates[i]; }
92 
96  B2Vector3<DataType>& operator = (const TVector3& b);
98  B2Vector3<DataType>& operator = (const ROOT::Math::XYZVector& b);
99 
101  operator TVector3() const { return GetTVector3(); }
103  operator ROOT::Math::XYZVector() const { return GetXYZVector(); }
104 
106  bool operator == (const B2Vector3<DataType>& b) const { return X() == b.X() && Y() == b.Y() && Z() == b.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(); }
112  bool operator != (const B2Vector3<DataType>& b) const { return !(*this == b); }
114  bool operator != (const TVector3& b) const { return !(*this == b); }
116  bool operator != (const ROOT::Math::XYZVector& b) const { return !(*this == b); }
117 
125  B2Vector3<DataType> operator - () const { return B2Vector3<DataType>(-X(), -Y(), -Z()); }
128  {
129  return B2Vector3<DataType>(X() + b.X(), Y() + b.Y(), Z() + b.Z());
130  }
133  {
134  return B2Vector3<DataType>(X() - b.X(), Y() - b.Y(), Z() - b.Z());
135  }
138  {
139  return B2Vector3<DataType>(a * X(), a * Y(), a * Z());
140  }
143  {
144  return B2Vector3<DataType>(X() / a, Y() / a, Z() / a);
145  }
147  DataType operator * (const B2Vector3<DataType>& b) const { return Dot(b); }
148 
149 
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()); }
160 
162  void SetPhi(DataType phi)
163  {
164  const double perp = Perp();
165  SetX(perp * cos((double)phi));
166  SetY(perp * sin((double)phi));
167  }
168 
170  void SetTheta(DataType theta)
171  {
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));
178  SetZ(ma * ctheta);
179  }
180 
182  void SetMag(DataType mag)
183  {
184  double factor = Mag();
185  if (factor == 0) {
186  B2WARNING(name() << "::SetMag: zero vector can't be stretched");
187  } else {
188  factor = mag / factor;
189  SetX(X()*factor);
190  SetY(Y()*factor);
191  SetZ(Z()*factor);
192  }
193  }
194 
196  DataType Perp2() const { return X() * X() + Y() * Y(); }
198  DataType Pt() const { return Perp(); }
200  DataType Perp() const { return std::hypot((double)X(), (double)Y()); }
201 
203  void SetPerp(DataType r)
204  {
205  const double p = Perp();
206  if (p != 0.0) {
207  m_coordinates[0] *= r / p;
208  m_coordinates[1] *= r / p;
209  }
210  }
211 
213  DataType Perp2(const B2Vector3<DataType>& axis) const
214  {
215  const double tot = axis.Mag2();
216  const double ss = Dot(axis);
217  double per = Mag2();
218  if (tot > 0.0) per -= ss * ss / tot;
219  if (per < 0) per = 0;
220  return per;
221  }
222 
224  DataType Pt(const B2Vector3<DataType>& axis) const { return Perp(axis); }
226  DataType Perp(const B2Vector3<DataType>& axis) const { return std::sqrt(Perp2(axis)); }
228  DataType DeltaPhi(const B2Vector3<DataType>& v) const { return Mpi_pi(Phi() - v.Phi()); }
229 
230 
232  static DataType Mpi_pi(DataType angle)
233  {
234  if (std::isnan(angle)) {
235  B2ERROR(name() << "::Mpi_pi: function called with NaN");
236  return angle;
237  }
238  angle = std::remainder(angle, 2 * M_PI);
239  //for compatibility with ROOT we flip the sign for exactly pi
240  if (angle == M_PI) angle = -M_PI;
241  return angle;
242  }
243 
245  DataType DeltaR(const B2Vector3<DataType>& v) const
246  {
247  const double deta = Eta() - v.Eta();
248  const double dphi = DeltaPhi(v);
249  return std::hypot(deta, dphi);
250  }
251 
253  DataType DrEtaPhi(const B2Vector3<DataType>& v) const
254  {
255  return DeltaR(v);
256  }
257 
259  void SetMagThetaPhi(DataType mag, DataType theta, DataType phi)
260  {
261  const double amag = std::abs(mag);
262  const double sinTheta = std::sin((double)theta);
263  m_coordinates[0] = amag * sinTheta * std::cos((double)phi);
264  m_coordinates[1] = amag * sinTheta * std::sin((double)phi);
265  m_coordinates[2] = amag * std::cos((double)theta);
266  }
267 
270  {
271  const double tot = Mag2();
272  B2Vector3<DataType> p(X(), Y(), Z());
273  return tot > 0.0 ? p *= (1.0 / std::sqrt(tot)) : p;
274  }
275 
278  {
279  const double xVal = std::abs((double)X());
280  const double yVal = std::abs((double)Y());
281  const double zVal = std::abs((double)Z());
282  if (xVal < yVal) {
283  return xVal < zVal ? B2Vector3<DataType>(0, Z(), -Y()) : B2Vector3<DataType>(Y(), -X(), 0);
284  } else {
285  return yVal < zVal ? B2Vector3<DataType>(-Z(), 0, X()) : B2Vector3<DataType>(Y(), -X(), 0);
286  }
287  }
288 
290  DataType Dot(const B2Vector3<DataType>& p) const
291  {
292  return X() * p.X() + Y() * p.Y() + Z() * p.Z();
293  }
294 
297  {
298  return B2Vector3<DataType>(Y() * p.Z() - p.Y() * Z(), Z() * p.X() - p.Z() * X(), X() * p.Y() - p.X() * Y());
299  }
300 
302  DataType Angle(const B2Vector3<DataType>& q) const
303  {
304  const double ptot2 = Mag2() * q.Mag2();
305  if (ptot2 <= 0) {
306  return 0.0;
307  } else {
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);
312  }
313  }
314 
319  DataType PseudoRapidity() const
320  {
321  const double cosTheta = CosTheta();
322  if (std::abs(cosTheta) < 1) return -0.5 * std::log((1.0 - cosTheta) / (1.0 + cosTheta));
323  if (Z() == 0) return 0;
324  //B2WARNING(name() << "::PseudoRapidity: transverse momentum = 0! return +/- 10e10");
325  if (Z() > 0) return 10e10;
326  else return -10e10;
327  }
328 
329 
331  DataType Eta() const { return PseudoRapidity(); }
332 
333 
335  void RotateX(DataType angle)
336  {
337  //rotate vector around X
338  const double s = std::sin((double)angle);
339  const double c = std::cos((double)angle);
340  const double yOld = Y();
341  m_coordinates[1] = c * yOld - s * Z();
342  m_coordinates[2] = s * yOld + c * Z();
343  }
344 
345 
347  void RotateY(DataType angle)
348  {
349  //rotate vector around Y
350  const double s = std::sin((double)angle);
351  const double c = std::cos((double)angle);
352  const double zOld = Z();
353  m_coordinates[0] = s * zOld + c * X();
354  m_coordinates[2] = c * zOld - s * X();
355  }
356 
357 
359  void RotateZ(DataType angle)
360  {
361  //rotate vector around Z
362  const double s = std::sin((double)angle);
363  const double c = std::cos((double)angle);
364  const double xOld = X();
365  m_coordinates[0] = c * xOld - s * Y();
366  m_coordinates[1] = s * xOld + c * Y();
367  }
368 
370  void RotateUz(const B2Vector3<DataType>& NewUzVector)
371  {
372  // NewUzVector must be normalized !
373 
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;
378 
379  if (up) {
380  up = std::sqrt(up);
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;
384  m_coordinates[2] = (u3 * u3 * px - px + u3 * up * pz) / up;
385  } else if (u3 < 0.) {
386  m_coordinates[0] = -m_coordinates[0];
387  m_coordinates[2] = -m_coordinates[2];
388  }
389  }
390 
399  void Rotate(DataType alpha, const B2Vector3<DataType>& v)
400  {
401  B2Vector3<DataType> n = v.Unit();
402  *this = (n * (n.Dot(*this)) + cos(alpha) * ((n.Cross(*this)).Cross(n)) + sin(alpha) * (n.Cross(*this)));
403  }
404 
406  void Abs()
407  {
408  m_coordinates[0] = std::abs(m_coordinates[0]);
409  m_coordinates[1] = std::abs(m_coordinates[1]);
410  m_coordinates[2] = std::abs(m_coordinates[2]);
411  }
412 
414  void Sqrt()
415  {
416  Abs();
420  }
421 
423  DataType at(unsigned i) const;
425  DataType x() const { return m_coordinates[0]; }
427  DataType y() const { return m_coordinates[1]; }
429  DataType z() const { return m_coordinates[2]; }
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(); }
442 
444  void GetXYZ(Double_t* carray) const;
446  void GetXYZ(Float_t* carray) const;
448  void GetXYZ(TVector3* tVec) const;
450  void GetXYZ(ROOT::Math::XYZVector* xyzVec) const;
452  TVector3 GetTVector3() const;
454  ROOT::Math::XYZVector GetXYZVector() const;
455 
457  void SetX(DataType x) { m_coordinates[0] = x; }
459  void SetY(DataType y) { m_coordinates[1] = y; }
461  void SetZ(DataType z) { m_coordinates[2] = z; }
462 
464  void SetXYZ(DataType x, DataType y, DataType z)
465  {
466  SetX(x); SetY(y); SetZ(z);
467  }
469  void SetXYZ(const TVector3& tVec);
471  void SetXYZ(const TVector3* tVec);
473  void SetXYZ(const ROOT::Math::XYZVector& xyzVec);
475  void SetXYZ(const ROOT::Math::XYZVector* xyzVec);
476 
478  static std::string name();
479 
481  std::string PrintString(unsigned precision = 4) const
482  {
483  return name() + " " + PrintStringXYZ(precision) + " " + PrintStringCyl(precision);
484  }
485 
487  std::string PrintStringXYZ(unsigned precision = 4) const
488  {
489  std::ostringstream output;
490  output << "(x,y,z)=("
491  << std::fixed << std::setprecision(precision)
492  << X() << "," << Y() << "," << Z() << ")";
493  return output.str();
494  }
495 
497  std::string PrintStringCyl(unsigned precision = 4) const
498  {
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 << ")";
503  return output.str();
504  }
505 
507  void Print()
508  {
509  //print vector parameters
510  Print(PrintString().c_str());
511  }
512 
513  };
514 
517 
520 
522  template <typename DataType>
523  Bool_t operator == (const TVector3& a, const B2Vector3<DataType>& b)
524  {
525  return (a.X() == b.X() && a.Y() == b.Y() && a.Z() == b.Z());
526  }
527 
529  template < typename DataType>
530  Bool_t operator != (const TVector3& a, const B2Vector3<DataType>& b)
531  {
532  return !(a == b);
533  }
534 
536  template < typename DataType>
538  {
539  return B2Vector3<DataType>(a * p.X(), a * p.Y(), a * p.Z());
540  }
541 
543  template < typename DataType>
545  {
546  return B2Vector3<DataType>(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
547  }
548 
550  template < typename DataType>
552  {
553  return B2Vector3<DataType>(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
554  }
555 
557  template < typename DataType>
559  {
560  return B2Vector3<DataType>(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
561  }
562 
564  template < typename DataType>
566  {
567  return B2Vector3<DataType>(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
568  }
569 
571  template < typename DataType>
572  B2Vector3<DataType> operator + (const ROOT::Math::XYZVector& a, const B2Vector3<DataType>& b)
573  {
574  return B2Vector3<DataType>(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
575  }
576 
578  template < typename DataType>
579  B2Vector3<DataType> operator - (const ROOT::Math::XYZVector& a, const B2Vector3<DataType>& b)
580  {
581  return B2Vector3<DataType>(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
582  }
583 
585  template < typename DataType>
586  B2Vector3<DataType> operator + (const B2Vector3<DataType>& a, const ROOT::Math::XYZVector& b)
587  {
588  return B2Vector3<DataType>(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
589  }
590 
592  template < typename DataType>
593  B2Vector3<DataType> operator - (const B2Vector3<DataType>& a, const ROOT::Math::XYZVector& b)
594  {
595  return B2Vector3<DataType>(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
596  }
597 
598 
600  template< typename DataType >
602  {
603  m_coordinates[0] = b.X();
604  m_coordinates[1] = b.Y();
605  m_coordinates[2] = b.Z();
606  return *this;
607  }
608 
610  template< typename DataType >
612  {
613  m_coordinates[0] = b.X();
614  m_coordinates[1] = b.Y();
615  m_coordinates[2] = b.Z();
616  return *this;
617  }
618 
620  template< typename DataType >
621  B2Vector3<DataType>& B2Vector3<DataType>::operator = (const ROOT::Math::XYZVector& b)
622  {
623  m_coordinates[0] = b.X();
624  m_coordinates[1] = b.Y();
625  m_coordinates[2] = b.Z();
626  return *this;
627  }
628 
630  template< typename DataType >
632  {
633  m_coordinates[0] += b.X();
634  m_coordinates[1] += b.Y();
635  m_coordinates[2] += b.Z();
636  return *this;
637  }
638 
639 
641  template< typename DataType >
643  {
644  m_coordinates[0] -= b.X();
645  m_coordinates[1] -= b.Y();
646  m_coordinates[2] -= b.Z();
647  return *this;
648  }
649 
651  template< typename DataType >
653  {
654  m_coordinates[0] *= a;
655  m_coordinates[1] *= a;
656  m_coordinates[2] *= a;
657  return *this;
658  }
659 
661  template< typename DataType >
662  void B2Vector3<DataType>::SetXYZ(const TVector3& tVec)
663  {
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());
667  }
668 
670  template< typename DataType >
671  void B2Vector3<DataType>::SetXYZ(const TVector3* tVec)
672  {
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());
676  }
677 
679  template< typename DataType >
680  void B2Vector3<DataType>::SetXYZ(const ROOT::Math::XYZVector& xyzVec)
681  {
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());
685  }
686 
688  template< typename DataType >
689  void B2Vector3<DataType>::SetXYZ(const ROOT::Math::XYZVector* xyzVec)
690  {
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());
694  }
695 
696  template< typename DataType >
697  void B2Vector3<DataType>::GetXYZ(double* carray) const
698  {
699  carray[0] = X();
700  carray[1] = Y();
701  carray[2] = Z();
702  }
703 
705  template< typename DataType >
706  void B2Vector3<DataType>::GetXYZ(TVector3* tVec) const
707  {
708  tVec->SetXYZ(static_cast<Double_t>(X()),
709  static_cast<Double_t>(Y()),
710  static_cast<Double_t>(Z()));
711  }
712 
714  template< typename DataType >
715  void B2Vector3<DataType>::GetXYZ(ROOT::Math::XYZVector* xyzVec) const
716  {
717  xyzVec->SetXYZ(static_cast<Double_t>(X()),
718  static_cast<Double_t>(Y()),
719  static_cast<Double_t>(Z()));
720  }
721 
722 
724  template< typename DataType >
726  {
727  return
728  TVector3(
729  static_cast<Double_t>(X()),
730  static_cast<Double_t>(Y()),
731  static_cast<Double_t>(Z())
732  );
733  }
734 
735 
737  template< typename DataType >
738  ROOT::Math::XYZVector B2Vector3<DataType>::GetXYZVector() const
739  {
740  return
741  ROOT::Math::XYZVector(
742  static_cast<Double_t>(X()),
743  static_cast<Double_t>(Y()),
744  static_cast<Double_t>(Z())
745  );
746  }
747 
748 
750  template < typename DataType>
751  DataType B2Vector3<DataType>::at(unsigned i) const
752  {
753  switch (i) {
754  case 0:
756  case 1:
758  case 2:
760  }
761  B2FATAL(this->name() << "::access operator: given index (i=" << i << ") is out of bounds!");
762  return 0.;
763  }
764 
766  template < typename DataType>
768  {
769  return std::string("B2Vector3<") + typeid(DataType).name() + std::string(">");
770  }
771 
773 } // end namespace Belle2
A fast and root compatible alternative to TVector3.
Definition: B2Vector3.h:42
DataType Phi() const
The azimuth angle.
Definition: B2Vector3.h:151
B2Vector3< DataType > operator+(const B2Vector3< DataType > &b) const
Addition of 3-vectors.
Definition: B2Vector3.h:127
B2Vector3< DataType > Cross(const B2Vector3< DataType > &p) const
Cross product.
Definition: B2Vector3.h:296
void SetMag(DataType mag)
Set magnitude keeping theta and phi constant.
Definition: B2Vector3.h:182
DataType Pz() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:441
DataType m_coordinates[3]
Make sure that we only have floating point vectors.
Definition: B2Vector3.h:45
B2Vector3(const B2Vector3< DataType > &b2Vec3)
Constructor expecting a B2Vector3 of same type.
Definition: B2Vector3.h:67
B2Vector3< DataType > operator/(DataType a) const
Scaling of 3-vectors with a real number.
Definition: B2Vector3.h:142
DataType & operator()(unsigned i)
member access without boundary check
Definition: B2Vector3.h:89
DataType operator()(unsigned i) const
member access without boundary check
Definition: B2Vector3.h:85
B2Vector3< DataType > operator*(DataType a) const
Scaling of 3-vectors with a real number.
Definition: B2Vector3.h:137
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:435
DataType Perp2() const
The transverse component squared (R^2 in cylindrical coordinate system).
Definition: B2Vector3.h:196
void SetPerp(DataType r)
Set the transverse component keeping phi and z constant.
Definition: B2Vector3.h:203
DataType y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:427
void SetX(DataType x)
set X/1st-coordinate
Definition: B2Vector3.h:457
DataType Theta() const
The polar angle.
Definition: B2Vector3.h:153
DataType CosTheta() const
Cosine of the polar angle.
Definition: B2Vector3.h:155
DataType & operator[](unsigned i)
member access without boundary check
Definition: B2Vector3.h:91
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)
Definition: B2Vector3.h:429
B2Vector3(const B2Vector3< DataType > *b2Vec3)
Constructor expecting a pointer to a B2Vector3.
Definition: B2Vector3.h:69
void SetMagThetaPhi(DataType mag, DataType theta, DataType phi)
setter with mag, theta, phi
Definition: B2Vector3.h:259
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:431
std::string PrintStringXYZ(unsigned precision=4) const
create a string containing vector in cartesian coordinates
Definition: B2Vector3.h:487
DataType Eta() const
Returns the pseudo-rapidity.
Definition: B2Vector3.h:331
DataType DeltaPhi(const B2Vector3< DataType > &v) const
returns phi in the interval [-PI,PI)
Definition: B2Vector3.h:228
void RotateY(DataType angle)
Rotates the B2Vector3 around the y-axis.
Definition: B2Vector3.h:347
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:433
DataType operator[](unsigned i) const
member access without boundary check
Definition: B2Vector3.h:87
std::string PrintString(unsigned precision=4) const
create a string containing vector in cartesian and spherical coordinates
Definition: B2Vector3.h:481
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:159
DataType value_type
storage type of the vector
Definition: B2Vector3.h:50
B2Vector3(const ROOT::Math::XYZVector &xyzVec)
Constructor expecting a XYZVector.
Definition: B2Vector3.h:79
bool operator==(const B2Vector3< DataType > &b) const
Comparison for equality with a B2Vector3.
Definition: B2Vector3.h:106
DataType x() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:425
B2Vector3(const DataType(*coords)[3])
Constructor using a pointer.
Definition: B2Vector3.h:59
B2Vector3< DataType > operator-() const
unary minus
Definition: B2Vector3.h:125
DataType DeltaR(const B2Vector3< DataType > &v) const
return deltaR with respect to input-vector
Definition: B2Vector3.h:245
DataType Perp(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
Definition: B2Vector3.h:226
void RotateX(DataType angle)
Rotates the B2Vector3 around the x-axis.
Definition: B2Vector3.h:335
void Abs()
calculates the absolute value of the coordinates element-wise
Definition: B2Vector3.h:406
B2Vector3(const DataType(&coords)[3])
Constructor using a reference.
Definition: B2Vector3.h:57
B2Vector3(const ROOT::Math::XYZVector *xyzVec)
Constructor expecting a pointer to a XYZVector.
Definition: B2Vector3.h:82
DataType Pt(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
Definition: B2Vector3.h:224
DataType Mag2() const
The magnitude squared (rho^2 in spherical coordinate system).
Definition: B2Vector3.h:157
void SetTheta(DataType theta)
Set theta keeping mag and phi constant.
Definition: B2Vector3.h:170
B2Vector3(const TVector3 &tVec3)
Constructor expecting a TVector3.
Definition: B2Vector3.h:62
void RotateZ(DataType angle)
Rotates the B2Vector3 around the z-axis.
Definition: B2Vector3.h:359
void Print()
just for backward compatibility, should not be used with new code
Definition: B2Vector3.h:507
void Sqrt()
calculates the square root of the absolute values of the coordinates element-wise
Definition: B2Vector3.h:414
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
Definition: B2Vector3.h:269
void SetZ(DataType z)
set Z/3rd-coordinate
Definition: B2Vector3.h:461
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
Definition: B2Vector3.h:290
B2Vector3< DataType > Orthogonal() const
Vector orthogonal to this one.
Definition: B2Vector3.h:277
void SetY(DataType y)
set Y/2nd-coordinate
Definition: B2Vector3.h:459
B2Vector3(const B2Vector3< OtherType > *b2Vec3)
Constructor expecting a pointer to a B2Vector3 of different type.
Definition: B2Vector3.h:75
DataType PseudoRapidity() const
Returns the pseudo-rapidity, i.e.
Definition: B2Vector3.h:319
DataType DrEtaPhi(const B2Vector3< DataType > &v) const
return DrEtaPhi with respect to input-vector
Definition: B2Vector3.h:253
bool operator!=(const B2Vector3< DataType > &b) const
Comparison != with a B2Vector3.
Definition: B2Vector3.h:112
DataType Perp() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:200
B2Vector3(void)
empty Constructor sets everything to 0
Definition: B2Vector3.h:53
DataType Perp2(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
Definition: B2Vector3.h:213
B2Vector3(const B2Vector3< OtherType > &b2Vec3)
Constructor expecting a B2Vector3 of different type.
Definition: B2Vector3.h:72
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).
Definition: B2Vector3.h:370
DataType Py() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:439
static DataType Mpi_pi(DataType angle)
returns given angle in the interval [-PI,PI)
Definition: B2Vector3.h:232
DataType Pt() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:198
B2Vector3(const TVector3 *tVec3)
Constructor expecting a pointer to a TVector3.
Definition: B2Vector3.h:65
DataType Px() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:437
void SetXYZ(DataType x, DataType y, DataType z)
set all coordinates using data type
Definition: B2Vector3.h:464
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
Definition: B2Vector3.h:302
std::string PrintStringCyl(unsigned precision=4) const
create a string containing vector in spherical coordinates
Definition: B2Vector3.h:497
void Rotate(DataType alpha, const B2Vector3< DataType > &v)
Rotation around an arbitrary axis v with angle alpha.
Definition: B2Vector3.h:399
void SetPhi(DataType phi)
Set phi keeping mag and theta constant.
Definition: B2Vector3.h:162
B2Vector3(const DataType xVal, const DataType yVal, const DataType zVal)
Constructor expecting 3 coordinates.
Definition: B2Vector3.h:55
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 ...
Definition: DecayNode.cc:48
bool operator!=(const DecayNode &node1, const DecayNode &node2)
Not equal: See operator==.
Definition: DecayNode.cc:65
void SetXYZ(const ROOT::Math::XYZVector *xyzVec)
set all coordinates using a pointer to XYZVector
Definition: B2Vector3.h:689
void GetXYZ(ROOT::Math::XYZVector *xyzVec) const
directly copies coordinates to a XYZVector
Definition: B2Vector3.h:715
void SetXYZ(const TVector3 &tVec)
set all coordinates using a reference to TVector3
Definition: B2Vector3.h:662
B2Vector3< DataType > & operator-=(const B2Vector3< DataType > &b)
subtraction
Definition: B2Vector3.h:642
TVector3 GetTVector3() const
returns a TVector3 containing the same coordinates
Definition: B2Vector3.h:725
B2Vector3< float > B2Vector3F
typedef for common usage with float
Definition: B2Vector3.h:519
void GetXYZ(TVector3 *tVec) const
directly copies coordinates to a TVector3
Definition: B2Vector3.h:706
void SetXYZ(const TVector3 *tVec)
set all coordinates using a pointer to TVector3
Definition: B2Vector3.h:671
ROOT::Math::XYZVector GetXYZVector() const
returns a XYZVector containing the same coordinates
Definition: B2Vector3.h:738
B2Vector3< DataType > & operator+=(const B2Vector3< DataType > &b)
addition
Definition: B2Vector3.h:631
B2Vector3< DataType > operator*(DataType a, const B2Vector3< DataType > &p)
non-memberfunction Scaling of 3-vectors with a real number
Definition: B2Vector3.h:537
B2Vector3< DataType > & operator*=(DataType a)
scaling with real numbers
Definition: B2Vector3.h:652
void SetXYZ(const ROOT::Math::XYZVector &xyzVec)
set all coordinates using a reference to XYZVector
Definition: B2Vector3.h:680
B2Vector3< DataType > & operator=(const B2Vector3< DataType > &b)
Assignment via B2Vector3.
Definition: B2Vector3.h:601
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
DataType at(unsigned i) const
safe member access (with boundary check!)
Definition: B2Vector3.h:751
static std::string name()
Returns the name of the B2Vector.
Definition: B2Vector3.h:767
B2Vector3< DataType > operator-(const TVector3 &a, const B2Vector3< DataType > &b)
non-memberfunction for substracting a TVector3 from a B2Vector3
Definition: B2Vector3.h:551
B2Vector3< DataType > operator+(const TVector3 &a, const B2Vector3< DataType > &b)
non-memberfunction for adding a TVector3 to a B2Vector3
Definition: B2Vector3.h:544
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.