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 
11 // Author: L. Vanelderen, S. Kurz
12 // Date: 29 May 2017
14 
15 namespace fastsim {
16 
18 
23  public:
25 
30 
33 
36 
38 
41  const double getRadius() const { return geomProperty_; }
42 
44 
50  const double getThickness(const math::XYZTLorentzVector &position) const override {
51  return thicknessHist_->GetBinContent(thicknessHist_->GetXaxis()->FindBin(fabs(position.Z())));
52  }
53 
55 
63  const math::XYZTLorentzVector &momentum) const override {
64  // Do projection of norm(layer) on momentum vector
65  // CosTheta = (momentum dot norm) / (length(momentum) / length(norm))
66  return getThickness(position) /
67  (fabs(momentum.X() * position.X() + momentum.Y() * position.Y()) /
68  (momentum.P() * std::sqrt(position.X() * position.X() + position.Y() * position.Y())));
69  }
70 
72 
77  const double getMagneticFieldZ(const math::XYZTLorentzVector &position) const override {
78  return magneticFieldHist_->GetBinContent(magneticFieldHist_->GetXaxis()->FindBin(fabs(position.z())));
79  }
80 
82 
86  bool isForward() const override { return false; }
87  };
88 
89 } // namespace fastsim
90 
91 #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).
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...
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:23
bool isForward() const override
Returns false since class for barrel layer.
~BarrelSimplifiedGeometry() override
Default destructor.
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:289
const double getRadius() const
Return radius of the barrel layer.
double geomProperty_
Geometric property of the layer: radius (barrel layer) / position z (forward layer) ...