CMS 3D CMS Logo

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