CMS 3D CMS Logo

BarrelSimplifiedGeometry.h
Go to the documentation of this file.
1 #ifndef FASTSIM_BARRELSIMPLIFIEDGEOMETRY_H
2 #define FASTSIM_BARRELSIMPLIFIEDGEOMETRY_H
3 
8 #include <TH1F.h>
9 
10 
12 // Author: L. Vanelderen, S. Kurz
13 // Date: 29 May 2017
15 
16 
17 namespace fastsim{
18 
20 
25  {
26  public:
28 
33  SimplifiedGeometry(radius) {}
34 
37 
40 
42 
45  const double getRadius() const { return geomProperty_; }
46 
48 
54  const double getThickness(const math::XYZTLorentzVector & position) const override
55  {
56  return thicknessHist_->GetBinContent(thicknessHist_->GetXaxis()->FindBin(fabs(position.Z())));
57  }
58 
60 
67  const double getThickness(const math::XYZTLorentzVector & position, const math::XYZTLorentzVector & momentum) const override
68  {
69  // Do projection of norm(layer) on momentum vector
70  // CosTheta = (momentum dot norm) / (length(momentum) / length(norm))
71  return getThickness(position) / (fabs(momentum.X() * position.X() + momentum.Y() * position.Y()) / (momentum.P() * std::sqrt(position.X() * position.X() + position.Y() * position.Y())));
72  }
73 
75 
80  const double getMagneticFieldZ(const math::XYZTLorentzVector & position) const override
81  {
82  return magneticFieldHist_->GetBinContent(magneticFieldHist_->GetXaxis()->FindBin(fabs(position.z())));
83  }
84 
86 
90  bool isForward() const override
91  {
92  return false;
93  }
94  };
95 
96 }
97 
98 #endif
const double getMagneticFieldZ(const math::XYZTLorentzVector &position) const override
Return magnetic field (field only has Z component!) on the barrel layer.
Implementation of a generic detector layer (base class for forward/barrel layers).
std::unique_ptr< TH1F > magneticFieldHist_
Histogram that stores the size of the magnetic field along the layer.
Implementation of a barrel detector layer (cylindrical).
BarrelSimplifiedGeometry(double radius)
Constructor.
std::unique_ptr< TH1F > thicknessHist_
Histogram that stores the tickness (radLengths) along the layer.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T sqrt(T t)
Definition: SSEVec.h:18
const double getThickness(const math::XYZTLorentzVector &position) const override
Return thickness of the barrel layer at a given position.
static int position[264][3]
Definition: ReadPGInfo.cc:509
bool isForward() const override
Returns false since class for barrel layer.
const double getRadius() const
Return radius of the barrel layer.
double geomProperty_
Geometric property of the layer: radius (barrel layer) / position z (forward layer) ...
const double getThickness(const math::XYZTLorentzVector &position, const math::XYZTLorentzVector &momentum) const override
Return thickness of the barrel layer at a given position, also considering the incident angle...