CMS 3D CMS Logo

Geometry.h
Go to the documentation of this file.
1 #ifndef FASTSIM_GEOMETRY_H
2 #define FASTSIM_GEOMETRY_H
3 
8 
9 
11 // Author: L. Vanelderen, S. Kurz
12 // Date: 29 May 2017
14 
15 
17 class MagneticField;
18 
19 #include <vector>
20 
21 namespace edm {
22  //class ParameterSet;
23  class EventSetup;
24 }
25 
26 namespace fastsim{
27  class InteractionModel;
28 
30 
34  class Geometry
35  {
36  public:
37 
40 
42  ~Geometry();
43 
45 
51  void update(const edm::EventSetup & iSetup, const std::map<std::string,InteractionModel*> & interactionModelMap);
52 
54 
59  double getMagneticFieldZ(const math::XYZTLorentzVector & position) const;
60 
62 
66  const std::vector<std::unique_ptr<BarrelSimplifiedGeometry> >& barrelLayers() const { return barrelLayers_; }
67 
69 
73  const std::vector<std::unique_ptr<ForwardSimplifiedGeometry> >& forwardLayers() const { return forwardLayers_; }
74 
76 
80  double getMaxRadius() { return maxRadius_;}
81 
83 
87  double getMaxZ() { return maxZ_;}
88 
90  friend std::ostream& operator << (std::ostream& o , const fastsim::Geometry & geometry);
91 
93 
99  {
100  if(layer == nullptr)
101  {
102  return nullptr;
103  }
104  unsigned nextLayerIndex = layer->index() + 1;
105  return nextLayerIndex < barrelLayers_.size() ? barrelLayers_[nextLayerIndex].get() : nullptr;
106  }
107 
109 
115  {
116  if(layer == nullptr)
117  {
118  return nullptr;
119  }
120  unsigned nextLayerIndex = layer->index() + 1;
121  return nextLayerIndex < forwardLayers_.size() ? forwardLayers_[nextLayerIndex].get() : nullptr;
122  }
123 
125 
131  {
132  if(layer == nullptr)
133  {
134  return barrelLayers_.back().get();
135  }
136  return layer->index() > 0 ? barrelLayers_[layer->index() -1].get() : nullptr;
137  }
138 
140 
146  {
147  if(layer == nullptr)
148  {
149  return forwardLayers_.back().get();
150  }
151  return layer->index() > 0 ? forwardLayers_[layer->index() -1].get() : nullptr;
152  }
153 
154  private:
155 
156  std::vector<std::unique_ptr<BarrelSimplifiedGeometry>> barrelLayers_;
157  std::vector<std::unique_ptr<ForwardSimplifiedGeometry>> forwardLayers_;
158  std::unique_ptr<MagneticField> ownedMagneticField_;
159 
162 
166  const double fixedMagneticFieldZ_;
169  const std::vector<edm::ParameterSet> barrelLayerCfg_;
170  const std::vector<edm::ParameterSet> forwardLayerCfg_;
171  const double maxRadius_;
172  const double maxZ_;
173 
174  const bool barrelBoundary_;
175  const bool forwardBoundary_;
178 
179  };
180  std::ostream& operator << (std::ostream& os , const fastsim::Geometry & geometry);
181 }
182 
183 
184 #endif
const bool useTrackerRecoGeometryRecord_
Use GeometricSearchTracker (active layers/reco geometry). Can be used to get position/radius of track...
Definition: Geometry.h:167
unsigned long long cacheIdentifierTrackerRecoGeometry_
Check interval of validity of the tracker geometry.
Definition: Geometry.h:160
Implementation of a forward detector layer (disk).
double getMaxRadius()
Upper bound of the radius of the whole tracker geometry.
Definition: Geometry.h:80
std::vector< std::unique_ptr< BarrelSimplifiedGeometry > > barrelLayers_
The vector of barrel layers (increasing radius)
Definition: Geometry.h:156
const GeometricSearchTracker * geometricSearchTracker_
Definition: Geometry.h:163
const std::vector< std::unique_ptr< ForwardSimplifiedGeometry > > & forwardLayers() const
Return the vector of forward layers (disks).
Definition: Geometry.h:73
const MagneticField * magneticField_
The tracker geometry.
Definition: Geometry.h:164
const BarrelSimplifiedGeometry * nextLayer(const BarrelSimplifiedGeometry *layer) const
Helps to navigate through the vector of barrel layers.
Definition: Geometry.h:98
const std::vector< edm::ParameterSet > forwardLayerCfg_
The config in which all parameters of the forward layers are defined.
Definition: Geometry.h:170
Implementation of a barrel detector layer (cylindrical).
const std::vector< edm::ParameterSet > barrelLayerCfg_
The config in which all parameters of the barrel layers are defined.
Definition: Geometry.h:169
const double fixedMagneticFieldZ_
Use a uniform magnetic field or non-uniform from MagneticFieldRecord.
Definition: Geometry.h:166
const double maxZ_
Upper bound of the radius of the whole tracker geometry.
Definition: Geometry.h:172
const edm::ParameterSet trackerForwardBoundaryCfg_
Hack to interface "old" calo to "new" tracking.
Definition: Geometry.h:177
const ForwardSimplifiedGeometry * nextLayer(const ForwardSimplifiedGeometry *layer) const
Helps to navigate through the vector of forward layers.
Definition: Geometry.h:114
const std::string trackerAlignmentLabel_
The tracker alignment label.
Definition: Geometry.h:168
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
std::unique_ptr< MagneticField > ownedMagneticField_
Needed to create a uniform magnetic field if speciefied in config.
Definition: Geometry.h:158
const bool forwardBoundary_
Hack to interface "old" calo to "new" tracking.
Definition: Geometry.h:175
const BarrelSimplifiedGeometry * previousLayer(const BarrelSimplifiedGeometry *layer) const
Helps to navigate through the vector of barrel layers.
Definition: Geometry.h:130
double getMaxZ()
Upper bound of the length/2 (0 to +Z) of the whole tracker geometry.
Definition: Geometry.h:87
const edm::ParameterSet trackerBarrelBoundaryCfg_
Hack to interface "old" calo to "new" tracking.
Definition: Geometry.h:176
Definition the tracker geometry (vectors of forward/barrel layers).
Definition: Geometry.h:34
const std::vector< std::unique_ptr< BarrelSimplifiedGeometry > > & barrelLayers() const
Return the vector of barrel layers.
Definition: Geometry.h:66
HLT enums.
const double maxRadius_
Definition: Geometry.h:171
int index() const
Return index of this layer (layers are numbered according to their position in the detector)...
#define update(a, b)
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::vector< std::unique_ptr< ForwardSimplifiedGeometry > > forwardLayers_
The vector of forward layers (increasing Z-position)
Definition: Geometry.h:157
const ForwardSimplifiedGeometry * previousLayer(const ForwardSimplifiedGeometry *layer) const
Helps to navigate through the vector of forward layers.
Definition: Geometry.h:145
const bool useFixedMagneticFieldZ_
Needed to create a uniform magnetic field if speciefied in config.
Definition: Geometry.h:165
std::ostream & operator<<(std::ostream &ost, const HLTGlobalStatus &hlt)
Formatted printout of trigger tbale.
const bool barrelBoundary_
Upper bound of the length/2 (0 to +Z) of the whole tracker geometry.
Definition: Geometry.h:174
unsigned long long cacheIdentifierIdealMagneticField_
Check interval of validity of the magnetic field.
Definition: Geometry.h:161