CMS 3D CMS Logo

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