CMS 3D CMS Logo

ME0ReDigiProducer.cc
Go to the documentation of this file.
1 #ifndef SimMuon_GEMDigitizer_ME0ReDigiProducer_h
2 #define SimMuon_GEMDigitizer_ME0ReDigiProducer_h
3 
4 /*
5  * This module smears and discretizes the timing and position of the
6  * ME0 pseudo digis.
7  */
8 
21 
24 
29 
30 #include "CLHEP/Random/RandGaussQ.h"
31 #include "CLHEP/Random/RandFlat.h"
32 #include "CLHEP/Units/PhysicalConstants.h"
33 
34 #include <sstream>
35 #include <string>
36 #include <map>
37 #include <vector>
38 
40 
42 private:
43  //Class used to define custom geometry, if required
44  //assume that all ME0 chambers have the same dimension
45  //and for the same layer have the same radial and Z position
46  //Good for now, can build in support for more varied geos later
47  //if necessary
49  public:
51  const unsigned int numberOfStrips,
52  const unsigned int numberOfPartitions);
54  unsigned int findEtaPartition(float locY) const;
55  const TrapezoidalStripTopology* getTopo(const unsigned int partIdx) const { return stripTopos[partIdx]; }
56  float getPartCenter(const unsigned int partIdx) const; //position of part. in chamber
57  float getCentralTOF(const ME0DetId& me0Id, unsigned int partIdx) const {
58  return tofs[me0Id.layer() - 1][partIdx];
59  } //in detId layer numbers stat at 1
60  unsigned int numLayers() const { return tofs.size(); }
61 
62  private:
63  TrapezoidalStripTopology* buildTopo(const std::vector<float>& _p) const;
64 
65  private:
66  float middleDistanceFromBeam; // radiusOfMainPartitionInCenter;
67  std::vector<TrapezoidalStripTopology*> stripTopos; // vector of Topos, one for each part
68  std::vector<std::vector<double>> tofs; //TOF to center of the partition: [layer][part]
69  std::vector<float> partitionTops; //Top of each partition in the chamber's local coords
70  };
71 
72 public:
73  explicit ME0ReDigiProducer(const edm::ParameterSet& ps);
74 
75  ~ME0ReDigiProducer() override;
76 
77  void beginRun(const edm::Run&, const edm::EventSetup&) override;
78 
79  void produce(edm::Event&, const edm::EventSetup&) override;
80 
84  CLHEP::HepRandomEngine* engine);
85 
86 private:
87  void fillCentralTOFs();
88  void getStripProperties(const ME0EtaPartition* etaPart,
89  const ME0DigiPreReco* inDigi,
90  float& tof,
91  int& strip,
92  LocalPoint& digiLocalPoint,
93  LocalError& digiLocalError) const;
95  const ME0DigiPreReco* inDigi,
96  float& tof,
97  int& strip,
98  LocalPoint& digiLocalPoint,
99  LocalError& digiLocalError) const;
100 
101  typedef std::tuple<unsigned int, unsigned int, unsigned int> DigiIndicies;
102  typedef std::map<DigiIndicies, unsigned int> ChamberDigiMap;
103  //fills map...returns -1 if digi is not already in the map
104  unsigned int fillDigiMap(
105  ChamberDigiMap& chDigiMap, unsigned int bx, unsigned int part, unsigned int strip, unsigned int currentIDX) const;
106 
107  //paramters
108  const float bxWidth; // ns
109  bool useCusGeoFor1PartGeo; //Use custom strips and partitions for digitization for single partition geometry
110  bool usePads; //sets strip granularity to x2 coarser
111  unsigned int numberOfStrips; // Custom number of strips per partition
112  unsigned int numberOfPartitions; // Custom number of partitions per chamber
113  double neutronAcceptance; // fraction of neutron events to keep in event (>= 1 means no filtering)
114  double timeResolution; // smear time by gaussian with this sigma (in ns)....negative for no smearing
115  int minBXReadout; // Minimum BX to readout
116  int maxBXReadout; // Maximum BX to readout
117  std::vector<int> layerReadout; // Don't readout layer if entry is 0 (Layer number 1 in the numbering scheme is idx 0)
118  bool mergeDigis; // Keep only one digi at the same chamber, strip, partition, and BX
121 
125  std::vector<std::vector<double>> tofs; //used for built in geo
126 };
127 
129  const unsigned int numberOfStrips,
130  const unsigned int numberOfPartitions) {
131  //First test geometry to make sure that it is compatible with our assumptions
132  const auto& chambers = geometry->chambers();
133  if (chambers.empty())
134  throw cms::Exception("Setup")
135  << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - No ME0Chambers in geometry.";
136  const auto* mainChamber = chambers.front();
137  const unsigned int nLayers = chambers.front()->nLayers();
138  if (!nLayers)
139  throw cms::Exception("Setup")
140  << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - ME0Chamber has no layers.";
141  const auto* mainLayer = mainChamber->layers()[0];
142  if (!mainLayer->nEtaPartitions())
143  throw cms::Exception("Setup")
144  << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - ME0Layer has no partitions.";
145  if (mainLayer->nEtaPartitions() != 1)
146  throw cms::Exception("Setup") << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - This module is only "
147  "compatitble with geometries that contain only one partition per ME0Layer.";
148 
149  const auto* mainPartition = mainLayer->etaPartitions()[0];
150  const TrapezoidalStripTopology* mainTopo = dynamic_cast<const TrapezoidalStripTopology*>(&mainPartition->topology());
151  if (!mainTopo)
152  throw cms::Exception("Setup") << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - ME0 strip topology "
153  "must be of type TrapezoidalStripTopology. This module cannot be used";
154 
155  for (const auto& chamber : geometry->chambers()) {
156  if (chamber->nLayers() != int(nLayers))
157  throw cms::Exception("Setup") << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - Not all "
158  "ME0Chambers have the same number of layers. This module cannot be used.";
159  for (unsigned int iL = 0; iL < nLayers; ++iL) {
160  if (chamber->layers()[iL]->nEtaPartitions() != mainLayer->nEtaPartitions())
161  throw cms::Exception("Setup") << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - Not all "
162  "ME0Layers have the same number of partitions. This module cannot be used.";
163  if (chamber->layers()[iL]->etaPartitions()[0]->specs()->parameters() != mainPartition->specs()->parameters())
164  throw cms::Exception("Setup") << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - Not all ME0 ETA "
165  "partitions have the same properties. This module cannot be used.";
166  if (std::fabs(chamber->layers()[iL]->etaPartitions()[0]->position().z()) !=
167  std::fabs(mainChamber->layers()[iL]->etaPartitions()[0]->position().z()))
168  throw cms::Exception("Setup")
169  << "ME0ReDigiProducer::TemporaryGeometry::TemporaryGeometry() - Not all ME0 ETA partitions in a single "
170  "layer have the same Z position. This module cannot be used.";
171  }
172  }
173 
174  //Calculate radius to center of partition
175  middleDistanceFromBeam = mainTopo->radius();
176 
177  //calculate the top of each eta partition, assuming equal distance in eta between partitions
178  const auto localTop = LocalPoint(0, mainTopo->stripLength() / 2);
179  const auto localBottom = LocalPoint(0, -1 * mainTopo->stripLength() / 2);
180  const auto globalTop = mainPartition->toGlobal(localTop);
181  const auto globalBottom = mainPartition->toGlobal(localBottom);
182  const double etaTop = globalTop.eta();
183  const double etaBottom = globalBottom.eta();
184  const double zBottom = globalBottom.z();
185 
186  //Build topologies
189  const auto& mainPars = mainPartition->specs()->parameters();
190  for (unsigned int iP = 0; iP < numberOfPartitions; ++iP) {
191  const double eta = (etaTop - etaBottom) * double(iP + 1) / double(numberOfPartitions) + etaBottom;
192  const double distFromBeam = std::fabs(zBottom / std::sinh(eta));
193  partitionTops.push_back(distFromBeam - middleDistanceFromBeam);
194  LogDebug("ME0ReDigiProducer::TemporaryGeometry") << "Top of new partition: " << partitionTops.back() << std::endl;
195 
196  std::vector<float> params(4, 0);
197 
198  //half width of trapezoid at local coordinate Y
199  auto getWidth = [&](float locY) -> float {
200  return (mainPars[2] * (mainPars[1] + mainPars[0]) + locY * (mainPars[1] - mainPars[0])) / (2 * mainPars[2]);
201  };
202 
203  params[0] = iP == 0 ? mainPars[0] : getWidth(partitionTops[iP - 1]); // Half width of bottom of chamber
204  params[1] =
205  iP + 1 == numberOfPartitions ? mainPars[1] : getWidth(partitionTops[iP]); // Half width of top of chamber
206  params[2] = ((iP + 1 == numberOfPartitions ? localTop.y() : partitionTops[iP]) -
207  (iP == 0 ? localBottom.y() : partitionTops[iP - 1])) /
208  2; // Half width of height of chamber
209  params[3] = numberOfStrips;
210 
211  stripTopos.push_back(buildTopo(params));
212  }
213 
214  //Get TOF at center of each partition
215  tofs.resize(nLayers);
216  LogDebug("ME0ReDigiProducer::TemporaryGeometry") << "TOF numbers [layer][partition]: ";
217  for (unsigned int iL = 0; iL < nLayers; ++iL) {
218  tofs[iL].resize(numberOfPartitions);
219  for (unsigned int iP = 0; iP < numberOfPartitions; ++iP) {
220  const LocalPoint partCenter(0., getPartCenter(iP), 0.);
221  const GlobalPoint centralGP(mainChamber->layers()[iL]->etaPartitions()[0]->toGlobal(partCenter));
222  tofs[iL][iP] = (centralGP.mag() / (CLHEP::c_light / CLHEP::cm)); //speed of light [cm/ns]
223  LogDebug("ME0ReDigiProducer::TemporaryGeometry")
224  << "[" << iL << "][" << iP << "]=" << tofs[iL][iP] << " " << std::endl;
225  }
226  }
227 }
228 
230  unsigned int etaPart = stripTopos.size() - 1;
231  for (unsigned int iP = 0; iP < stripTopos.size(); ++iP) {
232  if (locY < partitionTops[iP]) {
233  etaPart = iP;
234  break;
235  }
236  }
237  return etaPart;
238 }
239 
241  return stripTopos[partIdx]->radius() - middleDistanceFromBeam;
242 }
243 
245  for (auto* p : stripTopos) {
246  delete p;
247  }
248 }
249 
251  float b = _p[0];
252  float B = _p[1];
253  float h = _p[2];
254  float r0 = h * (B + b) / (B - b);
255  float striplength = h * 2;
256  float strips = _p[3];
257  float pitch = (b + B) / strips;
258  int nstrip = static_cast<int>(strips);
259 
260  LogDebug("ME0ReDigiProducer::TemporaryGeometry")
261  << "New partition parameters: "
262  << "bottom width(" << 2 * b << ") top width(" << 2 * B << ") height(" << 2 * h << ") radius to center(" << r0
263  << ") nStrips(" << strips << ") pitch(" << pitch << ")" << std::endl;
264 
265  return new TrapezoidalStripTopology(nstrip, pitch, striplength, r0);
266 }
267 
269  : bxWidth(25.0),
270  useCusGeoFor1PartGeo(ps.getParameter<bool>("useCusGeoFor1PartGeo")),
271  usePads(ps.getParameter<bool>("usePads")),
272  numberOfStrips(ps.getParameter<unsigned int>("numberOfStrips")),
273  numberOfPartitions(ps.getParameter<unsigned int>("numberOfPartitions")),
274  neutronAcceptance(ps.getParameter<double>("neutronAcceptance")),
275  timeResolution(ps.getParameter<double>("timeResolution")),
276  minBXReadout(ps.getParameter<int>("minBXReadout")),
277  maxBXReadout(ps.getParameter<int>("maxBXReadout")),
278  layerReadout(ps.getParameter<std::vector<int>>("layerReadout")),
279  mergeDigis(ps.getParameter<bool>("mergeDigis")),
280  token(consumes<ME0DigiPreRecoCollection>(edm::InputTag(ps.getParameter<std::string>("inputCollection")))),
282  produces<ME0DigiPreRecoCollection>();
283  produces<ME0DigiPreRecoMap>();
284 
286  if (!rng.isAvailable()) {
287  throw cms::Exception("Configuration")
288  << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - RandomNumberGeneratorService is not present in configuration "
289  "file.\n"
290  << "Add the service in the configuration file or remove the modules that require it.";
291  }
292  geometry = nullptr;
293  tempGeo = nullptr;
294  useBuiltinGeo = true;
295 
296  if (useCusGeoFor1PartGeo) {
297  if (usePads)
299  if (numberOfStrips == 0)
300  throw cms::Exception("Setup")
301  << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - Must have at least one strip if using custom geometry.";
302  if (numberOfPartitions == 0)
303  throw cms::Exception("Setup")
304  << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - Must have at least one partition if using custom geometry.";
305  }
306 
307  if (neutronAcceptance < 0)
308  throw cms::Exception("Setup") << "ME0ReDigiProducer::ME0PreRecoDigiProducer() - neutronAcceptance must be >= 0.";
309 }
310 
312  if (tempGeo)
313  delete tempGeo;
314 }
315 
317  // set geometry
319  geometry = &*hGeom;
320 
321  const auto& chambers = geometry->chambers();
322  if (chambers.empty())
323  throw cms::Exception("Setup") << "ME0ReDigiProducer::beginRun() - No ME0Chambers in geometry.";
324 
325  const unsigned int nLayers = chambers.front()->nLayers();
326  if (!nLayers)
327  throw cms::Exception("Setup") << "ME0ReDigiProducer::beginRun() - No layers in ME0 geometry.";
328 
329  const unsigned int nPartitions = chambers.front()->layers()[0]->nEtaPartitions();
330 
331  if (useCusGeoFor1PartGeo && nPartitions == 1) {
332  useBuiltinGeo = false;
333  }
334 
335  if (useBuiltinGeo) {
336  if (nLayers != layerReadout.size())
337  throw cms::Exception("Configuration")
338  << "ME0ReDigiProducer::beginRun() - The geometry has " << nLayers << " layers, but the readout of "
339  << layerReadout.size() << " were specified with the layerReadout parameter.";
340  fillCentralTOFs();
341  } else {
342  LogDebug("ME0ReDigiProducer") << "Building temporary geometry:" << std::endl;
344  LogDebug("ME0ReDigiProducer") << "Done building temporary geometry!" << std::endl;
345 
346  if (tempGeo->numLayers() != layerReadout.size())
347  throw cms::Exception("Configuration")
348  << "ME0ReDigiProducer::beginRun() - The geometry has " << tempGeo->numLayers()
349  << " layers, but the readout of " << layerReadout.size()
350  << " were specified with the layerReadout parameter.";
351  }
352 }
353 
356  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
357 
359  e.getByToken(token, input_digis);
360 
361  std::unique_ptr<ME0DigiPreRecoCollection> output_digis(new ME0DigiPreRecoCollection());
362  std::unique_ptr<ME0DigiPreRecoMap> output_digimap(new ME0DigiPreRecoMap());
363 
364  // build the digis
365  buildDigis(*(input_digis.product()), *output_digis, *output_digimap, engine);
366 
367  // store them in the event
368  e.put(std::move(output_digis));
369  e.put(std::move(output_digimap));
370 }
371 
373  ME0DigiPreRecoCollection& output_digis,
374  ME0DigiPreRecoMap& output_digimap,
375  CLHEP::HepRandomEngine* engine) {
376  /*
377  Starting form the incoming pseudo-digi, which has perfect time and position resolution:
378  1A. Smear time using sigma_t by some value
379  1B. Correct the smeared time with the central arrival time for partition
380  1C. Apply discretization: if the smeared time is outside the BX window (-12.5ns;+12.5ns),
381  the hit should be assigned to the next (or previous) BX
382 
383  2A. Find strip that the digi belongs to
384  2B. Get the center of this strip and the error on the position assuming the geometry
385 
386  3A. Filter event if a digi at this partition/strip/BX already exists
387  3B. Add to collection
388  */
389 
390  LogDebug("ME0ReDigiProducer::buildDigis") << "Begin building digis." << std::endl;
392  for (me0dgIt = input_digis.begin(); me0dgIt != input_digis.end(); ++me0dgIt) {
393  const auto& me0Id = (*me0dgIt).first;
394  LogTrace("ME0ReDigiProducer::buildDigis") << "Starting with roll: " << me0Id << std::endl;
395 
396  //setup map for this chamber/eta partition
397  ChamberDigiMap chDigiMap;
398 
399  int newDigiIdx = 0;
400  const ME0DigiPreRecoCollection::Range& range = (*me0dgIt).second;
401  for (ME0DigiPreRecoCollection::const_iterator digi = range.first; digi != range.second; digi++) {
402  LogTrace("ME0ReDigiProducer::buildDigis") << std::endl
403  << "(" << digi->x() << "," << digi->y() << "," << digi->tof() << ","
404  << digi->pdgid() << "," << digi->prompt() << ")-> ";
405 
406  //If we don't readout this layer skip
407  if (!layerReadout[me0Id.layer() - 1]) {
408  output_digimap.insertDigi(me0Id, -1);
409  continue;
410  }
411 
412  //if neutron and we are filtering skip
413  if (!digi->prompt() && neutronAcceptance < 1.0)
414  if (CLHEP::RandFlat::shoot(engine) > neutronAcceptance) {
415  output_digimap.insertDigi(me0Id, -1);
416  continue;
417  }
418 
419  //smear TOF if necessary
420  float tof = digi->tof() + (timeResolution < 0 ? 0.0 : CLHEP::RandGaussQ::shoot(engine, 0, timeResolution));
421 
422  //Values used to fill objet
423  int mapPartIDX = me0Id.roll() - 1;
424  int strip = 0;
425  LocalPoint digiLocalPoint;
426  LocalError digiLocalError;
427  if (useBuiltinGeo) {
428  getStripProperties(geometry->etaPartition(me0Id), &*digi, tof, strip, digiLocalPoint, digiLocalError);
429  } else {
430  mapPartIDX = getCustomStripProperties(me0Id, &*digi, tof, strip, digiLocalPoint, digiLocalError);
431  }
432 
433  //filter if outside of readout window
434  const int bxIdx = std::round(tof / bxWidth);
435  LogTrace("ME0ReDigiProducer::buildDigis") << tof << "(" << bxIdx << ") ";
436  if (bxIdx < minBXReadout) {
437  output_digimap.insertDigi(me0Id, -1);
438  continue;
439  }
440  if (bxIdx > maxBXReadout) {
441  output_digimap.insertDigi(me0Id, -1);
442  continue;
443  }
444  tof = bxIdx * bxWidth;
445 
446  //If we are merging check to see if it already exists
447  LogTrace("ME0ReDigiProducer::buildDigis") << "(" << bxIdx << "," << mapPartIDX << "," << strip << ") ";
448  if (mergeDigis) {
449  int matchIDX = fillDigiMap(chDigiMap, bxIdx, mapPartIDX, strip, newDigiIdx);
450  if (matchIDX >= 0) {
451  output_digimap.insertDigi(me0Id, matchIDX);
452  continue;
453  }
454  }
455 
456  //Digis store sigmaX,sigmaY, correlationCoef
457  const float sigmaX = std::sqrt(digiLocalError.xx());
458  const float sigmaY = std::sqrt(digiLocalError.yy());
459  const float corrCoef = digiLocalError.xy() / (sigmaX * sigmaY);
460 
461  //Fill in the new collection
462  ME0DigiPreReco out_digi(
463  digiLocalPoint.x(), digiLocalPoint.y(), sigmaX, sigmaY, corrCoef, tof, digi->pdgid(), digi->prompt());
464  output_digis.insertDigi(me0Id, out_digi);
465 
466  // store index of previous detid and digi
467  output_digimap.insertDigi(me0Id, newDigiIdx);
468  newDigiIdx++;
469 
470  LogTrace("ME0ReDigiProducer::buildDigis") << "(" << digiLocalPoint.x() << "," << digiLocalPoint.y() << ","
471  << sigmaX << "," << sigmaY << "," << tof << ") ";
472  }
473 
474  chDigiMap.clear();
475  }
476 }
477 
479  const auto* mainChamber = geometry->chambers().front();
480  const unsigned int nLayers = mainChamber->nLayers();
481  //Get TOF at center of each partition
482  tofs.clear();
483  tofs.resize(nLayers);
484  LogDebug("ME0ReDigiProducer::fillCentralTOFs()") << "TOF numbers [layer][partition]: ";
485  for (unsigned int iL = 0; iL < nLayers; ++iL) {
486  const auto* layer = mainChamber->layers()[iL];
487  const unsigned int mapLayIDX = layer->id().layer() - 1;
488  const unsigned int nPartitions = layer->nEtaPartitions();
489  if (!nPartitions)
490  throw cms::Exception("Setup") << "ME0ReDigiProducer::fillCentralTOFs() - ME0Layer has no partitions.";
491  tofs[mapLayIDX].resize(nPartitions);
492  for (unsigned int iP = 0; iP < nPartitions; ++iP) {
493  const unsigned int mapPartIDX = layer->etaPartitions()[iP]->id().roll() - 1;
494  const GlobalPoint centralGP(layer->etaPartitions()[iP]->position());
495  tofs[mapLayIDX][mapPartIDX] = (centralGP.mag() / (CLHEP::c_light / CLHEP::cm)); //speed of light [cm/ns]
496  LogDebug("ME0ReDigiProducer::fillCentralTOFs()")
497  << "[" << mapLayIDX << "][" << mapPartIDX << "]=" << tofs[mapLayIDX][mapPartIDX] << " " << std::endl;
498  }
499  }
500 }
502  const ME0DigiPreReco* inDigi,
503  float& tof,
504  int& strip,
505  LocalPoint& digiLocalPoint,
506  LocalError& digiLocalError) const {
507  const unsigned int partIdx = tempGeo->findEtaPartition(inDigi->y());
508  LogTrace("ME0ReDigiProducer::buildDigis") << partIdx << " ";
509  const float partMeanTof = tempGeo->getCentralTOF(detId, partIdx);
510 
511  //convert to relative to partition
512  tof -= partMeanTof;
513 
514  //get coordinates and errors
515  const float partCenter = tempGeo->getPartCenter(partIdx);
516  const auto* topo = tempGeo->getTopo(partIdx);
517 
518  //find channel
519  const LocalPoint partLocalPoint(inDigi->x(), inDigi->y() - partCenter, 0.);
520  strip = topo->channel(partLocalPoint);
521  const float stripF = float(strip) + 0.5;
522 
523  //get digitized location
524  LocalPoint digiPartLocalPoint = topo->localPosition(stripF);
525  digiLocalError = topo->localError(
526  stripF,
527  1. /
528  sqrt(
529  12.)); //std dev. flat distribution with length L is L/sqrt(12). The strip topology expects the error in units of strips.
530  digiLocalPoint = LocalPoint(digiPartLocalPoint.x(), digiPartLocalPoint.y() + partCenter, 0.0);
531  return partIdx;
532 }
534  const ME0DigiPreReco* inDigi,
535  float& tof,
536  int& strip,
537  LocalPoint& digiLocalPoint,
538  LocalError& digiLocalError) const {
539  //convert to relative to partition
540  tof -= tofs[etaPart->id().layer() - 1][etaPart->id().roll() - 1];
541 
542  const TrapezoidalStripTopology* origTopo = (const TrapezoidalStripTopology*)(&etaPart->specificTopology());
543  TrapezoidalStripTopology padTopo(
544  origTopo->nstrips() / 2, origTopo->pitch() * 2, origTopo->stripLength(), origTopo->radius());
545  const auto& topo = usePads ? padTopo : etaPart->specificTopology();
546 
547  //find channel
548  const LocalPoint partLocalPoint(inDigi->x(), inDigi->y(), 0.);
549  strip = topo.channel(partLocalPoint);
550  const float stripF = float(strip) + 0.5;
551 
552  //get digitized location
553  digiLocalPoint = topo.localPosition(stripF);
554  digiLocalError = topo.localError(stripF, 1. / sqrt(12.));
555 }
556 
558  ChamberDigiMap& chDigiMap, unsigned int bx, unsigned int part, unsigned int strip, unsigned int currentIDX) const {
559  DigiIndicies newIDX(bx, part, strip);
560  auto it1 = chDigiMap.find(newIDX);
561  if (it1 == chDigiMap.end()) {
562  chDigiMap[newIDX] = currentIDX;
563  return -1;
564  }
565  return it1->second;
566 }
567 
569 #endif
void getStripProperties(const ME0EtaPartition *etaPart, const ME0DigiPreReco *inDigi, float &tof, int &strip, LocalPoint &digiLocalPoint, LocalError &digiLocalError) const
edm::ESGetToken< ME0Geometry, MuonGeometryRecord > geom_token_
std::vector< std::vector< double > > tofs
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
float x() const
class Point3DBase< float, LocalTag > LocalPoint
static unsigned int partIdx(const InputGenJetsParticleSelector::ParticleVector &p, const reco::Candidate *particle)
const ME0Geometry * geometry
Definition: APVGainStruct.h:7
TrapezoidalStripTopology * buildTopo(const std::vector< float > &_p) const
edm::EDGetTokenT< ME0DigiPreRecoCollection > token
void buildDigis(const ME0DigiPreRecoCollection &, ME0DigiPreRecoCollection &, ME0DigiPreRecoMap &, CLHEP::HepRandomEngine *engine)
T eta() const
Definition: PV3DBase.h:73
T const * product() const
Definition: Handle.h:70
std::vector< int > layerReadout
~ME0ReDigiProducer() override
void insertDigi(const IndexType &index, const DigiType &digi)
insert a digi for a given DetUnit
std::tuple< unsigned int, unsigned int, unsigned int > DigiIndicies
TemporaryGeometry(const ME0Geometry *geometry, const unsigned int numberOfStrips, const unsigned int numberOfPartitions)
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
float getPartCenter(const unsigned int partIdx) const
float y() const
#define LogTrace(id)
float yy() const
Definition: LocalError.h:24
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
int layer() const
Layer id: each chamber has six layers of chambers: layer 1 is the inner layer and layer 6 is the oute...
Definition: ME0DetId.h:44
MuonDigiCollection< ME0DetId, ME0DigiPreReco > ME0DigiPreRecoCollection
T sqrt(T t)
Definition: SSEVec.h:19
float getCentralTOF(const ME0DetId &me0Id, unsigned int partIdx) const
const TrapezoidalStripTopology * getTopo(const unsigned int partIdx) const
Transition
Definition: Transition.h:12
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ME0DetId id() const
int getCustomStripProperties(const ME0DetId &detId, const ME0DigiPreReco *inDigi, float &tof, int &strip, LocalPoint &digiLocalPoint, LocalError &digiLocalError) const
std::map< DigiIndicies, unsigned int > ChamberDigiMap
float xy() const
Definition: LocalError.h:23
TemporaryGeometry * tempGeo
unsigned int fillDigiMap(ChamberDigiMap &chDigiMap, unsigned int bx, unsigned int part, unsigned int strip, unsigned int currentIDX) const
part
Definition: HCALResponse.h:20
unsigned int numberOfPartitions
unsigned int numberOfStrips
std::pair< const_iterator, const_iterator > Range
double b
Definition: hdecay.h:120
std::vector< DigiType >::const_iterator const_iterator
unsigned int findEtaPartition(float locY) const
MuonDigiCollection< ME0DetId, int > ME0DigiPreRecoMap
std::vector< TrapezoidalStripTopology * > stripTopos
ME0ReDigiProducer(const edm::ParameterSet &ps)
HLT enums.
const StripTopology & specificTopology() const
int roll() const
Definition: ME0DetId.h:48
std::vector< std::vector< double > > tofs
strips
#turn off noise in all subdetectors simHcalUnsuppressedDigis.doNoise = False mix.digitizers.hcal.doNoise = False simEcalUnsuppressedDigis.doNoise = False mix.digitizers.ecal.doNoise = False simEcalUnsuppressedDigis.doESNoise = False simSiPixelDigis.AddNoise = False mix.digitizers.pixel.AddNoise = False simSiStripDigis.Noise = False mix.digitizers.strip.AddNoise = False
Definition: DigiDM_cff.py:32
float stripLength() const override
det heigth (strip length in the middle)
static char chambers[264][20]
Definition: ReadPGInfo.cc:243
The Signals That Services Can Subscribe To This is based on ActivityRegistry h
Helper function to determine trigger accepts.
Definition: Activities.doc:4
void produce(edm::Event &, const edm::EventSetup &) override
float xx() const
Definition: LocalError.h:22
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
void beginRun(const edm::Run &, const edm::EventSetup &) override
#define LogDebug(id)
A container for a generic type of digis indexed by some index, implemented with a map<IndexType...