|
||||
File indexing completed on 2024-04-06 12:11:20
0001 #ifndef FASTSIM_SIMPLIFIEDGEOMETRY_H 0002 #define FASTSIM_SIMPLIFIEDGEOMETRY_H 0003 0004 #include "DataFormats/Math/interface/LorentzVector.h" 0005 0006 #include <memory> 0007 #include <vector> 0008 0009 /////////////////////////////////////////////// 0010 // Author: L. Vanelderen, S. Kurz 0011 // Date: 29 May 2017 0012 ////////////////////////////////////////////////////////// 0013 0014 class DetLayer; 0015 class MagneticField; 0016 class GeometricSearchTracker; 0017 class TH1F; 0018 0019 namespace edm { 0020 class ParameterSet; 0021 } 0022 0023 namespace fastsim { 0024 class InteractionModel; 0025 class SimplifiedGeometryFactory; 0026 0027 //! Implementation of a generic detector layer (base class for forward/barrel layers). 0028 /*! 0029 Each layer has a geometric attribute ('geomProperty') which depends on which kind of layer is actually created 0030 (radius for barrel layers, position z for forward layers). Furthermore, a thickness in radiation length is assigned 0031 which can vary throughout the layer. 0032 \sa BarrelSimplifiedGeometry 0033 \sa ForwardSimplifiedGeometry 0034 */ 0035 class SimplifiedGeometry { 0036 public: 0037 //! Default constructor. 0038 /*! 0039 'geomProperty' depends on which kind of layer is actually created 0040 - BarrelSimplifiedGeometry: radius 0041 - ForwardSimplifiedGeometry: z 0042 */ 0043 SimplifiedGeometry(double geomProperty); 0044 0045 SimplifiedGeometry(SimplifiedGeometry&&) = default; 0046 SimplifiedGeometry(SimplifiedGeometry const&) = delete; 0047 0048 SimplifiedGeometry& operator=(SimplifiedGeometry&&) = default; 0049 SimplifiedGeometry& operator=(SimplifiedGeometry const&) = delete; 0050 0051 //! Default destructor. 0052 virtual ~SimplifiedGeometry(); 0053 0054 //////// 0055 // HACK 0056 //////// 0057 0058 //! Hack to interface "old" Calorimetry with "new" Tracker 0059 enum CaloType { NONE, TRACKERBOUNDARY, PRESHOWER1, PRESHOWER2, ECAL, HCAL, VFCAL }; 0060 0061 //! Hack to interface "old" Calorimetry with "new" Tracker 0062 void setCaloType(CaloType type) { caloType_ = type; } 0063 0064 //! Hack to interface "old" Calorimetry with "new" Tracker 0065 CaloType getCaloType() const { return caloType_; } 0066 0067 //////// 0068 // END 0069 //////// 0070 0071 //! Set index of this layer (layers are numbered according to their position in the detector). 0072 /*! 0073 The (usual) order is increasing 'geomProperty': 0074 - BarrelLayers: innermost to outermost 0075 - ForwardLayers: -z to +z 0076 Makes it more convenient to navigate from one layer to the previous/next layer. 0077 */ 0078 void setIndex(int index) { index_ = index; } 0079 0080 //! Return index of this layer (layers are numbered according to their position in the detector). 0081 /*! 0082 The (usual) order is increasing 'geomProperty': 0083 - BarrelLayers: innermost to outermost 0084 - ForwardLayers: -z to +z 0085 Makes it more convenient to navigate from one layer to the previous/next layer. 0086 */ 0087 int index() const { return index_; } 0088 0089 //! Return geometric attribute of the layer. 0090 /*! 0091 'geomProperty' depends on which kind of layer is actually created 0092 - BarrelSimplifiedGeometry: radius 0093 - ForwardSimplifiedGeometry: z 0094 */ 0095 const double getGeomProperty() const { return geomProperty_; } 0096 0097 //! Return thickness of the layer at a given position. 0098 /*! 0099 Returns the thickness of the layer (in radiation length) at a specified position since the thickness can vary throughout the layer. 0100 \param position A position which has to be on the layer. 0101 \return The thickness of the layer (in radiation length). 0102 \sa getThickness(const math::XYZTLorentzVector & position, const math::XYZTLorentzVector & momentum) 0103 */ 0104 virtual const double getThickness(const math::XYZTLorentzVector& position) const = 0; 0105 0106 //! Return thickness of the layer at a given position, also considering the incident angle. 0107 /*! 0108 Returns the thickness of the layer (in radiation length) at a specified position and a given incident angle since the thickness can vary throughout the layer. 0109 \param position A position which has to be on the layer. 0110 \param momentum The momentum of the incident particle. 0111 \return The thickness of the layer (in radiation length). 0112 \sa getThickness(const math::XYZTLorentzVector & position) 0113 */ 0114 virtual const double getThickness(const math::XYZTLorentzVector& position, 0115 const math::XYZTLorentzVector& momentum) const = 0; 0116 0117 //! Some layers have a different thickness for nuclear interactions. 0118 /*! 0119 Right now only considered for TEC layers. 0120 */ 0121 const double getNuclearInteractionThicknessFactor() const { return nuclearInteractionThicknessFactor_; } 0122 0123 //! Return pointer to the assigned active layer (if any). 0124 /*! 0125 Necessary to create SimHits. 0126 */ 0127 const DetLayer* getDetLayer() const { return detLayer_; } 0128 0129 //! Return magnetic field (field only has Z component!) on the layer. 0130 /*! 0131 \param position A position which has to be on the layer. 0132 \return The magnetic field on the layer. 0133 */ 0134 virtual const double getMagneticFieldZ(const math::XYZTLorentzVector& position) const = 0; 0135 0136 //! Returns false/true depending if the object is a (non-abstract) barrel/forward layer. 0137 /*! 0138 Function to easily destinguish barrel from forward layers (which both inherit from SimplifiedGeometry). 0139 */ 0140 virtual bool isForward() const = 0; 0141 0142 //! Return the vector of all interaction models that are assigned with a layer. 0143 /*! 0144 This makes it easy to switch on/off some interactions for some layers. 0145 */ 0146 const std::vector<InteractionModel*>& getInteractionModels() const { return interactionModels_; } 0147 0148 //! Some basic output. 0149 friend std::ostream& operator<<(std::ostream& os, const SimplifiedGeometry& layer); 0150 0151 friend class fastsim::SimplifiedGeometryFactory; 0152 0153 protected: 0154 double geomProperty_; //!< Geometric property of the layer: radius (barrel layer) / position z (forward layer) 0155 int index_; //!< Return index of this layer (layers are numbered according to their position in the detector). The (usual) order is increasing 'geomProperty'. 0156 const DetLayer* detLayer_; //!< Return pointer to the assigned active layer (if any). 0157 std::unique_ptr<TH1F> magneticFieldHist_; //!< Histogram that stores the size of the magnetic field along the layer. 0158 std::unique_ptr<TH1F> thicknessHist_; //!< Histogram that stores the tickness (radLengths) along the layer. 0159 double nuclearInteractionThicknessFactor_; //!< Some layers have a different thickness for nuclear interactions. 0160 std::vector<InteractionModel*> 0161 interactionModels_; //!< Vector of all interaction models that are assigned with a layer. 0162 CaloType caloType_; //!< Hack to interface "old" Calorimetry with "new" Tracker 0163 }; 0164 0165 std::ostream& operator<<(std::ostream& os, const SimplifiedGeometry& layer); 0166 0167 } // namespace fastsim 0168 0169 #endif
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.2.1 LXR engine. The LXR team |