CMS 3D CMS Logo

CSCStripElectronicsSim.cc
Go to the documentation of this file.
12 
13 #include "CLHEP/Random/RandGaussQ.h"
14 #include "CLHEP/Units/GlobalPhysicalConstants.h"
15 #include "CLHEP/Units/GlobalSystemOfUnits.h"
16 
17 #include "boost/bind.hpp"
18 #include <cassert>
19 #include <list>
20 
21 // This is CSCStripElectronicsSim.cc
22 
25  theAmpResponse(theShapingTime, CSCStripAmpResponse::RADICAL),
26  theComparatorThreshold(20.),
27  theComparatorNoise(0.),
28  theComparatorRMSOffset(2.),
29  theComparatorSaturation(1057.),
30  theComparatorWait(50.),
31  theComparatorDeadTime(100.),
32  theDaqDeadTime(200.),
33  theTimingOffset(0.),
34  nScaBins_(p.getParameter<int>("nScaBins")),
35  doSuppression_(p.getParameter<bool>("doSuppression")),
36  doCrosstalk_(p.getParameter<bool>("doCrosstalk")),
37  theStripConditions(nullptr),
38  theCrosstalkGenerator(nullptr),
39  theComparatorClockJump(2),
40  sca_time_bin_size(50.),
41  sca_peak_bin(p.getParameter<int>("scaPeakBin")),
42  theComparatorTimeBinOffset(p.getParameter<double>("comparatorTimeBinOffset")),
43  theComparatorTimeOffset(p.getParameter<double>("comparatorTimeOffset")),
44  theComparatorSamplingTime(p.getParameter<double>("comparatorSamplingTime")),
45  theSCATimingOffsets(p.getParameter<std::vector<double>>("scaTimingOffsets")) {
46  if (doCrosstalk_) {
48  }
49 
51 }
52 
54  if (doCrosstalk_) {
55  delete theCrosstalkGenerator;
56  }
57 }
58 
62  // selfTest();
63 
64  // calculate the offset to the peak
65  float averageDistance = theLayer->surface().position().mag();
66  theAverageTimeOfFlight = averageDistance * cm / c_light; // Units of c_light: mm/ns
67  int chamberType = theSpecs->chamberType();
69  // TODO make sure config gets overridden
72  theNumberOfSamples = nScaBins_ * static_cast<int>(sca_time_bin_size / theSamplingTime);
73 }
74 
76 
78 
79 CSCAnalogSignal CSCStripElectronicsSim::makeNoiseSignal(int element, CLHEP::HepRandomEngine *engine) {
80  std::vector<float> noiseBins(nScaBins_);
81  CSCAnalogSignal tmpSignal(element, sca_time_bin_size, noiseBins);
82  if (doNoise_) {
83  theStripConditions->noisify(layerId(), tmpSignal, engine);
84  }
85  // now rebin it
86  std::vector<float> binValues(theNumberOfSamples);
87  for (int ibin = 0; ibin < theNumberOfSamples; ++ibin) {
88  binValues[ibin] = tmpSignal.getValue(ibin * theSamplingTime);
89  }
90  CSCAnalogSignal finalSignal(element, theSamplingTime, binValues, 0., theSignalStartTime);
91  return finalSignal;
92 }
93 
95  float time,
96  CLHEP::HepRandomEngine *engine) const {
97  return std::min(signal.getValue(time), theComparatorSaturation) +
98  theComparatorRMSOffset * CLHEP::RandGaussQ::shoot(engine);
99 }
100 
101 void CSCStripElectronicsSim::runComparator(std::vector<CSCComparatorDigi> &result, CLHEP::HepRandomEngine *engine) {
102  // first, make a list of all the comparators we actually
103  // need to run
104  std::list<int> comparatorsWithSignal;
105  CSCSignalMap::iterator signalMapItr;
106  for (signalMapItr = theSignalMap.begin(); signalMapItr != theSignalMap.end(); ++signalMapItr) {
107  // Elements in signal map count from 1
108  // 1,2->0, 3,4->1, 5,6->2, ...
109  comparatorsWithSignal.push_back(((*signalMapItr).first - 1) / 2);
110  }
111  // no need to sort
112  comparatorsWithSignal.unique();
113  for (std::list<int>::iterator listItr = comparatorsWithSignal.begin(); listItr != comparatorsWithSignal.end();
114  ++listItr) {
115  int iComparator = *listItr;
116  // find signal1 and signal2
117  // iComparator counts from 0
118  // icomp =0->1,2, =1->3,4, =2->5,6, ...
119  const CSCAnalogSignal &signal1 = find(readoutElement(iComparator * 2 + 1), engine);
120  const CSCAnalogSignal &signal2 = find(readoutElement(iComparator * 2 + 2), engine);
123  if (comparatorReading(signal1, time, engine) > theComparatorThreshold ||
124  comparatorReading(signal2, time, engine) > theComparatorThreshold) {
125  // wait a bit, so we can run the comparator at the signal peak
126  float comparatorTime = time;
128 
129  float height1 = comparatorReading(signal1, time, engine);
130  float height2 = comparatorReading(signal2, time, engine);
131  int output = 0;
132  int strip = 0;
133  // distrip logic; comparator output is for pairs of strips:
134  // hit bin dec
135  // x--- 100 4
136  // -x-- 101 5
137  // --x- 110 6
138  // ---x 111 7
139  // just to prevent a copy
140  const CSCAnalogSignal *mainSignal = nullptr;
141  // pick the higher of the two strips in the pair
142  if (height1 > height2) {
143  mainSignal = &signal1;
144  float leftStrip = 0.;
145  if (iComparator > 0) {
146  leftStrip = comparatorReading(find(readoutElement(iComparator * 2), engine), time, engine);
147  }
148  // if this strip is higher than either of its neighbors, make a
149  // comparator digi
150  if (leftStrip < height1 && height1 > theComparatorThreshold) {
151  output = (leftStrip < height2);
152  strip = iComparator * 2 + 1;
153  }
154  } else {
155  mainSignal = &signal2;
156  float rightStrip = 0.;
157  if (iComparator * 2 + 3 <= nElements) {
158  rightStrip = comparatorReading(find(readoutElement(iComparator * 2 + 3), engine), time, engine);
159  }
160  if (rightStrip < height2 && height2 > theComparatorThreshold) {
161  output = (height1 < rightStrip);
162  strip = iComparator * 2 + 2;
163  }
164  }
165  if (strip != 0) {
166  float bxFloat =
168 
169  // Comparator digi as of Nov-2006 adapted to real data: time word has
170  // 16 bits with set bit flagging appropriate bunch crossing, and bx 0
171  // corresponding to 9th bit i.e.
172 
173  // 1st bit set (bit 0) <-> bx -9
174  // 2nd 1 <-> bx -8
175  // ... ... ....
176  // 8th 9 <-> bx 0
177  // 9th 10 <-> bx +1
178  // ... ... ....
179  // 16th 15 <-> bx +6
180 
181  // Parameter theOffsetOfBxZero = 9 @@WARNING! This offset may be
182  // changed (hardware)!
183 
184  int timeWord = 0; // and this will remain if too early or late
185  if ((bxFloat >= 0) && (bxFloat < 16))
186  timeWord = (1 << static_cast<int>(bxFloat)); // set appropriate bit
187 
188  CSCComparatorDigi newDigi(strip, output, timeWord);
189  result.push_back(newDigi);
190  }
191 
192  // wait for the comparator to reset
194  // really should be zero, but strip signal doesn't go negative yet
195  float resetThreshold = 1;
196  while (time < theSignalStopTime && mainSignal->getValue(time) > resetThreshold) {
198  }
199 
200  } // if over threshold
201  } // loop over time samples
202  } // loop over comparators
203  // sort by time
204  sort(result.begin(), result.end());
205 }
206 
207 std::list<int> CSCStripElectronicsSim::getKeyStrips(const std::vector<CSCComparatorDigi> &comparators) const {
208  std::list<int> result;
209  for (std::vector<CSCComparatorDigi>::const_iterator compItr = comparators.begin(); compItr != comparators.end();
210  ++compItr) {
211  if (std::abs(compItr->getTimeBin() - theOffsetOfBxZero) <= 2) {
212  result.push_back(compItr->getStrip());
213  }
214  }
215  // need sort for unique to work.
216  result.sort();
217  result.unique();
218  return result;
219 }
220 
222  // assumes the detector hit map is filled
223  std::list<int> result;
225  theDetectorHitMap.end(),
226  // back_inserter(result),
227  // boost::bind(&DetectorHitMap::value_type::first,_1));
228  // suggested code from Chris Jones
229  back_inserter(result),
230  std::bind(&DetectorHitMap::value_type::first, std::placeholders::_1));
231  // back_inserter(result), [](DetectorHitMap::value_type const& iValue) {
232  // return iValue.first; } );
233  result.sort();
234  result.unique();
235  return result;
236 }
237 
238 std::list<int> CSCStripElectronicsSim::channelsToRead(const std::list<int> &keyStrips, int window) const {
239  std::list<int> result;
240  std::list<int>::const_iterator keyStripItr = keyStrips.begin();
241  if (doSuppression_) {
242  for (; keyStripItr != keyStrips.end(); ++keyStripItr) {
243  // pick the five strips around the comparator
244  for (int istrip = (*keyStripItr) - window; istrip <= (*keyStripItr) + window; ++istrip) {
245  if (istrip > 0 && istrip <= nElements) {
246  result.push_back(readoutElement(istrip));
247  }
248  }
249  }
250  result.sort();
251  result.unique();
252  } else {
253  // read the whole CFEB, 16 strips
254  std::list<int> cfebsToRead;
255  for (; keyStripItr != keyStrips.end(); ++keyStripItr) {
256  int cfeb = (readoutElement(*keyStripItr) - 1) / 16;
257  cfebsToRead.push_back(cfeb);
258  int remainder = (readoutElement(*keyStripItr) - 1) % 16;
259  // if we're within 3 strips of an edge, take neighboring CFEB, too
260  if (remainder < window && cfeb != 0) {
261  cfebsToRead.push_back(cfeb - 1);
262  }
263  // the 'readouElement' makes it so that ME1/1 has just one CFEB
264  int maxCFEBs = readoutElement(nElements) / 16 - 1;
265  if (remainder >= 16 - window && cfeb != maxCFEBs) {
266  cfebsToRead.push_back(cfeb + 1);
267  }
268  }
269  cfebsToRead.sort();
270  cfebsToRead.unique();
271 
272  // now convert the CFEBS to strips
273  for (std::list<int>::const_iterator cfebItr = cfebsToRead.begin(); cfebItr != cfebsToRead.end(); ++cfebItr) {
274  for (int i = 1; i <= 16; ++i) {
275  result.push_back((*cfebItr) * 16 + i);
276  }
277  }
278  }
279  return result;
280 }
281 
283  return (s1.getTotal() > s2.getTotal());
284 }
285 
287  CSCComparatorDigiCollection &comparators,
288  CLHEP::HepRandomEngine *engine) {
289  if (doCrosstalk_) {
290  addCrosstalk(engine);
291  }
292 
293  std::vector<CSCComparatorDigi> comparatorOutputs;
294  runComparator(comparatorOutputs, engine);
295  // copy these to the result
296  if (!comparatorOutputs.empty()) {
297  CSCComparatorDigiCollection::Range range(comparatorOutputs.begin(), comparatorOutputs.end());
298  comparators.put(range, layerId());
299  }
300 
301  // std::list<int> keyStrips = getKeyStrips(comparatorOutputs);
302  std::list<int> keyStrips = getKeyStripsFromMC();
303  fillStripDigis(keyStrips, digis, engine);
304 }
305 
306 void CSCStripElectronicsSim::fillStripDigis(const std::list<int> &keyStrips,
307  CSCStripDigiCollection &digis,
308  CLHEP::HepRandomEngine *engine) {
309  std::list<int> stripsToDo = channelsToRead(keyStrips, 3);
310  std::vector<CSCStripDigi> stripDigis;
311  stripDigis.reserve(stripsToDo.size());
312  for (std::list<int>::const_iterator stripItr = stripsToDo.begin(); stripItr != stripsToDo.end(); ++stripItr) {
313  createDigi(*stripItr, find(*stripItr, engine), stripDigis, engine);
314  }
315 
316  CSCStripDigiCollection::Range stripRange(stripDigis.begin(), stripDigis.end());
317  digis.put(stripRange, layerId());
318 }
319 
320 void CSCStripElectronicsSim::addCrosstalk(CLHEP::HepRandomEngine *engine) {
321  // this is needed so we can add a noise signal to the map
322  // without messing up any iterators
323  std::vector<CSCAnalogSignal> realSignals;
324  realSignals.reserve(theSignalMap.size());
325  CSCSignalMap::iterator mapI = theSignalMap.begin(), mapEnd = theSignalMap.end();
326  for (; mapI != mapEnd; ++mapI) {
327  realSignals.push_back((*mapI).second);
328  }
329  sort(realSignals.begin(), realSignals.end(), SortSignalsByTotal);
330  std::vector<CSCAnalogSignal>::iterator realSignalItr = realSignals.begin(), realSignalsEnd = realSignals.end();
331  for (; realSignalItr != realSignalsEnd; ++realSignalItr) {
332  int thisStrip = (*realSignalItr).getElement();
333  // add it to each neighbor
334  if (thisStrip > 1) {
335  int otherStrip = thisStrip - 1;
336  addCrosstalk(*realSignalItr, thisStrip, otherStrip, engine);
337  }
338  if (thisStrip < nElements) {
339  int otherStrip = thisStrip + 1;
340  addCrosstalk(*realSignalItr, thisStrip, otherStrip, engine);
341  }
342  }
343 }
344 
346  int thisStrip,
347  int otherStrip,
348  CLHEP::HepRandomEngine *engine) {
349  float capacitiveCrosstalk, resistiveCrosstalk;
350  bool leftRight = (otherStrip > thisStrip);
352  layerId(), thisStrip, theLayerGeometry->length(), leftRight, capacitiveCrosstalk, resistiveCrosstalk);
353  theCrosstalkGenerator->setParameters(capacitiveCrosstalk, 0., resistiveCrosstalk);
354  CSCAnalogSignal crosstalkSignal(theCrosstalkGenerator->getCrosstalk(signal));
355  find(readoutElement(otherStrip), engine).superimpose(crosstalkSignal);
356 
357  // Now subtract the crosstalk signal from the original signal
358  crosstalkSignal *= -1.;
359  find(thisStrip, engine).superimpose(crosstalkSignal);
360 }
361 
363  const CSCAnalogSignal &signal,
364  std::vector<CSCStripDigi> &result,
365  CLHEP::HepRandomEngine *engine) {
366  // fill in the sca information
367  std::vector<int> scaCounts(nScaBins_);
368 
369  float pedestal = theStripConditions->pedestal(layerId(), channel);
370  float gain = theStripConditions->smearedGain(layerId(), channel, engine);
371  int chamberType = theSpecs->chamberType();
372  float timeSmearing = CLHEP::RandGaussQ::shoot(engine) * theTimingCalibrationError[chamberType];
373  // undo the correction for TOF, instead, using some nominal
374  // value from ME2/1
375  float t0 = theSignalStartTime + theSCATimingOffsets[chamberType] + timeSmearing + 29. - theAverageTimeOfFlight;
376  for (int scaBin = 0; scaBin < nScaBins_; ++scaBin) {
377  float t = t0 + scaBin * sca_time_bin_size;
378  scaCounts[scaBin] = static_cast<int>(pedestal + signal.getValue(t) * gain);
379  }
380  CSCStripDigi newDigi(channel, scaCounts);
381 
382  // do saturation of 12-bit ADC
383  doSaturation(newDigi);
384 
385  result.push_back(newDigi);
386  addLinks(channelIndex(channel));
387  LogTrace("CSCStripElectronicsSim") << newDigi;
388 }
389 
391  std::vector<int> scaCounts(digi.getADCCounts());
392  for (unsigned scaBin = 0; scaBin < scaCounts.size(); ++scaBin) {
393  scaCounts[scaBin] = std::min(scaCounts[scaBin], 4095);
394  }
395  digi.setADCCounts(scaCounts);
396 }
397 
399  const CSCComparatorDigiCollection &comparators,
400  CSCStripDigiCollection &digis,
401  CLHEP::HepRandomEngine *engine) {
402  theSignalMap.clear();
403  setLayer(layer);
404  CSCDetId chamberId(theLayerId.chamberId());
405  // find all comparator key strips in this chamber
406  std::list<int> chamberKeyStrips;
407  for (CSCComparatorDigiCollection::DigiRangeIterator comparatorItr = comparators.begin();
408  comparatorItr != comparators.end();
409  ++comparatorItr) {
410  // could be more efficient
411  if (CSCDetId((*comparatorItr).first).chamberId() == chamberId) {
412  std::vector<CSCComparatorDigi> layerComparators((*comparatorItr).second.first, (*comparatorItr).second.second);
413  std::list<int> layerKeyStrips = getKeyStrips(layerComparators);
414  chamberKeyStrips.insert(chamberKeyStrips.end(), layerKeyStrips.begin(), layerKeyStrips.end());
415  }
416  }
417  chamberKeyStrips.sort();
418  chamberKeyStrips.unique();
419  fillStripDigis(chamberKeyStrips, digis, engine);
420 }
421 
423  // make sure the zero suppression algorithms work
424  std::list<int> keyStrips, stripsRead;
425  //
426  bool isGanged = (readoutElement(nElements) == 16);
427  keyStrips.push_back(readoutElement(19));
428  keyStrips.push_back(readoutElement(30));
429  keyStrips.push_back(readoutElement(32));
430  stripsRead = channelsToRead(keyStrips, 3);
431  if (doSuppression_) {
432  unsigned int expectedSize = isGanged ? 10 : 12;
433  assert(stripsRead.size() == expectedSize);
434  assert(stripsRead.front() == readoutElement(17));
435  } else {
436  unsigned int expectedSize = isGanged ? 16 : 48;
437  assert(stripsRead.size() == expectedSize);
438  assert(stripsRead.front() == 1);
439  }
440 }
std::vector< double > theBunchTimingOffsets
int readoutElement(int strip) const override
void initParameters() override
initialization for each layer
const CSCChamberSpecs * theSpecs
const CSCLayerGeometry * theLayerGeometry
bool SortSignalsByTotal(const CSCAnalogSignal &s1, const CSCAnalogSignal &s2)
std::vector< int > const & getADCCounts() const
Get ADC readings.
Definition: CSCStripDigi.h:44
std::list< int > getKeyStripsFromMC() const
get ths strips that have detector hits
void addCrosstalk(CLHEP::HepRandomEngine *)
float length() const override
#define nullptr
void runComparator(std::vector< CSCComparatorDigi > &result, CLHEP::HepRandomEngine *)
int numberOfStrips() const
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
void setADCCounts(const std::vector< int > &ADCCounts)
Definition: CSCStripDigi.cc:22
CSCStripConditions * theStripConditions
void setLayer(const CSCLayer *layer)
float comparatorReading(const CSCAnalogSignal &signal, float time, CLHEP::HepRandomEngine *) const
calculates the comparator reading, including saturation and offsets
std::vector< double > theSCATimingOffsets
T mag() const
Definition: PV3DBase.h:67
std::list< int > channelsToRead(const std::list< int > &keyStrips, int window) const
float calculateAmpResponse(float t) const override
void noisify(const CSCDetId &detId, CSCAnalogSignal &signal, CLHEP::HepRandomEngine *)
superimposes noise, in fC, on the signal
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CSCAnalogSignal makeNoiseSignal(int element, CLHEP::HepRandomEngine *) override
CSCDetId chamberId() const
Definition: CSCDetId.h:53
def window(xmin, xmax, ymin, ymax, x=0, y=0, width=100, height=100, xlogbase=None, ylogbase=None, minusInfinity=-1000, flipx=False, flipy=True)
Definition: svgfig.py:643
void fillMissingLayer(const CSCLayer *layer, const CSCComparatorDigiCollection &comparators, CSCStripDigiCollection &digis, CLHEP::HepRandomEngine *)
CSCCrosstalkGenerator * theCrosstalkGenerator
T min(T a, T b)
Definition: MathUtil.h:58
void superimpose(const CSCAnalogSignal &signal2)
void doSaturation(CSCStripDigi &digi)
#define LogTrace(id)
int channel(int strip) const
float calculateAmpResponse(float t) const
virtual void addLinks(int channelIndex)
CSCStripElectronicsSim(const edm::ParameterSet &p)
configurable parameters
void setParameters(float crosstalk, float delay, float resistiveFraction)
int chamberType() const
std::list< int > getKeyStrips(const std::vector< CSCComparatorDigi > &comparators) const
finds the key strips from these comparators
float getTotal() const
CSCAnalogSignal getCrosstalk(const CSCAnalogSignal &inputSignal) const
virtual float pedestal(const CSCDetId &detId, int channel) const =0
in ADC counts
CSCDetId layerId() const
the CSCDetId corresponding to the current layer
CSCAnalogSignal & find(int element, CLHEP::HepRandomEngine *)
void fillDigis(CSCStripDigiCollection &digis, CSCComparatorDigiCollection &comparators, CLHEP::HepRandomEngine *)
virtual void crosstalk(const CSCDetId &detId, int channel, double stripLength, bool leftRight, float &capacitive, float &resistive) const =0
float getValue(float t) const
virtual float smearedGain(const CSCDetId &detId, int channel, CLHEP::HepRandomEngine *) const
virtual int channelIndex(int channel) const
lets users map channels to different indices for links
const JetExtendedData & getValue(const Container &, const reco::JetBaseRef &)
get value for the association. Throw exception if no association found
std::pair< const_iterator, const_iterator > Range
const PositionType & position() const
CSCStripAmpResponse theAmpResponse
void createDigi(int istrip, const CSCAnalogSignal &signal, std::vector< CSCStripDigi > &result, CLHEP::HepRandomEngine *)
void fillStripDigis(const std::list< int > &keyStrips, CSCStripDigiCollection &digis, CLHEP::HepRandomEngine *)
std::vector< double > theTimingCalibrationError