Belle II Software  release-08-01-10
PCmsLabTransform.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 #include <mdst/dbobjects/CollisionBoostVector.h>
12 #include <mdst/dbobjects/CollisionInvariantMass.h>
13 #include <mdst/dbobjects/CollisionAxisCMS.h>
14 #include <framework/database/DBObjPtr.h>
15 #include <framework/geometry/B2Vector3.h>
16 
17 #include <Math/Boost.h>
18 #include <Math/AxisAngle.h>
19 #include <Math/LorentzRotation.h>
20 #include <Math/Vector4D.h>
21 
22 namespace Belle2 {
32 
33  public:
34 
39 
44  {
45  return m_boostVectorDB->getBoost();
46  }
47 
51  double getCMSEnergy() const
52  {
53  return m_invariantMassDB->getMass();
54  }
55 
59  ROOT::Math::PxPyPzEVector getBeamFourMomentum() const
60  {
61  return rotateCmsToLab() * ROOT::Math::PxPyPzEVector(0, 0, 0, getCMSEnergy());
62  }
63 
69  const ROOT::Math::LorentzRotation rotateLabToCms() const
70  {
71  //boost to CM frame
72  ROOT::Math::LorentzRotation boost(ROOT::Math::Boost(-1.*getBoostVector()));
73 
74 
75  //rotation such that the collision axis is aligned with z-axis
76  ROOT::Math::XYZVector zaxis(0., 0., 1.); //target collision axis
77 
78  double tanAngleXZ = tan(m_axisCmsDB->getAngleXZ());
79  double tanAngleYZ = tan(m_axisCmsDB->getAngleYZ());
80  double Norm = 1 / sqrt(1 + pow(tanAngleXZ, 2) + pow(tanAngleYZ, 2));
81  ROOT::Math::XYZVector electronCMS(Norm * tanAngleXZ, Norm * tanAngleYZ, Norm); //current collision axis
82 
83  ROOT::Math::XYZVector rotAxis = zaxis.Cross(electronCMS);
84  double rotangle = asin(rotAxis.R());
85 
86  ROOT::Math::LorentzRotation rotation(ROOT::Math::AxisAngle(rotAxis, -rotangle));
87 
88 
89  ROOT::Math::LorentzRotation trans = rotation * boost;
90  return trans;
91  }
92 
97  const ROOT::Math::LorentzRotation rotateCmsToLab() const
98  {
99  return rotateLabToCms().Inverse();
100  }
101 
107  static ROOT::Math::PxPyPzMVector labToCms(const ROOT::Math::PxPyPzMVector& vec);
108 
114  static ROOT::Math::PxPyPzMVector cmsToLab(const ROOT::Math::PxPyPzMVector& vec);
115 
121  static ROOT::Math::PxPyPzEVector labToCms(const ROOT::Math::PxPyPzEVector& vec);
122 
128  static ROOT::Math::PxPyPzEVector cmsToLab(const ROOT::Math::PxPyPzEVector& vec);
129 
130  private:
134  };
135 
137 } // Belle2 namespace
138 
139 
140 
141 
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Class to hold Lorentz transformations from/to CMS and boost vector.
ROOT::Math::PxPyPzEVector getBeamFourMomentum() const
Returns LAB four-momentum of e+e-, i.e.
double getCMSEnergy() const
Returns CMS energy of e+e- (aka.
static ROOT::Math::PxPyPzMVector labToCms(const ROOT::Math::PxPyPzMVector &vec)
Transforms Lorentz vector into CM System.
const ROOT::Math::LorentzRotation rotateLabToCms() const
Returns Lorentz transformation from Lab to CMS.
PCmsLabTransform()
Constructor.
static ROOT::Math::PxPyPzMVector cmsToLab(const ROOT::Math::PxPyPzMVector &vec)
Transforms Lorentz vector into Laboratory System.
const DBObjPtr< CollisionBoostVector > m_boostVectorDB
db object for boost vector.
const ROOT::Math::LorentzRotation rotateCmsToLab() const
Returns Lorentz transformation from CMS to Lab.
const DBObjPtr< CollisionInvariantMass > m_invariantMassDB
db object for invariant mass.
B2Vector3D getBoostVector() const
Returns boost vector (beta=p/E)
const DBObjPtr< CollisionAxisCMS > m_axisCmsDB
db object for collision axis in CM system from boost.
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
double tan(double a)
tan for double
Definition: beamHelpers.h:31
Abstract base class for different kinds of events.