CMS 3D CMS Logo

ForwardSimplifiedGeometry.h
Go to the documentation of this file.
1 #ifndef FASTSIM_FORWARDSIMPLIFIEDGEOMETRY_H
2 #define FASTSIM_FORWARDSIMPLIFIEDGEOMETRY_H
3 
6 #include <TH1F.h>
8 
9 
11 // Author: L. Vanelderen, S. Kurz
12 // Date: 29 May 2017
14 
15 
16 namespace fastsim{
17 
19 
24  {
25  public:
27 
32  SimplifiedGeometry(z) {}
33 
36 
39 
41 
44  const double getZ() const { return geomProperty_; }
45 
47 
53  const double getThickness(const math::XYZTLorentzVector & position) const override
54  {
55  return thicknessHist_->GetBinContent(thicknessHist_->GetXaxis()->FindBin(position.Pt()));
56  }
57 
59 
66  const double getThickness(const math::XYZTLorentzVector & position, const math::XYZTLorentzVector & momentum) const override
67  {
68  return getThickness(position) / fabs(momentum.Pz()) * momentum.P();
69  }
70 
72 
77  const double getMagneticFieldZ(const math::XYZTLorentzVector & position) const override
78  {
79  return magneticFieldHist_->GetBinContent(magneticFieldHist_->GetXaxis()->FindBin(position.Pt()));
80  }
81 
83 
87  bool isForward() const override {return true;}
88  };
89 
90 }
91 
92 #endif
Implementation of a forward detector layer (disk).
Implementation of a generic detector layer (base class for forward/barrel layers).
const double getMagneticFieldZ(const math::XYZTLorentzVector &position) const override
Return magnetic field (field only has Z component!) on the forward layer.
std::unique_ptr< TH1F > magneticFieldHist_
Histogram that stores the size of the magnetic field along the layer.
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
const double getZ() const
Return z-position of the forward layer.
bool isForward() const override
Returns true since class for forward layer.
const double getThickness(const math::XYZTLorentzVector &position, const math::XYZTLorentzVector &momentum) const override
Return thickness of the forward layer at a given position, also considering the incident angle...
const double getThickness(const math::XYZTLorentzVector &position) const override
Return thickness of the forward layer at a given position.
static int position[264][3]
Definition: ReadPGInfo.cc:509
double geomProperty_
Geometric property of the layer: radius (barrel layer) / position z (forward layer) ...