CMS 3D CMS Logo

SimplifiedGeometryFactory.cc
Go to the documentation of this file.
2 
9 
10 #include <TH1F.h>
11 #include <cctype>
12 
14  const GeometricSearchTracker *geometricSearchTracker,
16  const std::map<std::string, fastsim::InteractionModel *> &interactionModelMap,
17  double magneticFieldHistMaxR,
18  double magneticFieldHistMaxZ)
19  : geometricSearchTracker_(geometricSearchTracker),
20  magneticField_(&magneticField),
21  interactionModelMap_(&interactionModelMap),
22  magneticFieldHistMaxR_(magneticFieldHistMaxR),
23  magneticFieldHistMaxZ_(magneticFieldHistMaxZ) {
24  // naming convention for barrel DetLayer lists
28 
29  // naming convention for forward DetLayer lists
36 }
37 
38 std::unique_ptr<fastsim::BarrelSimplifiedGeometry> fastsim::SimplifiedGeometryFactory::createBarrelSimplifiedGeometry(
39  const edm::ParameterSet &cfg) const {
40  std::unique_ptr<fastsim::SimplifiedGeometry> layer = createSimplifiedGeometry(BARREL, cfg);
41  return std::unique_ptr<fastsim::BarrelSimplifiedGeometry>(
42  static_cast<fastsim::BarrelSimplifiedGeometry *>(layer.release()));
43 }
44 
45 std::unique_ptr<fastsim::ForwardSimplifiedGeometry> fastsim::SimplifiedGeometryFactory::createForwardSimplifiedGeometry(
47  if (layerType != NEGFWD && layerType != POSFWD) {
48  throw cms::Exception("fastsim::SimplifiedGeometry::createForwardLayer")
49  << " called with forbidden layerType. Allowed layerTypes are NEGFWD and POSFWD";
50  }
51  std::unique_ptr<fastsim::SimplifiedGeometry> layer = createSimplifiedGeometry(layerType, cfg);
52  return std::unique_ptr<fastsim::ForwardSimplifiedGeometry>(
53  static_cast<fastsim::ForwardSimplifiedGeometry *>(layer.release()));
54 }
55 
56 std::unique_ptr<fastsim::SimplifiedGeometry> fastsim::SimplifiedGeometryFactory::createSimplifiedGeometry(
58  // some flags for internal usage
59  bool isForward = true;
60  bool isOnPositiveSide = false;
61  if (layerType == BARREL) {
62  isForward = false;
63  } else if (layerType == POSFWD) {
64  isOnPositiveSide = true;
65  }
66 
67  // -------------------------------
68  // extract DetLayer (i.e. full geometry of tracker modules)
69  // -------------------------------
70 
71  std::string detLayerName = cfg.getUntrackedParameter<std::string>("activeLayer", "");
72  const DetLayer *detLayer = nullptr;
73 
74  if (!detLayerName.empty() && geometricSearchTracker_) {
75  if (isForward) {
76  detLayerName = (isOnPositiveSide ? "pos" : "neg") + detLayerName;
77  }
78  detLayer = getDetLayer(detLayerName, *geometricSearchTracker_);
79  }
80 
81  // ------------------------------
82  // radius / z of layers
83  // ------------------------------
84 
85  // first try to get it from the configuration
86  double position = 0;
87  std::string positionParameterName = (isForward ? "z" : "radius");
88  if (cfg.exists(positionParameterName)) {
89  position = fabs(cfg.getUntrackedParameter<double>(positionParameterName));
90  if (isForward && !isOnPositiveSide) {
91  position = -position;
92  }
93  }
94  // then try extracting from detLayer
95  else if (detLayer) {
96  if (isForward) {
97  position = static_cast<ForwardDetLayer const *>(detLayer)->surface().position().z();
98  } else {
99  position = static_cast<BarrelDetLayer const *>(detLayer)->specificSurface().radius();
100  }
101  }
102  // then throw error
103  else {
104  std::string cfgString;
105  cfg.allToString(cfgString);
106  throw cms::Exception("fastsim::SimplifiedGeometry")
107  << "Cannot extract a " << (isForward ? "position" : "radius") << " for this "
108  << (isForward ? "forward" : "barrel") << " layer:\n"
109  << cfgString;
110  }
111 
112  // -----------------------------
113  // create the layers
114  // -----------------------------
115 
116  std::unique_ptr<fastsim::SimplifiedGeometry> layer;
117  if (isForward) {
118  layer.reset(new fastsim::ForwardSimplifiedGeometry(position));
119  } else {
120  layer.reset(new fastsim::BarrelSimplifiedGeometry(position));
121  }
122  layer->detLayer_ = detLayer;
123 
124  // -----------------------------
125  // thickness histogram
126  // -----------------------------
127 
128  // Get limits
129  const std::vector<double> &limits = cfg.getUntrackedParameter<std::vector<double> >("limits");
130  // and check order.
131  for (unsigned index = 1; index < limits.size(); index++) {
132  if (limits[index] < limits[index - 1]) {
133  std::string cfgString;
134  cfg.allToString(cfgString);
135  throw cms::Exception("fastsim::SimplifiedGeometryFactory")
136  << "limits must be provided in increasing order. error in:\n"
137  << cfgString;
138  }
139  }
140  // Get thickness values
141  const std::vector<double> &thickness = cfg.getUntrackedParameter<std::vector<double> >("thickness");
142  // and check compatibility with limits
143  if (limits.size() < 2 || thickness.size() != limits.size() - 1) {
144  std::string cfgString;
145  cfg.allToString(cfgString);
146  throw cms::Exception("fastim::SimplifiedGeometryFactory")
147  << "layer thickness and limits not configured properly! error in:" << cfgString;
148  }
149  // create the histogram
150  layer->thicknessHist_.reset(new TH1F("h", "h", limits.size() - 1, &limits[0]));
151  layer->thicknessHist_->SetDirectory(nullptr);
152  for (unsigned i = 1; i < limits.size(); ++i) {
153  layer->thicknessHist_->SetBinContent(i, thickness[i - 1]);
154  }
155 
156  // -----------------------------
157  // nuclear interaction thickness factor
158  // -----------------------------
159 
160  layer->nuclearInteractionThicknessFactor_ =
161  cfg.getUntrackedParameter<double>("nuclearInteractionThicknessFactor", 1.);
162 
163  // -----------------------------
164  // magnetic field
165  // -----------------------------
166 
167  layer->magneticFieldHist_.reset(
168  new TH1F("h", "h", 100, 0., isForward ? magneticFieldHistMaxR_ : magneticFieldHistMaxZ_));
169  layer->magneticFieldHist_->SetDirectory(nullptr);
170  for (int i = 1; i <= 101; i++) {
171  GlobalPoint point = isForward ? GlobalPoint(layer->magneticFieldHist_->GetXaxis()->GetBinCenter(i), 0., position)
172  : GlobalPoint(position, 0., layer->magneticFieldHist_->GetXaxis()->GetBinCenter(i));
173  layer->magneticFieldHist_->SetBinContent(i, magneticField_->inTesla(point).z());
174  }
175 
176  // -----------------------------
177  // list of interaction models
178  // -----------------------------
179 
180  std::vector<std::string> interactionModelLabels =
181  cfg.getUntrackedParameter<std::vector<std::string> >("interactionModels");
182  for (const auto &label : interactionModelLabels) {
183  std::map<std::string, fastsim::InteractionModel *>::const_iterator interactionModel =
185  if (interactionModel == interactionModelMap_->end()) {
186  throw cms::Exception("fastsim::SimplifiedGeometryFactory") << "unknown interaction model '" << label << "'";
187  }
188  layer->interactionModels_.push_back(interactionModel->second);
189  }
190 
191  // -----------------------------
192  // Hack to interface "old" calorimetry with "new" propagation in tracker
193  // -----------------------------
194 
195  if (cfg.exists("caloType")) {
197 
198  if (caloType == "PRESHOWER1") {
199  layer->setCaloType(SimplifiedGeometry::PRESHOWER1);
200  } else if (caloType == "PRESHOWER2") {
201  layer->setCaloType(SimplifiedGeometry::PRESHOWER2);
202  } else if (caloType == "ECAL") {
203  layer->setCaloType(SimplifiedGeometry::ECAL);
204  } else if (caloType == "HCAL") {
205  layer->setCaloType(SimplifiedGeometry::HCAL);
206  } else if (caloType == "VFCAL") {
207  layer->setCaloType(SimplifiedGeometry::VFCAL);
208  } else {
209  throw cms::Exception("fastsim::SimplifiedGeometryFactory")
210  << "unknown caloType '" << caloType << "' (defined PRESHOWER1, PRESHOWER2, ECAL, HCAL, VFCAL)";
211  }
212  }
213 
214  // -----------------------------
215  // and return the layer!
216  // -----------------------------
217 
218  return layer;
219 }
220 
222  const std::string &detLayerName, const GeometricSearchTracker &geometricSearchTracker) const {
223  if (detLayerName.empty()) {
224  return nullptr;
225  }
226 
227  // obtain the index from the detLayerName
228  unsigned pos = detLayerName.size();
229  while (isdigit(detLayerName[pos - 1])) {
230  pos -= 1;
231  }
232  if (pos == detLayerName.size()) {
233  throw cms::Exception("fastsim::SimplifiedGeometry::getDetLayer")
234  << "last part of detLayerName must be index of DetLayer in list. Error in detLayerName" << detLayerName
235  << std::endl;
236  }
237  int index = atoi(detLayerName.substr(pos).c_str());
238  std::string detLayerListName = detLayerName.substr(0, pos);
239 
240  try {
241  // try to find the detLayer in the barrel map
242  if (barrelDetLayersMap_.find(detLayerListName) != barrelDetLayersMap_.end()) {
243  auto detLayerList = barrelDetLayersMap_.find(detLayerListName)->second;
244  return detLayerList->at(index - 1); // use at, to provoce the throwing of an error in case of a bad index
245  }
246 
247  // try to find the detLayer in the forward map
248  else if (forwardDetLayersMap_.find(detLayerListName) != forwardDetLayersMap_.end()) {
249  auto detLayerList = forwardDetLayersMap_.find(detLayerListName)->second;
250  return detLayerList->at(index - 1); // use at, to provoce the throwing of an error in case of a bad index
251  }
252  // throw an error
253  else {
254  throw cms::Exception("fastsim::SimplifiedGeometry::getDetLayer")
255  << " could not find list of detLayers corresponding to detLayerName " << detLayerName << std::endl;
256  }
257  } catch (const std::out_of_range &error) {
258  throw cms::Exception("fastsim::SimplifiedGeometry::getDetLayer")
259  << " index out of range for detLayerName: " << detLayerName << " " << error.what() << std::endl;
260  }
261 }
T getUntrackedParameter(std::string const &, T const &) const
Implementation of a forward detector layer (disk).
std::map< std::string, const std::vector< BarrelDetLayer const * > * > barrelDetLayersMap_
A map of strings and pointers to detLayers.
const double magneticFieldHistMaxZ_
Limit in +-Z for histogram of magnetic field.
const DetLayer * getDetLayer(const std::string &detLayerName, const GeometricSearchTracker &geometricSearchTracker) const
Method returns a pointer to a DetLayer according to the string that was passed.
std::vector< ForwardDetLayer const * > const & posPixelForwardLayers() const
const MagneticField *const magneticField_
The full magnetic field.
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
bool isForward(DetId const &)
bool exists(std::string const &parameterName) const
checks if a parameter exists
Implementation of a barrel detector layer (cylindrical).
std::vector< BarrelDetLayer const * > const & tobLayers() const
SimplifiedGeometryFactory(const GeometricSearchTracker *geometricSearchTracker, const MagneticField &magneticField, const std::map< std::string, fastsim::InteractionModel * > &interactionModelMap, double magneticFieldHistMaxR, double magneticFieldHistMaxZ)
Constructor.
std::vector< ForwardDetLayer const * > const & negPixelForwardLayers() const
char const * label
T z() const
Definition: PV3DBase.h:61
std::vector< BarrelDetLayer const * > const & pixelBarrelLayers() const
const GeometricSearchTracker *const geometricSearchTracker_
The full tracker geometry.
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
std::vector< BarrelDetLayer const * > const & tibLayers() const
virtual const Surface::PositionType & position() const
Returns position of the surface.
std::vector< ForwardDetLayer const * > const & posTecLayers() const
std::unique_ptr< SimplifiedGeometry > createSimplifiedGeometry(LayerType type, const edm::ParameterSet &cfg) const
Main method of this class. Creates a new detector layer (SimplifiedGeometry).
std::vector< ForwardDetLayer const * > const & negTidLayers() const
std::map< std::string, const std::vector< ForwardDetLayer const * > * > forwardDetLayersMap_
A map of strings and pointers to detLayers.
void allToString(std::string &result) const
std::vector< ForwardDetLayer const * > const & posTidLayers() const
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::unique_ptr< ForwardSimplifiedGeometry > createForwardSimplifiedGeometry(LayerType type, const edm::ParameterSet &cfg) const
Helper method for createSimplifiedGeometry(..) to create a forward layer (ForwardSimplifiedGeometry)...
std::unique_ptr< BarrelSimplifiedGeometry > createBarrelSimplifiedGeometry(const edm::ParameterSet &cfg) const
Helper method for createSimplifiedGeometry(..) to create a barrel layer (BarrelSimplifiedGeometry).
const std::map< std::string, fastsim::InteractionModel * > * interactionModelMap_
Map of interaction models.
std::vector< ForwardDetLayer const * > const & negTecLayers() const
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
LayerType
Each layer is either a barrel layer, or a forward layer (either at ppositive or negative Z)...
const double magneticFieldHistMaxR_
Limit in R for histogram of magnetic field.