CMS 3D CMS Logo

List of all members | Public Member Functions
fastsim::BarrelSimplifiedGeometry Class Reference

Implementation of a barrel detector layer (cylindrical). More...

#include <BarrelSimplifiedGeometry.h>

Inheritance diagram for fastsim::BarrelSimplifiedGeometry:
fastsim::SimplifiedGeometry

Public Member Functions

 BarrelSimplifiedGeometry (double radius)
 Constructor. More...
 
 BarrelSimplifiedGeometry (BarrelSimplifiedGeometry &&)=default
 Move constructor. More...
 
const double getMagneticFieldZ (const math::XYZTLorentzVector &position) const override
 Return magnetic field (field only has Z component!) on the barrel layer. More...
 
const double getRadius () const
 Return radius of the barrel layer. More...
 
const double getThickness (const math::XYZTLorentzVector &position) const override
 Return thickness of the barrel layer at a given position. More...
 
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. More...
 
bool isForward () const override
 Returns false since class for barrel layer. More...
 
 ~BarrelSimplifiedGeometry ()
 Default destructor. More...
 
- Public Member Functions inherited from fastsim::SimplifiedGeometry
CaloType getCaloType () const
 Hack to interface "old" Calorimetry with "new" Tracker. More...
 
const DetLayergetDetLayer () const
 Return pointer to the assigned active layer (if any). More...
 
const double getGeomProperty () const
 Return geometric attribute of the layer. More...
 
const std::vector< InteractionModel * > & getInteractionModels () const
 Return the vector of all interaction models that are assigned with a layer. More...
 
const double getNuclearInteractionThicknessFactor () const
 Some layers have a different thickness for nuclear interactions. More...
 
int index () const
 Return index of this layer (layers are numbered according to their position in the detector). More...
 
void setCaloType (CaloType type)
 Hack to interface "old" Calorimetry with "new" Tracker. More...
 
void setIndex (int index)
 Set index of this layer (layers are numbered according to their position in the detector). More...
 
 SimplifiedGeometry (double geomProperty)
 Default constructor. More...
 
 ~SimplifiedGeometry ()
 Default destructor. More...
 

Additional Inherited Members

- Public Types inherited from fastsim::SimplifiedGeometry
enum  CaloType {
  NONE, TRACKERBOUNDARY, PRESHOWER1, PRESHOWER2,
  ECAL, HCAL, VFCAL
}
 Hack to interface "old" Calorimetry with "new" Tracker. More...
 
- Protected Attributes inherited from fastsim::SimplifiedGeometry
CaloType caloType_
 Hack to interface "old" Calorimetry with "new" Tracker. More...
 
const DetLayerdetLayer_
 Return pointer to the assigned active layer (if any). More...
 
double geomProperty_
 Geometric property of the layer: radius (barrel layer) / position z (forward layer) More...
 
int index_
 Return index of this layer (layers are numbered according to their position in the detector). The (usual) order is increasing 'geomProperty'. More...
 
std::vector< InteractionModel * > interactionModels_
 Vector of all interaction models that are assigned with a layer. More...
 
std::unique_ptr< TH1F > magneticFieldHist_
 Histogram that stores the size of the magnetic field along the layer. More...
 
double nuclearInteractionThicknessFactor_
 Some layers have a different thickness for nuclear interactions. More...
 
std::unique_ptr< TH1F > thicknessHist_
 Histogram that stores the tickness (radLengths) along the layer. More...
 

Detailed Description

Implementation of a barrel detector layer (cylindrical).

A cylindrical layer with a given radius and a thickness (in radiation length). The layer is regarded infinitely long in Z-direction, however the thickness can vary (as a function of Z) and also be 0.

Definition at line 24 of file BarrelSimplifiedGeometry.h.

Constructor & Destructor Documentation

fastsim::BarrelSimplifiedGeometry::BarrelSimplifiedGeometry ( double  radius)
inline

Constructor.

Create a barrel layer with a given radius.

Parameters
radiusThe radius of the layer (in cm).

Definition at line 32 of file BarrelSimplifiedGeometry.h.

32  :
SimplifiedGeometry(double geomProperty)
Default constructor.
fastsim::BarrelSimplifiedGeometry::BarrelSimplifiedGeometry ( BarrelSimplifiedGeometry &&  )
default

Move constructor.

fastsim::BarrelSimplifiedGeometry::~BarrelSimplifiedGeometry ( )
inline

Default destructor.

Definition at line 39 of file BarrelSimplifiedGeometry.h.

39 {};

Member Function Documentation

const double fastsim::BarrelSimplifiedGeometry::getMagneticFieldZ ( const math::XYZTLorentzVector position) const
inlineoverridevirtual

Return magnetic field (field only has Z component!) on the barrel layer.

Returns the magnetic field along the barrel layer at a specified position (radial symmetric).

Parameters
positionA position which has to be on the barrel layer.
Returns
The magnetic field on the layer.

Implements fastsim::SimplifiedGeometry.

Definition at line 80 of file BarrelSimplifiedGeometry.h.

References fastsim::SimplifiedGeometry::magneticFieldHist_.

81  {
82  return magneticFieldHist_->GetBinContent(magneticFieldHist_->GetXaxis()->FindBin(fabs(position.z())));
83  }
std::unique_ptr< TH1F > magneticFieldHist_
Histogram that stores the size of the magnetic field along the layer.
static int position[264][3]
Definition: ReadPGInfo.cc:509
const double fastsim::BarrelSimplifiedGeometry::getRadius ( ) const
inline

Return radius of the barrel layer.

Returns
The radius of the layer (in cm).

Definition at line 45 of file BarrelSimplifiedGeometry.h.

References fastsim::SimplifiedGeometry::geomProperty_.

Referenced by fastsim::HelixTrajectory::crosses(), fastsim::HelixTrajectory::nextCrossingTimeC(), and fastsim::StraightTrajectory::nextCrossingTimeC().

45 { return geomProperty_; }
double geomProperty_
Geometric property of the layer: radius (barrel layer) / position z (forward layer) ...
const double fastsim::BarrelSimplifiedGeometry::getThickness ( const math::XYZTLorentzVector position) const
inlineoverridevirtual

Return thickness of the barrel layer at a given position.

Returns the thickness of the barrel layer (in radiation length) at a specified position since the thickness can vary as a function of Z.

Parameters
positionA position which has to be on the barrel layer.
Returns
The thickness of the layer (in radiation length).
See also
getThickness(const math::XYZTLorentzVector & position, const math::XYZTLorentzVector & momentum)

Implements fastsim::SimplifiedGeometry.

Definition at line 54 of file BarrelSimplifiedGeometry.h.

References fastsim::SimplifiedGeometry::thicknessHist_.

Referenced by getThickness().

55  {
56  return thicknessHist_->GetBinContent(thicknessHist_->GetXaxis()->FindBin(fabs(position.Z())));
57  }
std::unique_ptr< TH1F > thicknessHist_
Histogram that stores the tickness (radLengths) along the layer.
static int position[264][3]
Definition: ReadPGInfo.cc:509
const double fastsim::BarrelSimplifiedGeometry::getThickness ( const math::XYZTLorentzVector position,
const math::XYZTLorentzVector momentum 
) const
inlineoverridevirtual

Return thickness of the barrel layer at a given position, also considering the incident angle.

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.

Parameters
positionA position which has to be on the barrel layer.
momentumThe momentum of the incident particle.
Returns
The thickness of the layer (in radiation length).
See also
getThickness(const math::XYZTLorentzVector & position)

Implements fastsim::SimplifiedGeometry.

Definition at line 67 of file BarrelSimplifiedGeometry.h.

References getThickness(), and mathSSE::sqrt().

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  }
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 fastsim::BarrelSimplifiedGeometry::isForward ( ) const
inlineoverridevirtual

Returns false since class for barrel layer.

Function to easily destinguish barrel from forward layers (which both inherit from SimplifiedGeometry).

Returns
false

Implements fastsim::SimplifiedGeometry.

Definition at line 90 of file BarrelSimplifiedGeometry.h.

91  {
92  return false;
93  }