Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:46:07

0001 #ifndef FASTSIM_BARRELSIMPLIFIEDGEOMETRY_H
0002 #define FASTSIM_BARRELSIMPLIFIEDGEOMETRY_H
0003 
0004 #include "FastSimulation/SimplifiedGeometryPropagator/interface/SimplifiedGeometry.h"
0005 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0006 #include "DataFormats/Math/interface/Vector3D.h"
0007 #include "FWCore/Utilities/interface/Exception.h"
0008 #include <TH1F.h>
0009 
0010 ///////////////////////////////////////////////
0011 // Author: L. Vanelderen, S. Kurz
0012 // Date: 29 May 2017
0013 //////////////////////////////////////////////////////////
0014 
0015 namespace fastsim {
0016 
0017   //! Implementation of a barrel detector layer (cylindrical).
0018   /*!
0019         A cylindrical layer with a given radius and a thickness (in radiation length).
0020         The layer is regarded infinitely long in Z-direction, however the thickness can vary (as a function of Z) and also be 0.
0021     */
0022   class BarrelSimplifiedGeometry : public SimplifiedGeometry {
0023   public:
0024     //! Constructor.
0025     /*!
0026             Create a barrel layer with a given radius.
0027             \param radius The radius of the layer (in cm).
0028         */
0029     BarrelSimplifiedGeometry(double radius) : SimplifiedGeometry(radius) {}
0030 
0031     //! Move constructor.
0032     BarrelSimplifiedGeometry(BarrelSimplifiedGeometry &&) = default;
0033 
0034     //! Default destructor.
0035     ~BarrelSimplifiedGeometry() override{};
0036 
0037     //! Return radius of the barrel layer.
0038     /*!
0039             \return The radius of the layer (in cm).
0040         */
0041     const double getRadius() const { return geomProperty_; }
0042 
0043     //! Return thickness of the barrel layer at a given position.
0044     /*!
0045             Returns the thickness of the barrel layer (in radiation length) at a specified position since the thickness can vary as a function of Z.
0046             \param position A position which has to be on the barrel layer.
0047             \return The thickness of the layer (in radiation length).
0048             \sa getThickness(const math::XYZTLorentzVector & position, const math::XYZTLorentzVector & momentum)
0049         */
0050     const double getThickness(const math::XYZTLorentzVector &position) const override {
0051       return thicknessHist_->GetBinContent(thicknessHist_->GetXaxis()->FindBin(fabs(position.Z())));
0052     }
0053 
0054     //! Return thickness of the barrel layer at a given position, also considering the incident angle.
0055     /*!
0056             Returns the thickness of the barrel layer (in radiation length) at a specified position and a given incident angle since the thickness can vary as a function of Z.
0057             \param position A position which has to be on the barrel layer.
0058             \param momentum The momentum of the incident particle.
0059             \return The thickness of the layer (in radiation length).
0060             \sa getThickness(const math::XYZTLorentzVector & position)
0061         */
0062     const double getThickness(const math::XYZTLorentzVector &position,
0063                               const math::XYZTLorentzVector &momentum) const override {
0064       // Do projection of norm(layer) on momentum vector
0065       // CosTheta = (momentum dot norm) / (length(momentum) / length(norm))
0066       return getThickness(position) /
0067              (fabs(momentum.X() * position.X() + momentum.Y() * position.Y()) /
0068               (momentum.P() * std::sqrt(position.X() * position.X() + position.Y() * position.Y())));
0069     }
0070 
0071     //! Return magnetic field (field only has Z component!) on the barrel layer.
0072     /*!
0073             Returns the magnetic field along the barrel layer at a specified position (radial symmetric).
0074             \param position A position which has to be on the barrel layer.
0075             \return The magnetic field on the layer.
0076         */
0077     const double getMagneticFieldZ(const math::XYZTLorentzVector &position) const override {
0078       return magneticFieldHist_->GetBinContent(magneticFieldHist_->GetXaxis()->FindBin(fabs(position.z())));
0079     }
0080 
0081     //! Returns false since class for barrel layer.
0082     /*!
0083             Function to easily destinguish barrel from forward layers (which both inherit from SimplifiedGeometry).
0084             \return false
0085         */
0086     bool isForward() const override { return false; }
0087   };
0088 
0089 }  // namespace fastsim
0090 
0091 #endif