CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
10 // Author: L. Vanelderen, S. Kurz
11 // Date: 29 May 2017
13 
14 namespace fastsim {
15 
17 
22  public:
24 
29 
32 
35 
37 
40  const double getZ() const { return geomProperty_; }
41 
43 
49  const double getThickness(const math::XYZTLorentzVector &position) const override {
50  return thicknessHist_->GetBinContent(thicknessHist_->GetXaxis()->FindBin(position.Pt()));
51  }
52 
54 
62  const math::XYZTLorentzVector &momentum) const override {
63  return getThickness(position) / fabs(momentum.Pz()) * momentum.P();
64  }
65 
67 
72  const double getMagneticFieldZ(const math::XYZTLorentzVector &position) const override {
73  return magneticFieldHist_->GetBinContent(magneticFieldHist_->GetXaxis()->FindBin(position.Pt()));
74  }
75 
77 
81  bool isForward() const override { return true; }
82  };
83 
84 } // namespace fastsim
85 
86 #endif
Implementation of a forward detector layer (disk).
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.
~ForwardSimplifiedGeometry() override
Default destructor.
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 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 getZ() const
Return z-position of the forward layer.
const double getMagneticFieldZ(const math::XYZTLorentzVector &position) const override
Return magnetic field (field only has Z component!) on the forward layer.
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:289
bool isForward() const override
Returns true since class for forward layer.
double geomProperty_
Geometric property of the layer: radius (barrel layer) / position z (forward layer) ...