Belle II Software  release-06-01-15
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 <string>
16 #include <iostream> // std::cout, std::fixed
17 #include <iomanip> // std::setprecision
18 #include <typeinfo>
19 #include <cmath>
20 
21 
22 namespace Belle2 {
40  template<typename DataType>
41  class B2Vector3 {
42  protected:
44  static_assert(std::is_floating_point<DataType>::value, "B2Vector3 only works with floating point types");
46  DataType m_coordinates[3];
47  public:
49  typedef DataType value_type;
50 
52  B2Vector3(void) : m_coordinates {static_cast<DataType>(0), static_cast<DataType>(0), static_cast<DataType>(0)} {};
54  B2Vector3(const DataType xVal, const DataType yVal, const DataType zVal): m_coordinates {xVal, yVal, zVal} {};
56  explicit B2Vector3(const DataType(& coords)[3]): m_coordinates {coords[0], coords[1], coords[2]} {};
58  explicit B2Vector3(const DataType(* coords)[3]): m_coordinates {(*coords)[0], (*coords)[1], (*coords)[2]} {};
60  // cppcheck-suppress noExplicitConstructor
61  B2Vector3(const TVector3& tVec3): m_coordinates {static_cast<DataType>(tVec3.X()), static_cast<DataType>(tVec3.Y()), static_cast<DataType>(tVec3.Z())} {};
63  // cppcheck-suppress noExplicitConstructor
64  B2Vector3(const TVector3* tVec3): m_coordinates {static_cast<DataType>(tVec3->X()), static_cast<DataType>(tVec3->Y()), static_cast<DataType>(tVec3->Z())} {};
66  explicit B2Vector3(const B2Vector3<DataType>& b2Vec3): m_coordinates {b2Vec3.X(), b2Vec3.Y(), b2Vec3.Z()} {};
68  explicit B2Vector3(const B2Vector3<DataType>* b2Vec3): m_coordinates {b2Vec3->X(), b2Vec3->Y(), b2Vec3->Z()} {};
70  // cppcheck-suppress noExplicitConstructor
71  template <typename OtherType> B2Vector3(const B2Vector3<OtherType>& b2Vec3):
72  m_coordinates {static_cast<DataType>(b2Vec3.X()), static_cast<DataType>(b2Vec3.Y()), static_cast<DataType>(b2Vec3.Z())} {};
74  template <typename OtherType> explicit B2Vector3(const B2Vector3<OtherType>* b2Vec3):
75  m_coordinates {static_cast<DataType>(b2Vec3->X()), static_cast<DataType>(b2Vec3->Y()), static_cast<DataType>(b2Vec3->Z())} {};
76 
78  DataType operator()(unsigned i) const { return m_coordinates[i]; }
80  DataType operator[](unsigned i) const { return m_coordinates[i]; }
82  DataType& operator()(unsigned i) { return m_coordinates[i]; }
84  DataType& operator[](unsigned i) { return m_coordinates[i]; }
85 
89  B2Vector3<DataType>& operator = (const TVector3& b);
90 
92  operator TVector3() const { return GetTVector3(); }
93 
95  bool operator == (const B2Vector3<DataType>& b) const { return X() == b.X() && Y() == b.Y() && Z() == b.Z(); }
97  bool operator == (const TVector3& b) const { return X() == b.X() && Y() == b.Y() && Z() == b.Z(); }
99  bool operator != (const B2Vector3<DataType>& b) const { return !(*this == b); }
101  bool operator != (const TVector3& b) const { return !(*this == b); }
102 
110  B2Vector3<DataType> operator - () const { return B2Vector3<DataType>(-X(), -Y(), -Z()); }
113  {
114  return B2Vector3<DataType>(X() + b.X(), Y() + b.Y(), Z() + b.Z());
115  }
118  {
119  return B2Vector3<DataType>(X() - b.X(), Y() - b.Y(), Z() - b.Z());
120  }
123  {
124  return B2Vector3<DataType>(a * X(), a * Y(), a * Z());
125  }
128  {
129  return B2Vector3<DataType>(X() / a, Y() / a, Z() / a);
130  }
132  DataType operator * (const B2Vector3<DataType>& b) const { return Dot(b); }
133 
134 
136  DataType Phi() const { return X() == 0 && Y() == 0 ? 0 : atan2(Y(), X()); }
138  DataType Theta() const { return X() == 0 && Y() == 0 && Z() == 0 ? 0 : atan2(Perp(), Z()); }
140  DataType CosTheta() const { const double pTot = Mag(); return pTot == 0 ? 1 : Z() / pTot; }
142  DataType Mag2() const { return X() * X() + Y() * Y() + Z() * Z(); }
144  DataType Mag() const { return std::hypot((double)Perp(), (double)Z()); }
145 
147  void SetPhi(DataType phi)
148  {
149  const double perp = Perp();
150  SetX(perp * cos((double)phi));
151  SetY(perp * sin((double)phi));
152  }
153 
155  void SetTheta(DataType theta)
156  {
157  const double ma = Mag();
158  const double ph = Phi();
159  const double ctheta = std::cos((double) theta);
160  const double stheta = std::sin((double) theta);
161  SetX(ma * stheta * std::cos(ph));
162  SetY(ma * stheta * std::cos(ph));
163  SetZ(ma * ctheta);
164  }
165 
167  void SetMag(DataType mag)
168  {
169  double factor = Mag();
170  if (factor == 0) {
171  B2WARNING(name() << "::SetMag: zero vector can't be stretched");
172  } else {
173  factor = mag / factor;
174  SetX(X()*factor);
175  SetY(Y()*factor);
176  SetZ(Z()*factor);
177  }
178  }
179 
181  DataType Perp2() const { return X() * X() + Y() * Y(); }
183  DataType Pt() const { return Perp(); }
185  DataType Perp() const { return std::hypot((double)X(), (double)Y()); }
186 
188  void SetPerp(DataType r)
189  {
190  const double p = Perp();
191  if (p != 0.0) {
192  m_coordinates[0] *= r / p;
193  m_coordinates[1] *= r / p;
194  }
195  }
196 
198  DataType Perp2(const B2Vector3<DataType>& axis) const
199  {
200  const double tot = axis.Mag2();
201  const double ss = Dot(axis);
202  double per = Mag2();
203  if (tot > 0.0) per -= ss * ss / tot;
204  if (per < 0) per = 0;
205  return per;
206  }
207 
209  DataType Pt(const B2Vector3<DataType>& axis) const { return Perp(axis); }
211  DataType Perp(const B2Vector3<DataType>& axis) const { return std::sqrt(Perp2(axis)); }
213  DataType DeltaPhi(const B2Vector3<DataType>& v) const { return Mpi_pi(Phi() - v.Phi()); }
214 
215 
217  static DataType Mpi_pi(DataType angle)
218  {
219  if (std::isnan(angle)) {
220  B2ERROR(name() << "::Mpi_pi: function called with NaN");
221  return angle;
222  }
223  angle = std::remainder(angle, 2 * M_PI);
224  //for compatibility with ROOT we flip the sign for exactly pi
225  if (angle == M_PI) angle = -M_PI;
226  return angle;
227  }
228 
230  DataType DeltaR(const B2Vector3<DataType>& v) const
231  {
232  const double deta = Eta() - v.Eta();
233  const double dphi = DeltaPhi(v);
234  return std::hypot(deta, dphi);
235  }
236 
238  DataType DrEtaPhi(const B2Vector3<DataType>& v) const
239  {
240  return DeltaR(v);
241  }
242 
244  void SetMagThetaPhi(DataType mag, DataType theta, DataType phi)
245  {
246  const double amag = std::abs(mag);
247  const double sinTheta = std::sin((double)theta);
248  m_coordinates[0] = amag * sinTheta * std::cos((double)phi);
249  m_coordinates[1] = amag * sinTheta * std::sin((double)phi);
250  m_coordinates[2] = amag * std::cos((double)theta);
251  }
252 
255  {
256  const double tot = Mag2();
257  B2Vector3<DataType> p(X(), Y(), Z());
258  return tot > 0.0 ? p *= (1.0 / std::sqrt(tot)) : p;
259  }
260 
263  {
264  const double xVal = std::abs((double)X());
265  const double yVal = std::abs((double)Y());
266  const double zVal = std::abs((double)Z());
267  if (xVal < yVal) {
268  return xVal < zVal ? B2Vector3<DataType>(0, Z(), -Y()) : B2Vector3<DataType>(Y(), -X(), 0);
269  } else {
270  return yVal < zVal ? B2Vector3<DataType>(-Z(), 0, X()) : B2Vector3<DataType>(Y(), -X(), 0);
271  }
272  }
273 
275  DataType Dot(const B2Vector3<DataType>& p) const
276  {
277  return X() * p.X() + Y() * p.Y() + Z() * p.Z();
278  }
279 
282  {
283  return B2Vector3<DataType>(Y() * p.Z() - p.Y() * Z(), Z() * p.X() - p.Z() * X(), X() * p.Y() - p.X() * Y());
284  }
285 
287  DataType Angle(const B2Vector3<DataType>& q) const
288  {
289  const double ptot2 = Mag2() * q.Mag2();
290  if (ptot2 <= 0) {
291  return 0.0;
292  } else {
293  double arg = Dot(q) / std::sqrt(ptot2);
294  if (arg > 1.0) arg = 1.0;
295  if (arg < -1.0) arg = -1.0;
296  return std::acos(arg);
297  }
298  }
299 
304  DataType PseudoRapidity() const
305  {
306  const double cosTheta = CosTheta();
307  if (std::abs(cosTheta) < 1) return -0.5 * std::log((1.0 - cosTheta) / (1.0 + cosTheta));
308  if (Z() == 0) return 0;
309  //B2WARNING(name() << "::PseudoRapidity: transverse momentum = 0! return +/- 10e10");
310  if (Z() > 0) return 10e10;
311  else return -10e10;
312  }
313 
314 
316  DataType Eta() const { return PseudoRapidity(); }
317 
318 
320  void RotateX(DataType angle)
321  {
322  //rotate vector around X
323  const double s = std::sin((double)angle);
324  const double c = std::cos((double)angle);
325  const double yOld = Y();
326  m_coordinates[1] = c * yOld - s * Z();
327  m_coordinates[2] = s * yOld + c * Z();
328  }
329 
330 
332  void RotateY(DataType angle)
333  {
334  //rotate vector around Y
335  const double s = std::sin((double)angle);
336  const double c = std::cos((double)angle);
337  const double zOld = Z();
338  m_coordinates[0] = s * zOld + c * X();
339  m_coordinates[2] = c * zOld - s * X();
340  }
341 
342 
344  void RotateZ(DataType angle)
345  {
346  //rotate vector around Z
347  const double s = std::sin((double)angle);
348  const double c = std::cos((double)angle);
349  const double xOld = X();
350  m_coordinates[0] = c * xOld - s * Y();
351  m_coordinates[1] = s * xOld + c * Y();
352  }
353 
355  void RotateUz(const B2Vector3<DataType>& NewUzVector)
356  {
357  // NewUzVector must be normalized !
358 
359  const double u1 = NewUzVector.X();
360  const double u2 = NewUzVector.Y();
361  const double u3 = NewUzVector.Z();
362  double up = u1 * u1 + u2 * u2;
363 
364  if (up) {
365  up = std::sqrt(up);
366  DataType px = X(), py = Y(), pz = Z();
367  m_coordinates[0] = (u1 * u3 * px - u2 * py + u1 * up * pz) / up;
368  m_coordinates[1] = (u2 * u3 * px + u1 * py + u2 * up * pz) / up;
369  m_coordinates[2] = (u3 * u3 * px - px + u3 * up * pz) / up;
370  } else if (u3 < 0.) {
371  m_coordinates[0] = -m_coordinates[0];
372  m_coordinates[2] = -m_coordinates[2];
373  }
374  }
375 
384  void Rotate(DataType alpha, const B2Vector3<DataType>& v)
385  {
386  B2Vector3<DataType> n = v.Unit();
387  *this = (n * (n.Dot(*this)) + cos(alpha) * ((n.Cross(*this)).Cross(n)) + sin(alpha) * (n.Cross(*this)));
388  }
389 
391  void Abs()
392  {
393  m_coordinates[0] = std::abs(m_coordinates[0]);
394  m_coordinates[1] = std::abs(m_coordinates[1]);
395  m_coordinates[2] = std::abs(m_coordinates[2]);
396  }
397 
399  void Sqrt()
400  {
401  Abs();
402  m_coordinates[0] = std::sqrt(m_coordinates[0]);
403  m_coordinates[1] = std::sqrt(m_coordinates[1]);
404  m_coordinates[2] = std::sqrt(m_coordinates[2]);
405  }
406 
408  DataType at(unsigned i) const;
410  DataType x() const { return m_coordinates[0]; }
412  DataType y() const { return m_coordinates[1]; }
414  DataType z() const { return m_coordinates[2]; }
416  DataType X() const { return x(); }
418  DataType Y() const { return y(); }
420  DataType Z() const { return z(); }
422  DataType Px() const { return x(); }
424  DataType Py() const { return y(); }
426  DataType Pz() const { return z(); }
427 
429  void GetXYZ(Double_t* carray) const;
431  void GetXYZ(Float_t* carray) const;
433  void GetXYZ(TVector3* tVec) const;
435  TVector3 GetTVector3() const;
436 
438  void SetX(DataType x) { m_coordinates[0] = x; }
440  void SetY(DataType y) { m_coordinates[1] = y; }
442  void SetZ(DataType z) { m_coordinates[2] = z; }
443 
445  void SetXYZ(DataType x, DataType y, DataType z)
446  {
447  SetX(x); SetY(y); SetZ(z);
448  }
450  void SetXYZ(const TVector3& tVec);
452  void SetXYZ(const TVector3* tVec);
453 
455  static std::string name();
456 
458  std::string PrintString(unsigned precision = 4) const
459  {
460  return name() + " " + PrintStringXYZ(precision) + " " + PrintStringCyl(precision);
461  }
462 
464  std::string PrintStringXYZ(unsigned precision = 4) const
465  {
466  std::ostringstream output;
467  output << "(x,y,z)=("
468  << std::fixed << std::setprecision(precision)
469  << X() << "," << Y() << "," << Z() << ")";
470  return output.str();
471  }
472 
474  std::string PrintStringCyl(unsigned precision = 4) const
475  {
476  std::ostringstream output;
477  output << "(rho, theta, phi)=("
478  << std::fixed << std::setprecision(precision)
479  << Mag() << "," << Theta() * 180. / M_PI << "," << Phi() * 180. / M_PI << ")";
480  return output.str();
481  }
482 
484  void Print()
485  {
486  //print vector parameters
487  Print(PrintString().c_str());
488  }
489 
490  };
491 
494 
497 
499  template <typename DataType>
500  Bool_t operator == (const TVector3& a, const B2Vector3<DataType>& b)
501  {
502  return (a.X() == b.X() && a.Y() == b.Y() && a.Z() == b.Z());
503  }
504 
506  template < typename DataType>
507  Bool_t operator != (const TVector3& a, const B2Vector3<DataType>& b)
508  {
509  return !(a == b);
510  }
511 
513  template < typename DataType>
515  {
516  return B2Vector3<DataType>(a * p.X(), a * p.Y(), a * p.Z());
517  }
518 
520  template < typename DataType>
522  {
523  return B2Vector3<DataType>(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
524  }
525 
527  template < typename DataType>
529  {
530  return B2Vector3<DataType>(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
531  }
532 
534  template < typename DataType>
536  {
537  return B2Vector3<DataType>(a.X() + b.X(), a.Y() + b.Y(), a.Z() + b.Z());
538  }
539 
541  template < typename DataType>
543  {
544  return B2Vector3<DataType>(a.X() - b.X(), a.Y() - b.Y(), a.Z() - b.Z());
545  }
546 
547 
549  template< typename DataType >
551  {
552  m_coordinates[0] = b.X();
553  m_coordinates[1] = b.Y();
554  m_coordinates[2] = b.Z();
555  return *this;
556  }
557 
559  template< typename DataType >
561  {
562  m_coordinates[0] = b.X();
563  m_coordinates[1] = b.Y();
564  m_coordinates[2] = b.Z();
565  return *this;
566  }
567 
569  template< typename DataType >
571  {
572  m_coordinates[0] += b.X();
573  m_coordinates[1] += b.Y();
574  m_coordinates[2] += b.Z();
575  return *this;
576  }
577 
578 
580  template< typename DataType >
582  {
583  m_coordinates[0] -= b.X();
584  m_coordinates[1] -= b.Y();
585  m_coordinates[2] -= b.Z();
586  return *this;
587  }
588 
590  template< typename DataType >
592  {
593  m_coordinates[0] *= a;
594  m_coordinates[1] *= a;
595  m_coordinates[2] *= a;
596  return *this;
597  }
598 
600  template< typename DataType >
601  void B2Vector3<DataType>::SetXYZ(const TVector3& tVec)
602  {
603  m_coordinates[0] = static_cast<Double_t>(tVec.X());
604  m_coordinates[1] = static_cast<Double_t>(tVec.Y());
605  m_coordinates[2] = static_cast<Double_t>(tVec.Z());
606  }
607 
609  template< typename DataType >
610  void B2Vector3<DataType>::SetXYZ(const TVector3* tVec)
611  {
612  m_coordinates[0] = static_cast<Double_t>(tVec->X());
613  m_coordinates[1] = static_cast<Double_t>(tVec->Y());
614  m_coordinates[2] = static_cast<Double_t>(tVec->Z());
615  }
616 
617  template< typename DataType >
618  void B2Vector3<DataType>::GetXYZ(double* carray) const
619  {
620  carray[0] = X();
621  carray[1] = Y();
622  carray[2] = Z();
623  }
624 
626  template< typename DataType >
627  void B2Vector3<DataType>::GetXYZ(TVector3* tVec) const
628  {
629  tVec->SetXYZ(static_cast<Double_t>(X()),
630  static_cast<Double_t>(Y()),
631  static_cast<Double_t>(Z()));
632  }
633 
634 
636  template< typename DataType >
638  {
639  return
640  TVector3(
641  static_cast<Double_t>(X()),
642  static_cast<Double_t>(Y()),
643  static_cast<Double_t>(Z())
644  );
645  }
646 
647 
649  template < typename DataType>
650  DataType B2Vector3<DataType>::at(unsigned i) const
651  {
652  switch (i) {
653  case 0:
655  case 1:
657  case 2:
659  }
660  B2FATAL(this->name() << "::access operator: given index (i=" << i << ") is out of bounds!");
661  return 0.;
662  }
663 
665  template < typename DataType>
667  {
668  return std::string("B2Vector3<") + typeid(DataType).name() + std::string(">");
669  }
670 
672 } // end namespace Belle2
A fast and root compatible alternative to TVector3.
Definition: B2Vector3.h:41
DataType Phi() const
The azimuth angle.
Definition: B2Vector3.h:136
B2Vector3< DataType > operator+(const B2Vector3< DataType > &b) const
Addition of 3-vectors.
Definition: B2Vector3.h:112
B2Vector3< DataType > Cross(const B2Vector3< DataType > &p) const
Cross product.
Definition: B2Vector3.h:281
void SetMag(DataType mag)
Set magnitude keeping theta and phi constant.
Definition: B2Vector3.h:167
DataType Pz() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:426
DataType m_coordinates[3]
Make sure that we only have floating point vectors.
Definition: B2Vector3.h:44
B2Vector3(const B2Vector3< DataType > &b2Vec3)
Constructor expecting a B2Vector3 of same type.
Definition: B2Vector3.h:66
B2Vector3< DataType > operator/(DataType a) const
Scaling of 3-vectors with a real number.
Definition: B2Vector3.h:127
DataType & operator()(unsigned i)
member access without boundary check
Definition: B2Vector3.h:82
DataType operator()(unsigned i) const
member access without boundary check
Definition: B2Vector3.h:78
B2Vector3< DataType > operator*(DataType a) const
Scaling of 3-vectors with a real number.
Definition: B2Vector3.h:122
DataType Z() const
access variable Z (= .at(2) without boundary check)
Definition: B2Vector3.h:420
DataType Perp2() const
The transverse component squared (R^2 in cylindrical coordinate system).
Definition: B2Vector3.h:181
void SetPerp(DataType r)
Set the transverse component keeping phi and z constant.
Definition: B2Vector3.h:188
DataType y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:412
void SetX(DataType x)
set X/1st-coordinate
Definition: B2Vector3.h:438
DataType Theta() const
The polar angle.
Definition: B2Vector3.h:138
DataType CosTheta() const
Cosine of the polar angle.
Definition: B2Vector3.h:140
DataType & operator[](unsigned i)
member access without boundary check
Definition: B2Vector3.h:84
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:414
B2Vector3(const B2Vector3< DataType > *b2Vec3)
Constructor expecting a pointer to a B2Vector3.
Definition: B2Vector3.h:68
void SetMagThetaPhi(DataType mag, DataType theta, DataType phi)
setter with mag, theta, phi
Definition: B2Vector3.h:244
DataType X() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:416
std::string PrintStringXYZ(unsigned precision=4) const
create a string containing vector in cartesian coordinates
Definition: B2Vector3.h:464
DataType Eta() const
Returns the pseudo-rapidity.
Definition: B2Vector3.h:316
DataType DeltaPhi(const B2Vector3< DataType > &v) const
returns phi in the interval [-PI,PI)
Definition: B2Vector3.h:213
void RotateY(DataType angle)
Rotates the B2Vector3 around the y-axis.
Definition: B2Vector3.h:332
DataType Y() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:418
DataType operator[](unsigned i) const
member access without boundary check
Definition: B2Vector3.h:80
std::string PrintString(unsigned precision=4) const
create a string containing vector in cartesian and spherical coordinates
Definition: B2Vector3.h:458
DataType Mag() const
The magnitude (rho in spherical coordinate system).
Definition: B2Vector3.h:144
DataType value_type
storage type of the vector
Definition: B2Vector3.h:49
bool operator==(const B2Vector3< DataType > &b) const
Comparison for equality with a B2Vector3.
Definition: B2Vector3.h:95
DataType x() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:410
B2Vector3(const DataType(*coords)[3])
Constructor using a pointer.
Definition: B2Vector3.h:58
B2Vector3< DataType > operator-() const
unary minus
Definition: B2Vector3.h:110
DataType DeltaR(const B2Vector3< DataType > &v) const
return deltaR with respect to input-vector
Definition: B2Vector3.h:230
DataType Perp(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
Definition: B2Vector3.h:211
void RotateX(DataType angle)
Rotates the B2Vector3 around the x-axis.
Definition: B2Vector3.h:320
void Abs()
calculates the absolute value of the coordinates element-wise
Definition: B2Vector3.h:391
B2Vector3(const DataType(&coords)[3])
Constructor using a reference.
Definition: B2Vector3.h:56
DataType Pt(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
Definition: B2Vector3.h:209
DataType Mag2() const
The magnitude squared (rho^2 in spherical coordinate system).
Definition: B2Vector3.h:142
void SetTheta(DataType theta)
Set theta keeping mag and phi constant.
Definition: B2Vector3.h:155
B2Vector3(const TVector3 &tVec3)
Constructor expecting a TVector3.
Definition: B2Vector3.h:61
void RotateZ(DataType angle)
Rotates the B2Vector3 around the z-axis.
Definition: B2Vector3.h:344
void Print()
just for backward compatibility, should not be used with new code
Definition: B2Vector3.h:484
void Sqrt()
calculates the square root of the absolute values of the coordinates element-wise
Definition: B2Vector3.h:399
B2Vector3< DataType > Unit() const
Unit vector parallel to this.
Definition: B2Vector3.h:254
void SetZ(DataType z)
set Z/3rd-coordinate
Definition: B2Vector3.h:442
DataType Dot(const B2Vector3< DataType > &p) const
Scalar product.
Definition: B2Vector3.h:275
B2Vector3< DataType > Orthogonal() const
Vector orthogonal to this one.
Definition: B2Vector3.h:262
void SetY(DataType y)
set Y/2nd-coordinate
Definition: B2Vector3.h:440
B2Vector3(const B2Vector3< OtherType > *b2Vec3)
Constructor expecting a pointer to a B2Vector3 of different type.
Definition: B2Vector3.h:74
DataType PseudoRapidity() const
Returns the pseudo-rapidity, i.e.
Definition: B2Vector3.h:304
DataType DrEtaPhi(const B2Vector3< DataType > &v) const
return DrEtaPhi with respect to input-vector
Definition: B2Vector3.h:238
bool operator!=(const B2Vector3< DataType > &b) const
Comparison != with a B2Vector3.
Definition: B2Vector3.h:99
DataType Perp() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:185
B2Vector3(void)
empty Constructor sets everything to 0
Definition: B2Vector3.h:52
DataType Perp2(const B2Vector3< DataType > &axis) const
The transverse component w.r.t.
Definition: B2Vector3.h:198
B2Vector3(const B2Vector3< OtherType > &b2Vec3)
Constructor expecting a B2Vector3 of different type.
Definition: B2Vector3.h:71
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:355
DataType Py() const
access variable Y (= .at(1) without boundary check)
Definition: B2Vector3.h:424
static DataType Mpi_pi(DataType angle)
returns given angle in the interval [-PI,PI)
Definition: B2Vector3.h:217
DataType Pt() const
The transverse component (R in cylindrical coordinate system).
Definition: B2Vector3.h:183
B2Vector3(const TVector3 *tVec3)
Constructor expecting a pointer to a TVector3.
Definition: B2Vector3.h:64
DataType Px() const
access variable X (= .at(0) without boundary check)
Definition: B2Vector3.h:422
void SetXYZ(DataType x, DataType y, DataType z)
set all coordinates using data type
Definition: B2Vector3.h:445
DataType Angle(const B2Vector3< DataType > &q) const
The angle w.r.t.
Definition: B2Vector3.h:287
std::string PrintStringCyl(unsigned precision=4) const
create a string containing vector in spherical coordinates
Definition: B2Vector3.h:474
void Rotate(DataType alpha, const B2Vector3< DataType > &v)
Rotation around an arbitrary axis v with angle alpha.
Definition: B2Vector3.h:384
void SetPhi(DataType phi)
Set phi keeping mag and theta constant.
Definition: B2Vector3.h:147
B2Vector3(const DataType xVal, const DataType yVal, const DataType zVal)
Constructor expecting 3 coordinates.
Definition: B2Vector3.h:54
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 TVector3 &tVec)
set all coordinates using a reference to TVector3
Definition: B2Vector3.h:601
B2Vector3< DataType > & operator-=(const B2Vector3< DataType > &b)
subtraction
Definition: B2Vector3.h:581
TVector3 GetTVector3() const
returns a TVector3 containing the same coordinates
Definition: B2Vector3.h:637
B2Vector3< float > B2Vector3F
typedef for common usage with float
Definition: B2Vector3.h:496
void GetXYZ(TVector3 *tVec) const
directly copies coordinates to a TVector3
Definition: B2Vector3.h:627
void SetXYZ(const TVector3 *tVec)
set all coordinates using a pointer to TVector3
Definition: B2Vector3.h:610
B2Vector3< DataType > & operator+=(const B2Vector3< DataType > &b)
addition
Definition: B2Vector3.h:570
B2Vector3< DataType > operator*(DataType a, const B2Vector3< DataType > &p)
non-memberfunction Scaling of 3-vectors with a real number
Definition: B2Vector3.h:514
B2Vector3< DataType > & operator*=(DataType a)
scaling with real numbers
Definition: B2Vector3.h:591
B2Vector3< DataType > & operator=(const B2Vector3< DataType > &b)
Assignment via B2Vector3.
Definition: B2Vector3.h:550
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:493
DataType at(unsigned i) const
safe member access (with boundary check!)
Definition: B2Vector3.h:650
static std::string name()
Returns the name of the B2Vector.
Definition: B2Vector3.h:666
B2Vector3< DataType > operator-(const TVector3 &a, const B2Vector3< DataType > &b)
non-memberfunction for substracting a TVector3 from a B2Vector3
Definition: B2Vector3.h:528
B2Vector3< DataType > operator+(const TVector3 &a, const B2Vector3< DataType > &b)
non-memberfunction for adding a TVector3 to a B2Vector3
Definition: B2Vector3.h:521
Abstract base class for different kinds of events.