CMS 3D CMS Logo

DigiSimLinkAlgorithm.cc
Go to the documentation of this file.
1 // File: DigiSimLinkAlgorithm.cc
2 // Description: Steering class for digitization.
3 
4 #include <vector>
5 #include <algorithm>
6 #include <iostream>
9 #include "DigiSimLinkAlgorithm.h"
14 
15 #include "CLHEP/Random/RandFlat.h"
16 
17 #define CBOLTZ (1.38E-23)
18 #define e_SI (1.6E-19)
19 
21  theThreshold = conf_.getParameter<double>("NoiseSigmaThreshold");
22  theFedAlgo = conf_.getParameter<int>("FedAlgorithm");
23  peakMode = conf_.getParameter<bool>("APVpeakmode");
24  theElectronPerADC = conf_.getParameter<double>(peakMode ? "electronPerAdcPeak" : "electronPerAdcDec");
25  noise = conf_.getParameter<bool>("Noise");
26  zeroSuppression = conf_.getParameter<bool>("ZeroSuppression");
27  theTOFCutForPeak = conf_.getParameter<double>("TOFCutForPeak");
28  theTOFCutForDeconvolution = conf_.getParameter<double>("TOFCutForDeconvolution");
29  cosmicShift = conf_.getUntrackedParameter<double>("CosmicDelayShift");
30  inefficiency = conf_.getParameter<double>("Inefficiency");
31  RealPedestals = conf_.getParameter<bool>("RealPedestals");
32  SingleStripNoise = conf_.getParameter<bool>("SingleStripNoise");
33  CommonModeNoise = conf_.getParameter<bool>("CommonModeNoise");
34  BaselineShift = conf_.getParameter<bool>("BaselineShift");
35  APVSaturationFromHIP = conf_.getParameter<bool>("APVSaturationFromHIP");
36  APVSaturationProb = conf_.getParameter<double>("APVSaturationProb");
37  cmnRMStib = conf_.getParameter<double>("cmnRMStib");
38  cmnRMStob = conf_.getParameter<double>("cmnRMStob");
39  cmnRMStid = conf_.getParameter<double>("cmnRMStid");
40  cmnRMStec = conf_.getParameter<double>("cmnRMStec");
41  pedOffset = (unsigned int)conf_.getParameter<double>("PedestalsOffset");
42  PreMixing_ = conf_.getParameter<bool>("PreMixingMode");
43  if (peakMode) {
45  LogDebug("StripDigiInfo") << "APVs running in peak mode (poor time resolution)";
46  } else {
48  LogDebug("StripDigiInfo") << "APVs running in deconvolution mode (good time resolution)";
49  };
50 
56 }
57 
59  delete theSiHitDigitizer;
61  delete theSiNoiseAdder;
62  delete theSiDigitalConverter;
63  delete theSiZeroSuppress;
64  //delete particle;
65  //delete pdt;
66 }
67 
68 // Run the algorithm for a given module
69 // ------------------------------------
70 
72  edm::DetSet<SiStripRawDigi>& outrawdigi,
73  const std::vector<std::pair<const PSimHit*, int> >& input,
74  StripGeomDetUnit const* det,
75  GlobalVector bfield,
76  float langle,
77  edm::ESHandle<SiStripGain>& gainHandle,
78  edm::ESHandle<SiStripThreshold>& thresholdHandle,
79  edm::ESHandle<SiStripNoises>& noiseHandle,
80  edm::ESHandle<SiStripPedestals>& pedestalHandle,
81  edm::ESHandle<SiStripBadStrip>& deadChannelHandle,
82  const TrackerTopology* tTopo,
83  CLHEP::HepRandomEngine* engine) {
85  unsigned int detID = det->geographicalId().rawId();
86  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
87  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
88  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
89  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detID);
90  numStrips = (det->specificTopology()).nstrips();
91 
92  //stroing the bad stip of the the module. the module is not removed but just signal put to 0
93  std::vector<bool> badChannels;
94  badChannels.clear();
95  badChannels.insert(badChannels.begin(), numStrips, false);
97  for (SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
98  fs = deadChannelHandle->decode(*it);
99  for (int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip)
100  badChannels[strip] = true;
101  }
102 
103  // local amplitude of detector channels (from processed PSimHit)
104  // locAmpl.clear();
105  detAmpl.clear();
106  // locAmpl.insert(locAmpl.begin(),numStrips,0.);
107  // total amplitude of detector channels
108  detAmpl.insert(detAmpl.begin(), numStrips, 0.);
109 
112 
113  // First: loop on the SimHits
114  std::vector<std::pair<const PSimHit*, int> >::const_iterator simHitIter = input.begin();
115  std::vector<std::pair<const PSimHit*, int> >::const_iterator simHitIterEnd = input.end();
116  if (CLHEP::RandFlat::shoot(engine) > inefficiency) {
117  for (; simHitIter != simHitIterEnd; ++simHitIter) {
118  locAmpl.clear();
119  locAmpl.insert(locAmpl.begin(), numStrips, 0.);
120  // check TOF
121  if (std::fabs(((*simHitIter).first)->tof() - cosmicShift -
122  det->surface().toGlobal(((*simHitIter).first)->localPosition()).mag() / 30.) < tofCut &&
123  ((*simHitIter).first)->energyLoss() > 0) {
125  localLastChannel = 0;
126  // process the hit
128  ((*simHitIter).first), *det, bfield, langle, locAmpl, localFirstChannel, localLastChannel, tTopo, engine);
129 
130  //APV Killer to simulate HIP effect
131  //------------------------------------------------------
132 
134  int pdg_id = ((*simHitIter).first)->particleType();
135  particle = pdt->particle(pdg_id);
136  if (particle != nullptr) {
137  float charge = particle->charge();
138  bool isHadron = particle->isHadron();
139  if (charge != 0 && isHadron) {
140  if (CLHEP::RandFlat::shoot(engine) < APVSaturationProb) {
141  int FirstAPV = localFirstChannel / 128;
142  int LastAPV = localLastChannel / 128;
143  std::cout << "-------------------HIP--------------" << std::endl;
144  std::cout << "Killing APVs " << FirstAPV << " - " << LastAPV << " " << detID << std::endl;
145  for (int strip = FirstAPV * 128; strip < LastAPV * 128 + 128; ++strip)
146  badChannels[strip] = true; //doing like that I remove the signal information only after the
147  //stip that got the HIP but it remains the signal of the previous
148  //one. I'll make a further loop to remove all signal
149  }
150  }
151  }
152  }
153 
155  locAmpl, localFirstChannel, localLastChannel, ((*simHitIter).first), (*simHitIter).second);
156 
157  // sum signal on strips
158  for (size_t iChannel = localFirstChannel; iChannel < localLastChannel; iChannel++) {
159  if (locAmpl[iChannel] > 0.) {
160  //if(!badChannels[iChannel]) detAmpl[iChannel]+=locAmpl[iChannel];
161  //locAmpl[iChannel]=0;
162  detAmpl[iChannel] += locAmpl[iChannel];
163  }
164  }
167  if (lastChannelWithSignal < localLastChannel)
169  }
170  }
171  }
172 
173  //removing signal from the dead (and HIP effected) strips
174  for (int strip = 0; strip < numStrips; ++strip)
175  if (badChannels[strip])
176  detAmpl[strip] = 0;
177 
181 
182  if (zeroSuppression) {
183  if (noise) {
184  int RefStrip = int(numStrips / 2.);
185  while (RefStrip < numStrips &&
186  badChannels[RefStrip]) { //if the refstrip is bad, I move up to when I don't find it
187  RefStrip++;
188  }
189  if (RefStrip < numStrips) {
190  float noiseRMS = noiseHandle->getNoise(RefStrip, detNoiseRange);
191  float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
195  numStrips,
196  noiseRMS * theElectronPerADC / gainValue,
197  engine);
198  }
199  }
200  digis.clear();
202  theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID, noiseHandle, thresholdHandle);
203  push_link(digis, theLink, theCounterLink, detAmpl, detID);
204  outdigi.data = digis;
205  }
206 
207  if (!zeroSuppression) {
208  //if(noise){
209  // the constant pedestal offset is needed because
210  // negative adc counts are not allowed in case
211  // Pedestal and CMN subtraction is performed.
212  // The pedestal value read from the conditions
213  // is pedValue and after the pedestal subtraction
214  // the baseline is zero. The Common Mode Noise
215  // is not subtracted from the negative adc counts
216  // channels. Adding pedOffset the baseline is set
217  // to pedOffset after pedestal subtraction and CMN
218  // is subtracted to all the channels since none of
219  // them has negative adc value. The pedOffset is
220  // treated as a constant component in the CMN
221  // estimation and subtracted as CMN.
222 
223  //calculating the charge deposited on each APV and subtracting the shift
224  //------------------------------------------------------
225  if (BaselineShift) {
227  }
228 
229  //Adding the strip noise
230  //------------------------------------------------------
231  if (noise) {
232  std::vector<float> noiseRMSv;
233  noiseRMSv.clear();
234  noiseRMSv.insert(noiseRMSv.begin(), numStrips, 0.);
235 
236  if (SingleStripNoise) {
237  for (int strip = 0; strip < numStrips; ++strip) {
238  if (!badChannels[strip])
239  noiseRMSv[strip] = (noiseHandle->getNoise(strip, detNoiseRange)) * theElectronPerADC;
240  }
241 
242  } else {
243  int RefStrip = 0; //int(numStrips/2.);
244  while (RefStrip < numStrips &&
245  badChannels[RefStrip]) { //if the refstrip is bad, I move up to when I don't find it
246  RefStrip++;
247  }
248  if (RefStrip < numStrips) {
249  float noiseRMS = noiseHandle->getNoise(RefStrip, detNoiseRange) * theElectronPerADC;
250  for (int strip = 0; strip < numStrips; ++strip) {
251  if (!badChannels[strip])
252  noiseRMSv[strip] = noiseRMS;
253  }
254  }
255  }
256 
257  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
258  }
259 
260  //adding the CMN
261  //------------------------------------------------------
262  if (CommonModeNoise) {
263  float cmnRMS = 0.;
264  DetId detId(detID);
265  uint32_t SubDet = detId.subdetId();
266  if (SubDet == 3) {
267  cmnRMS = cmnRMStib;
268  } else if (SubDet == 4) {
269  cmnRMS = cmnRMStid;
270  } else if (SubDet == 5) {
271  cmnRMS = cmnRMStob;
272  } else if (SubDet == 6) {
273  cmnRMS = cmnRMStec;
274  }
275  cmnRMS *= theElectronPerADC;
276  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels, engine);
277  }
278 
279  //Adding the pedestals
280  //------------------------------------------------------
281 
282  std::vector<float> vPeds;
283  vPeds.clear();
284  vPeds.insert(vPeds.begin(), numStrips, 0.);
285 
286  if (RealPedestals) {
287  for (int strip = 0; strip < numStrips; ++strip) {
288  if (!badChannels[strip])
289  vPeds[strip] = (pedestalHandle->getPed(strip, detPedestalRange) + pedOffset) * theElectronPerADC;
290  }
291  } else {
292  for (int strip = 0; strip < numStrips; ++strip) {
293  if (!badChannels[strip])
294  vPeds[strip] = pedOffset * theElectronPerADC;
295  }
296  }
297 
299 
300  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
301  // edm::LogWarning("DigiSimLinkAlgorithm")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
302  //}else{
303 
304  rawdigis.clear();
305  rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
306  push_link_raw(rawdigis, theLink, theCounterLink, detAmpl, detID);
307  outrawdigi.data = rawdigis;
308 
309  //}
310  }
311 }
312 
314  const HitToDigisMapType& htd,
315  const HitCounterToDigisMapType& hctd,
316  const std::vector<float>& afterNoise,
317  unsigned int detID) {
318  link_coll.clear();
319  for (DigitalVecType::const_iterator i = digis.begin(); i != digis.end(); i++) {
320  // Instead of checking the validity of the links against the digis,
321  // let's loop over digis and push the corresponding link
322  HitToDigisMapType::const_iterator mi(htd.find(i->strip()));
323  if (mi == htd.end())
324  continue;
325  HitCounterToDigisMapType::const_iterator cmi(hctd.find(i->strip()));
326  std::map<const PSimHit*, Amplitude> totalAmplitudePerSimHit;
327  for (std::vector<std::pair<const PSimHit*, Amplitude> >::const_iterator simul = (*mi).second.begin();
328  simul != (*mi).second.end();
329  simul++) {
330  totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
331  }
332 
333  //--- include the noise as well
334  double totalAmplitude1 = afterNoise[(*mi).first];
335 
336  //--- digisimlink
337  int sim_counter = 0;
338  for (std::map<const PSimHit*, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin();
339  iter != totalAmplitudePerSimHit.end();
340  iter++) {
341  float threshold = 0.;
342  float fraction = (*iter).second / totalAmplitude1;
343  if (fraction >= threshold) {
344  // Noise fluctuation could make fraction>1. Unphysical, set it by hand = 1.
345  if (fraction > 1.)
346  fraction = 1.;
347  for (std::vector<std::pair<const PSimHit*, int> >::const_iterator simcount = (*cmi).second.begin();
348  simcount != (*cmi).second.end();
349  ++simcount) {
350  if ((*iter).first == (*simcount).first)
351  sim_counter = (*simcount).second;
352  }
353  link_coll.push_back(StripDigiSimLink((*mi).first, //channel
354  ((*iter).first)->trackId(), //simhit trackId
355  sim_counter, //simhit counter
356  ((*iter).first)->eventId(), //simhit eventId
357  fraction)); //fraction
358  }
359  }
360  }
361 }
362 
364  const HitToDigisMapType& htd,
365  const HitCounterToDigisMapType& hctd,
366  const std::vector<float>& afterNoise,
367  unsigned int detID) {
368  link_coll.clear();
369  int nstrip = -1;
370  for (DigitalRawVecType::const_iterator i = digis.begin(); i != digis.end(); i++) {
371  nstrip++;
372  // Instead of checking the validity of the links against the digis,
373  // let's loop over digis and push the corresponding link
374  HitToDigisMapType::const_iterator mi(htd.find(nstrip));
375  HitCounterToDigisMapType::const_iterator cmi(hctd.find(nstrip));
376  if (mi == htd.end())
377  continue;
378  std::map<const PSimHit*, Amplitude> totalAmplitudePerSimHit;
379  for (std::vector<std::pair<const PSimHit*, Amplitude> >::const_iterator simul = (*mi).second.begin();
380  simul != (*mi).second.end();
381  simul++) {
382  totalAmplitudePerSimHit[(*simul).first] += (*simul).second;
383  }
384 
385  //--- include the noise as well
386  double totalAmplitude1 = afterNoise[(*mi).first];
387 
388  //--- digisimlink
389  int sim_counter_raw = 0;
390  for (std::map<const PSimHit*, Amplitude>::const_iterator iter = totalAmplitudePerSimHit.begin();
391  iter != totalAmplitudePerSimHit.end();
392  iter++) {
393  float threshold = 0.;
394  float fraction = (*iter).second / totalAmplitude1;
395  if (fraction >= threshold) {
396  //Noise fluctuation could make fraction>1. Unphysical, set it by hand.
397  if (fraction > 1.)
398  fraction = 1.;
399  //add counter information
400  for (std::vector<std::pair<const PSimHit*, int> >::const_iterator simcount = (*cmi).second.begin();
401  simcount != (*cmi).second.end();
402  ++simcount) {
403  if ((*iter).first == (*simcount).first)
404  sim_counter_raw = (*simcount).second;
405  }
406  link_coll.push_back(StripDigiSimLink((*mi).first, //channel
407  ((*iter).first)->trackId(), //simhit trackId
408  sim_counter_raw, //simhit counter
409  ((*iter).first)->eventId(), //simhit eventId
410  fraction)); //fraction
411  }
412  }
413  }
414 }
#define LogDebug(id)
unsigned short range
void push_link(const DigitalVecType &, const HitToDigisMapType &, const HitCounterToDigisMapType &, const std::vector< float > &, unsigned int)
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
void addNoise(std::vector< float > &, size_t &, size_t &, int, float, CLHEP::HepRandomEngine *) const override
SiHitDigitizer * theSiHitDigitizer
DigitalRawVecType rawdigis
edm::ParameterSet conf_
SiStripFedZeroSuppression * theSiZeroSuppress
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< unsigned int >::const_iterator ContainerIterator
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
std::pair< ContainerIterator, ContainerIterator > Range
void suppress(const std::vector< SiStripDigi > &in, std::vector< SiStripDigi > &selectedSignal, uint32_t detId, edm::ESHandle< SiStripNoises > &, edm::ESHandle< SiStripThreshold > &)
const ParticleDataTable * pdt
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
void push_link_raw(const DigitalRawVecType &, const HitToDigisMapType &, const HitCounterToDigisMapType &, const std::vector< float > &, unsigned int)
SiDigitalConverter::DigitalRawVecType DigitalRawVecType
SiDigitalConverter::DigitalVecType DigitalVecType
static std::string const input
Definition: EdmProvDump.cc:48
std::vector< StripDigiSimLink > link_coll
float getPed(const uint16_t &strip, const Range &range) const
void addNoiseVR(std::vector< float > &, std::vector< float > &, CLHEP::HepRandomEngine *) const override
SiTrivialDigitalConverter * theSiDigitalConverter
static float getNoise(uint16_t strip, const Range &range)
Definition: SiStripNoises.h:74
std::map< int, std::vector< std::pair< const PSimHit *, Amplitude > >, std::less< int > > HitToDigisMapType
virtual void add(const std::vector< float > &locAmpl, const size_t &firstChannelWithSignal, const size_t &lastChannelWithSignal, const PSimHit *hit, const int &counter)
void run(edm::DetSet< SiStripDigi > &, edm::DetSet< SiStripRawDigi > &, const std::vector< std::pair< const PSimHit *, int > > &, StripGeomDetUnit const *, GlobalVector, float, edm::ESHandle< SiStripGain > &, edm::ESHandle< SiStripThreshold > &, edm::ESHandle< SiStripNoises > &, edm::ESHandle< SiStripPedestals > &, edm::ESHandle< SiStripBadStrip > &, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
static float getStripGain(const uint16_t &strip, const SiStripApvGain::Range &range)
Definition: SiStripGain.h:73
std::vector< float > detAmpl
void addCMNoise(std::vector< float > &, float, std::vector< bool > &, CLHEP::HepRandomEngine *) const override
SiGaussianTailNoiseAdder * theSiNoiseAdder
const HitToDigisMapType & dumpLink() const
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
std::pair< ContainerIterator, ContainerIterator > Range
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
const HitCounterToDigisMapType & dumpCounterLink() const
void addPedestals(std::vector< float > &, std::vector< float > &) const override
DigitalVecType convert(const std::vector< float > &, edm::ESHandle< SiStripGain > &, unsigned int detid) override
DigiSimLinkPileUpSignals::HitCounterToDigisMapType HitCounterToDigisMapType
void processHit(const PSimHit *, const StripGeomDetUnit &, GlobalVector, float, std::vector< float > &, size_t &, size_t &, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
DigiSimLinkAlgorithm(const edm::ParameterSet &conf)
Definition: DetId.h:18
const ParticleData * particle
unsigned short firstStrip
collection_type data
Definition: DetSet.h:80
const Range getRange(const uint32_t detID) const
DigitalRawVecType convertRaw(const std::vector< float > &, edm::ESHandle< SiStripGain > &, unsigned int detid) override
const Range getRange(const uint32_t detID) const
std::pair< ContainerIterator, ContainerIterator > Range
DigiSimLinkPileUpSignals::HitToDigisMapType HitToDigisMapType
DigiSimLinkPileUpSignals * theDigiSimLinkPileUpSignals
const Range getRange(const uint32_t &detID) const
std::vector< float > locAmpl
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:50
void addBaselineShift(std::vector< float > &, std::vector< bool > &) const override
data decode(const unsigned int &value) const
std::map< int, std::vector< std::pair< const PSimHit *, int > >, std::less< int > > HitCounterToDigisMapType
const SiStripApvGain::Range getRange(uint32_t detID) const
Definition: SiStripGain.h:71