CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
< BarrelSimplifiedGeometry
createBarrelSimplifiedGeometry (const edm::ParameterSet &cfg) const
 Helper method for createSimplifiedGeometry(..) to create a barrel layer (BarrelSimplifiedGeometry). More...
 
std::unique_ptr
< ForwardSimplifiedGeometry
createForwardSimplifiedGeometry (LayerType type, const edm::ParameterSet &cfg) const
 Helper method for createSimplifiedGeometry(..) to create a forward layer (ForwardSimplifiedGeometry). More...
 
std::unique_ptr
< SimplifiedGeometry
createSimplifiedGeometry (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 14 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().

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 }
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 39 of file SimplifiedGeometryFactory.cc.

References BARREL, and phase1PixelTopology::layer.

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

40  {
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 }
Implementation of a barrel detector layer (cylindrical).
constexpr std::array< uint8_t, layerIndexSize > layer
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 46 of file SimplifiedGeometryFactory.cc.

References Exception, and phase1PixelTopology::layer.

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

47  {
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 }
Implementation of a forward detector layer (disk).
constexpr std::array< uint8_t, layerIndexSize > layer
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 57 of file SimplifiedGeometryFactory.cc.

References edm::ParameterSet::allToString(), BARREL, fastsim::SimplifiedGeometry::ECAL, Exception, edm::ParameterSet::exists(), getDetLayer(), edm::ParameterSet::getUntrackedParameter(), fastsim::SimplifiedGeometry::HCAL, mps_fire::i, ecaldqm::isForward(), label, phase1PixelTopology::layer, TH2PolyOfflineMaps::limits, point, GeometricSearchDet::position(), position, fastsim::SimplifiedGeometry::PRESHOWER1, fastsim::SimplifiedGeometry::PRESHOWER2, AlCaHLTBitMon_QueryRunRegistry::string, TrackerMaterial_cfi::thickness, fastsim::SimplifiedGeometry::VFCAL, and PV3DBase< T, PVType, FrameType >::z().

58  {
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 =
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") {
200  layer->setCaloType(SimplifiedGeometry::PRESHOWER1);
201  } else if (caloType == "PRESHOWER2") {
202  layer->setCaloType(SimplifiedGeometry::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 }
T getUntrackedParameter(std::string const &, T const &) const
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.
virtual GlobalVector inTesla(const GlobalPoint &gp) const =0
Field value ad specified global point, in Tesla.
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
constexpr std::array< uint8_t, layerIndexSize > layer
char const * label
T z() const
Definition: PV3DBase.h:61
const GeometricSearchTracker *const geometricSearchTracker_
The full tracker geometry.
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 222 of file SimplifiedGeometryFactory.cc.

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

223  {
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 }
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 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 SimplifiedGeometryFactory().

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

The full tracker geometry.

Definition at line 96 of file SimplifiedGeometryFactory.h.

Referenced by SimplifiedGeometryFactory().

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

Map of interaction models.

Definition at line 98 of file SimplifiedGeometryFactory.h.

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

The full magnetic field.

Definition at line 97 of file SimplifiedGeometryFactory.h.

const double fastsim::SimplifiedGeometryFactory::magneticFieldHistMaxR_
private

Limit in R for histogram of magnetic field.

Definition at line 99 of file SimplifiedGeometryFactory.h.

const double fastsim::SimplifiedGeometryFactory::magneticFieldHistMaxZ_
private

Limit in +-Z for histogram of magnetic field.

Definition at line 100 of file SimplifiedGeometryFactory.h.