CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Private Member Functions | Private Attributes
fastsim::SimplifiedGeometryFactory Class Reference

Constructs a tracker layer according to entry in python config (incl interaction models). More...

#include <SimplifiedGeometryFactory.h>

Public Types

enum  LayerType { BARREL, POSFWD, NEGFWD }
 Each layer is either a barrel layer, or a forward layer (either at ppositive or negative Z). More...
 

Public Member Functions

std::unique_ptr< BarrelSimplifiedGeometrycreateBarrelSimplifiedGeometry (const edm::ParameterSet &cfg) const
 Helper method for createSimplifiedGeometry(..) to create a barrel layer (BarrelSimplifiedGeometry). More...
 
std::unique_ptr< ForwardSimplifiedGeometrycreateForwardSimplifiedGeometry (LayerType type, const edm::ParameterSet &cfg) const
 Helper method for createSimplifiedGeometry(..) to create a forward layer (ForwardSimplifiedGeometry). More...
 
std::unique_ptr< SimplifiedGeometrycreateSimplifiedGeometry (LayerType type, const edm::ParameterSet &cfg) const
 Main method of this class. Creates a new detector layer (SimplifiedGeometry). More...
 
 SimplifiedGeometryFactory (const GeometricSearchTracker *geometricSearchTracker, const MagneticField &magneticField, const std::map< std::string, fastsim::InteractionModel * > &interactionModelMap, double magneticFieldHistMaxR, double magneticFieldHistMaxZ)
 Constructor. More...
 

Private Member Functions

const DetLayergetDetLayer (const std::string &detLayerName, const GeometricSearchTracker &geometricSearchTracker) const
 Method returns a pointer to a DetLayer according to the string that was passed. More...
 

Private Attributes

std::map< std::string, const std::vector< BarrelDetLayer const * > * > barrelDetLayersMap_
 A map of strings and pointers to detLayers. More...
 
std::map< std::string, const std::vector< ForwardDetLayer const * > * > forwardDetLayersMap_
 A map of strings and pointers to detLayers. More...
 
const GeometricSearchTracker *const geometricSearchTracker_
 The full tracker geometry. More...
 
const std::map< std::string, fastsim::InteractionModel * > * interactionModelMap_
 Map of interaction models. More...
 
const MagneticField *const magneticField_
 The full magnetic field. More...
 
const double magneticFieldHistMaxR_
 Limit in R for histogram of magnetic field. More...
 
const double magneticFieldHistMaxZ_
 Limit in +-Z for histogram of magnetic field. More...
 

Detailed Description

Constructs a tracker layer according to entry in python config (incl interaction models).

Also creates links to DetLayer (if active layer) and stores strength of magnetic field along the layer. If some parameters are not stored in the config file, tries to get them from the full detector geometry (GeometricSearchTracker). This is however only possible for active layers.

See also
Geometry()
SimplifiedGeometry()

Definition at line 40 of file SimplifiedGeometryFactory.h.

Member Enumeration Documentation

Each layer is either a barrel layer, or a forward layer (either at ppositive or negative Z).

Enumerator
BARREL 
POSFWD 
NEGFWD 

Definition at line 57 of file SimplifiedGeometryFactory.h.

Constructor & Destructor Documentation

fastsim::SimplifiedGeometryFactory::SimplifiedGeometryFactory ( const GeometricSearchTracker geometricSearchTracker,
const MagneticField magneticField,
const std::map< std::string, fastsim::InteractionModel * > &  interactionModelMap,
double  magneticFieldHistMaxR,
double  magneticFieldHistMaxZ 
)

Constructor.

Parameters
geometricSearchTrackerThe full tracker geometry (needed for links to active detLayers).
magneticFieldThe full magnetic field.
interactionModelMapMap of interaction models that should be assigned for that layer.
magneticFieldHistMaxRMax Radius for initialization of magnetic field histogram (TH1, limit of axis).
magneticFieldHistMaxZMax Z for initialization of magnetic field histogram (TH1, limit of axis).

Definition at line 13 of file SimplifiedGeometryFactory.cc.

References barrelDetLayersMap_, forwardDetLayersMap_, geometricSearchTracker_, GeometricSearchTracker::negPixelForwardLayers(), GeometricSearchTracker::negTecLayers(), GeometricSearchTracker::negTidLayers(), GeometricSearchTracker::pixelBarrelLayers(), GeometricSearchTracker::posPixelForwardLayers(), GeometricSearchTracker::posTecLayers(), GeometricSearchTracker::posTidLayers(), GeometricSearchTracker::tibLayers(), and GeometricSearchTracker::tobLayers().

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 }
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.
std::vector< ForwardDetLayer const * > const & posPixelForwardLayers() const
const MagneticField *const magneticField_
The full magnetic field.
std::vector< BarrelDetLayer const * > const & tobLayers() const
std::vector< ForwardDetLayer const * > const & negPixelForwardLayers() const
std::vector< BarrelDetLayer const * > const & pixelBarrelLayers() const
const GeometricSearchTracker *const geometricSearchTracker_
The full tracker geometry.
std::vector< BarrelDetLayer const * > const & tibLayers() const
std::vector< ForwardDetLayer const * > const & posTecLayers() const
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.
std::vector< ForwardDetLayer const * > const & posTidLayers() const
const std::map< std::string, fastsim::InteractionModel * > * interactionModelMap_
Map of interaction models.
std::vector< ForwardDetLayer const * > const & negTecLayers() const
const double magneticFieldHistMaxR_
Limit in R for histogram of magnetic field.

Member Function Documentation

std::unique_ptr< fastsim::BarrelSimplifiedGeometry > fastsim::SimplifiedGeometryFactory::createBarrelSimplifiedGeometry ( const edm::ParameterSet cfg) const

Helper method for createSimplifiedGeometry(..) to create a barrel layer (BarrelSimplifiedGeometry).

Returns
A BarrelSimplifiedGeometry
See also
createSimplifiedGeometry(LayerType type, const edm::ParameterSet & cfg)

Definition at line 38 of file SimplifiedGeometryFactory.cc.

References BARREL, and createSimplifiedGeometry().

Referenced by fastsim::Geometry::update().

39  {
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 }
Implementation of a barrel detector layer (cylindrical).
std::unique_ptr< SimplifiedGeometry > createSimplifiedGeometry(LayerType type, const edm::ParameterSet &cfg) const
Main method of this class. Creates a new detector layer (SimplifiedGeometry).
std::unique_ptr< fastsim::ForwardSimplifiedGeometry > fastsim::SimplifiedGeometryFactory::createForwardSimplifiedGeometry ( LayerType  type,
const edm::ParameterSet cfg 
) const

Helper method for createSimplifiedGeometry(..) to create a forward layer (ForwardSimplifiedGeometry).

Parameters
typeEither POSFWD or NEGFWD.
Returns
A ForwardSimplifiedGeometry
See also
createSimplifiedGeometry(LayerType type, const edm::ParameterSet & cfg)

Definition at line 45 of file SimplifiedGeometryFactory.cc.

References createSimplifiedGeometry(), Exception, NEGFWD, and POSFWD.

Referenced by fastsim::Geometry::update().

46  {
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 }
Implementation of a forward detector layer (disk).
std::unique_ptr< SimplifiedGeometry > createSimplifiedGeometry(LayerType type, const edm::ParameterSet &cfg) const
Main method of this class. Creates a new detector layer (SimplifiedGeometry).
std::unique_ptr< fastsim::SimplifiedGeometry > fastsim::SimplifiedGeometryFactory::createSimplifiedGeometry ( LayerType  type,
const edm::ParameterSet cfg 
) const

Main method of this class. Creates a new detector layer (SimplifiedGeometry).

Reads the config file, does all the initialization etc. and creates either a forward or a barrel layer (depends on LayerType type).

Parameters
typeEither BARREL, POSFWD or NEGFWD.
Returns
A SimplifiedGeometry (either ForwardSimplifiedGeometry or BarrelSimplifiedGeometry).

Definition at line 56 of file SimplifiedGeometryFactory.cc.

References edm::ParameterSet::allToString(), BARREL, CaloMaterial_cfi::caloType, fastsim::SimplifiedGeometry::ECAL, Exception, edm::ParameterSet::exists(), geometricSearchTracker_, getDetLayer(), edm::ParameterSet::getUntrackedParameter(), fastsim::SimplifiedGeometry::HCAL, mps_fire::i, interactionModelMap_, MagneticField::inTesla(), ecaldqm::isForward(), label, TH2PolyOfflineMaps::limits, magneticField_, magneticFieldHistMaxR_, magneticFieldHistMaxZ_, point, POSFWD, GeometricSearchDet::position(), position, fastsim::SimplifiedGeometry::PRESHOWER1, fastsim::SimplifiedGeometry::PRESHOWER2, AlCaHLTBitMon_QueryRunRegistry::string, Calorimetry_cff::thickness, fastsim::SimplifiedGeometry::VFCAL, and PV3DBase< T, PVType, FrameType >::z().

Referenced by createBarrelSimplifiedGeometry(), and createForwardSimplifiedGeometry().

57  {
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 }
T getUntrackedParameter(std::string const &, T const &) const
Implementation of a forward detector layer (disk).
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.
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).
char const * label
T z() const
Definition: PV3DBase.h:61
const GeometricSearchTracker *const geometricSearchTracker_
The full tracker geometry.
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
virtual const Surface::PositionType & position() const
Returns position of the surface.
void allToString(std::string &result) const
static int position[264][3]
Definition: ReadPGInfo.cc:289
const std::map< std::string, fastsim::InteractionModel * > * interactionModelMap_
Map of interaction models.
*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
const double magneticFieldHistMaxR_
Limit in R for histogram of magnetic field.
const DetLayer * fastsim::SimplifiedGeometryFactory::getDetLayer ( const std::string &  detLayerName,
const GeometricSearchTracker geometricSearchTracker 
) const
private

Method returns a pointer to a DetLayer according to the string that was passed.

A convention for the name of the layer is used. For barrel layers this is "XXX?" where XXX is a part of the tracker and ? is the index of the layer (starting at one). For forward layers one has to add neg/pos in front to distinguish between the disk at -Z and +Z spatial position, so the convention is "xxxXXX?" Valid names: BPix, TIB, TOB, negFPix, posFPix, negTID, posTID, negTEC, posTEC Accordingly, the innermost layer of the barrel pixel detector is "BPix1".

Parameters
detLayerNameA string following the naming convention.

Definition at line 221 of file SimplifiedGeometryFactory.cc.

References barrelDetLayersMap_, relativeConstraints::error, Exception, forwardDetLayersMap_, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by createSimplifiedGeometry().

222  {
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 }
std::map< std::string, const std::vector< BarrelDetLayer const * > * > barrelDetLayersMap_
A map of strings and pointers to detLayers.
std::map< std::string, const std::vector< ForwardDetLayer const * > * > forwardDetLayersMap_
A map of strings and pointers to detLayers.

Member Data Documentation

std::map<std::string, const std::vector<BarrelDetLayer const *> *> fastsim::SimplifiedGeometryFactory::barrelDetLayersMap_
private

A map of strings and pointers to detLayers.

Definition at line 102 of file SimplifiedGeometryFactory.h.

Referenced by getDetLayer(), and SimplifiedGeometryFactory().

std::map<std::string, const std::vector<ForwardDetLayer const *> *> fastsim::SimplifiedGeometryFactory::forwardDetLayersMap_
private

A map of strings and pointers to detLayers.

Definition at line 104 of file SimplifiedGeometryFactory.h.

Referenced by getDetLayer(), and SimplifiedGeometryFactory().

const GeometricSearchTracker* const fastsim::SimplifiedGeometryFactory::geometricSearchTracker_
private

The full tracker geometry.

Definition at line 96 of file SimplifiedGeometryFactory.h.

Referenced by createSimplifiedGeometry(), and SimplifiedGeometryFactory().

const std::map<std::string, fastsim::InteractionModel *>* fastsim::SimplifiedGeometryFactory::interactionModelMap_
private

Map of interaction models.

Definition at line 98 of file SimplifiedGeometryFactory.h.

Referenced by createSimplifiedGeometry().

const MagneticField* const fastsim::SimplifiedGeometryFactory::magneticField_
private

The full magnetic field.

Definition at line 97 of file SimplifiedGeometryFactory.h.

Referenced by createSimplifiedGeometry().

const double fastsim::SimplifiedGeometryFactory::magneticFieldHistMaxR_
private

Limit in R for histogram of magnetic field.

Definition at line 99 of file SimplifiedGeometryFactory.h.

Referenced by createSimplifiedGeometry().

const double fastsim::SimplifiedGeometryFactory::magneticFieldHistMaxZ_
private

Limit in +-Z for histogram of magnetic field.

Definition at line 100 of file SimplifiedGeometryFactory.h.

Referenced by createSimplifiedGeometry().