CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SiStripDigitizerAlgorithm.cc
Go to the documentation of this file.
1 // File: SiStripDigitizerAlgorithm.cc
2 // Description: Steering class for digitization.
3 
4 // Modified 15/May/2013 mark.grimes@bristol.ac.uk - Modified so that the digi-sim link has the correct
5 // index for the sim hits stored. It was previously always set to zero (I won't mention that it was
6 // me who originally wrote that).
7 
8 #include <vector>
9 #include <algorithm>
10 #include <iostream>
27 #include "CLHEP/Random/RandFlat.h"
28 
30  lorentzAngleName(conf.getParameter<std::string>("LorentzAngle")),
31  theThreshold(conf.getParameter<double>("NoiseSigmaThreshold")),
32  cmnRMStib(conf.getParameter<double>("cmnRMStib")),
33  cmnRMStob(conf.getParameter<double>("cmnRMStob")),
34  cmnRMStid(conf.getParameter<double>("cmnRMStid")),
35  cmnRMStec(conf.getParameter<double>("cmnRMStec")),
36  APVSaturationProb(conf.getParameter<double>("APVSaturationProb")),
37  makeDigiSimLinks_(conf.getUntrackedParameter<bool>("makeDigiSimLinks", false)),
38  peakMode(conf.getParameter<bool>("APVpeakmode")),
39  noise(conf.getParameter<bool>("Noise")),
40  RealPedestals(conf.getParameter<bool>("RealPedestals")),
41  SingleStripNoise(conf.getParameter<bool>("SingleStripNoise")),
42  CommonModeNoise(conf.getParameter<bool>("CommonModeNoise")),
43  BaselineShift(conf.getParameter<bool>("BaselineShift")),
44  APVSaturationFromHIP(conf.getParameter<bool>("APVSaturationFromHIP")),
45  theFedAlgo(conf.getParameter<int>("FedAlgorithm")),
46  zeroSuppression(conf.getParameter<bool>("ZeroSuppression")),
47  theElectronPerADC(conf.getParameter<double>( peakMode ? "electronPerAdcPeak" : "electronPerAdcDec" )),
48  theTOFCutForPeak(conf.getParameter<double>("TOFCutForPeak")),
49  theTOFCutForDeconvolution(conf.getParameter<double>("TOFCutForDeconvolution")),
50  tofCut(peakMode ? theTOFCutForPeak : theTOFCutForDeconvolution),
51  cosmicShift(conf.getUntrackedParameter<double>("CosmicDelayShift")),
52  inefficiency(conf.getParameter<double>("Inefficiency")),
53  pedOffset((unsigned int)conf.getParameter<double>("PedestalsOffset")),
54  theSiHitDigitizer(new SiHitDigitizer(conf)),
55  theSiPileUpSignals(new SiPileUpSignals()),
56  theSiNoiseAdder(new SiGaussianTailNoiseAdder(theThreshold)),
57  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC)),
58  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo)) {
59 
60  if (peakMode) {
61  LogDebug("StripDigiInfo")<<"APVs running in peak mode (poor time resolution)";
62  } else {
63  LogDebug("StripDigiInfo")<<"APVs running in deconvolution mode (good time resolution)";
64  };
65 }
66 
68 }
69 
70 void
72  edm::ESHandle<SiStripBadStrip> deadChannelHandle;
73  iSetup.get<SiStripBadChannelRcd>().get(deadChannelHandle);
74 
75  unsigned int detId = det->geographicalId().rawId();
76  int numStrips = (det->specificTopology()).nstrips();
77 
78  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detId);
79  //storing the bad strip of the the module. the module is not removed but just signal put to 0
80  std::vector<bool>& badChannels = allBadChannels[detId];
81  badChannels.clear();
82  badChannels.insert(badChannels.begin(), numStrips, false);
83  for(SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
84  SiStripBadStrip::data fs = deadChannelHandle->decode(*it);
85  for(int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip) badChannels[strip] = true;
86  }
87  firstChannelsWithSignal[detId] = numStrips;
88  lastChannelsWithSignal[detId]= 0;
89 }
90 
91 void
93  theSiPileUpSignals->reset();
94  // This should be clear by after all calls to digitize(), but I might as well make sure
96 
97  //get gain noise pedestal lorentzAngle from ES handle
99  iSetup.getData(pdt);
100  setParticleDataTable(&*pdt);
102 }
103 
104 // Run the algorithm for a given module
105 // ------------------------------------
106 
107 void
108 SiStripDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
109  std::vector<PSimHit>::const_iterator inputEnd,
110  size_t inputBeginGlobalIndex,
111  const StripGeomDetUnit* det,
112  const GlobalVector& bfield,
113  const TrackerTopology *tTopo,
114  CLHEP::HepRandomEngine* engine) {
115  // produce SignalPoints for all SimHits in detector
116  unsigned int detID = det->geographicalId().rawId();
117  int numStrips = (det->specificTopology()).nstrips();
118 
119  std::vector<bool>& badChannels = allBadChannels[detID];
120  size_t thisFirstChannelWithSignal = numStrips;
121  size_t thisLastChannelWithSignal = numStrips;
122 
123  float langle = (lorentzAngleHandle.isValid()) ? lorentzAngleHandle->getLorentzAngle(detID) : 0.;
124 
125  std::vector<float> locAmpl(numStrips, 0.);
126 
127  // Loop over hits
128 
129  uint32_t detId = det->geographicalId().rawId();
130  // First: loop on the SimHits
131  if(CLHEP::RandFlat::shoot(engine) > inefficiency) {
132  AssociationInfoForChannel* pDetIDAssociationInfo; // I only need this if makeDigiSimLinks_ is true...
133  if( makeDigiSimLinks_ ) pDetIDAssociationInfo=&(associationInfoForDetId_[detId]); // ...so only search the map if that is the case
134  std::vector<float> previousLocalAmplitude; // Only used if makeDigiSimLinks_ is true. Needed to work out the change in amplitude.
135 
136  size_t simHitGlobalIndex=inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
137  for (std::vector<PSimHit>::const_iterator simHitIter = inputBegin; simHitIter != inputEnd; ++simHitIter, ++simHitGlobalIndex ) {
138  // skip hits not in this detector.
139  if((*simHitIter).detUnitId() != detId) {
140  continue;
141  }
142  // check TOF
143  if (std::fabs(simHitIter->tof() - cosmicShift - det->surface().toGlobal(simHitIter->localPosition()).mag()/30.) < tofCut && simHitIter->energyLoss()>0) {
144  if( makeDigiSimLinks_ ) previousLocalAmplitude=locAmpl; // Not needed except to make the sim link association.
145  size_t localFirstChannel = numStrips;
146  size_t localLastChannel = 0;
147  // process the hit
148  theSiHitDigitizer->processHit(&*simHitIter, *det, bfield, langle, locAmpl, localFirstChannel, localLastChannel, tTopo, engine);
149 
150  //APV Killer to simulate HIP effect
151  //------------------------------------------------------
152 
154  int pdg_id = simHitIter->particleType();
155  particle = pdt->particle(pdg_id);
156  if(particle != NULL){
157  float charge = particle->charge();
158  bool isHadron = particle->isHadron();
159  if(charge!=0 && isHadron){
160  if(CLHEP::RandFlat::shoot(engine) < APVSaturationProb){
161  int FirstAPV = localFirstChannel/128;
162  int LastAPV = localLastChannel/128;
163  std::cout << "-------------------HIP--------------" << std::endl;
164  std::cout << "Killing APVs " << FirstAPV << " - " <<LastAPV << " " << detID <<std::endl;
165  for(int strip = FirstAPV*128; strip < LastAPV*128 +128; ++strip) {
166  badChannels[strip] = true;
167  }
168  //doing like that I remove the signal information only after the
169  //stip that got the HIP but it remains the signal of the previous
170  //one. I'll make a further loop to remove all signal
171  }
172  }
173  }
174  }
175 
176 
177  if(thisFirstChannelWithSignal > localFirstChannel) thisFirstChannelWithSignal = localFirstChannel;
178  if(thisLastChannelWithSignal < localLastChannel) thisLastChannelWithSignal = localLastChannel;
179 
180  if( makeDigiSimLinks_ ) { // No need to do any of this if truth association was turned off in the configuration
181  for( size_t stripIndex=0; stripIndex<locAmpl.size(); ++stripIndex ) {
182  // Work out the amplitude from this SimHit from the difference of what it was before and what it is now
183  float signalFromThisSimHit=locAmpl[stripIndex]-previousLocalAmplitude[stripIndex];
184  if( signalFromThisSimHit!=0 ) { // If this SimHit had any contribution I need to record it.
185  auto& associationVector=(*pDetIDAssociationInfo)[stripIndex];
186  bool addNewEntry=true;
187  // Make sure the hit isn't in already. I've seen this a few times, it always seems to happen in pairs so I think
188  // it's something to do with the stereo strips.
189  for( auto& associationInfo : associationVector ) {
190  if( associationInfo.trackID==simHitIter->trackId() && associationInfo.eventID==simHitIter->eventId() ) {
191  // The hit is already in, so add this second contribution and move on
192  associationInfo.contributionToADC+=signalFromThisSimHit;
193  addNewEntry=false;
194  break;
195  }
196  } // end of loop over associationVector
197  // If the hit wasn't already in create a new association info structure.
198  if( addNewEntry ) associationVector.push_back( AssociationInfo{ simHitIter->trackId(), simHitIter->eventId(), signalFromThisSimHit, simHitGlobalIndex } );
199  } // end of "if( signalFromThisSimHit!=0 )"
200  } // end of loop over locAmpl strips
201  } // end of "if( makeDigiSimLinks_ )"
202  } // end of TOF check
203  } // end for
204  }
205  theSiPileUpSignals->add(detID, locAmpl, thisFirstChannelWithSignal, thisLastChannelWithSignal);
206 
207  if(firstChannelsWithSignal[detID] > thisFirstChannelWithSignal) firstChannelsWithSignal[detID] = thisFirstChannelWithSignal;
208  if(lastChannelsWithSignal[detID] < thisLastChannelWithSignal) lastChannelsWithSignal[detID] = thisLastChannelWithSignal;
209 }
210 
211 void
213  edm::DetSet<SiStripDigi>& outdigi,
214  edm::DetSet<SiStripRawDigi>& outrawdigi,
216  const StripGeomDetUnit *det,
217  edm::ESHandle<SiStripGain> & gainHandle,
218  edm::ESHandle<SiStripThreshold> & thresholdHandle,
219  edm::ESHandle<SiStripNoises> & noiseHandle,
220  edm::ESHandle<SiStripPedestals> & pedestalHandle,
221  CLHEP::HepRandomEngine* engine) {
222  unsigned int detID = det->geographicalId().rawId();
223  int numStrips = (det->specificTopology()).nstrips();
224 
225  const SiPileUpSignals::SignalMapType* theSignal(theSiPileUpSignals->getSignal(detID));
226 
227  std::vector<float> detAmpl(numStrips, 0.);
228  if(theSignal) {
229  for(const auto& amp : *theSignal) {
230  detAmpl[amp.first] = amp.second;
231  }
232  }
233 
234  //removing signal from the dead (and HIP effected) strips
235  std::vector<bool>& badChannels = allBadChannels[detID];
236  for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] = 0.;
237 
238  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
239  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
240  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
241 
242 // -----------------------------------------------------------
243 
244  auto& firstChannelWithSignal = firstChannelsWithSignal[detID];
245  auto& lastChannelWithSignal = lastChannelsWithSignal[detID];
246  auto iAssociationInfoByChannel=associationInfoForDetId_.find(detID); // Use an iterator so that I can easily remove it once finished
247 
248  if(zeroSuppression){
249  if(noise){
250  int RefStrip = int(numStrips/2.);
251  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
252  RefStrip++;
253  }
254  if(RefStrip<numStrips){
255  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange);
256  float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
257  theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,noiseRMS*theElectronPerADC/gainValue, engine);
258  }
259  }
260  DigitalVecType digis;
261  theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
262  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
263  // will be empty and the iterator will always equal associationInfoForDetId_.end().
264  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
265  for( const auto& iDigi : digis ) {
266  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
267  const std::vector<AssociationInfo>& associationInfo=associationInfoByChannel[iDigi.channel()];
268 
269  // Need to find the total from all sim hits, because this might not be the same as the total
270  // digitised due to noise or whatever.
271  float totalSimADC=0;
272  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
273  // Now I know that I can loop again and create the links
274  for( const auto& iAssociationInfo : associationInfo ) {
275  // Note simHitGlobalIndex has +1 because TrackerHitAssociator (the only place I can find this value being used)
276  // expects counting to start at 1, not 0.
277  outLink.push_back( StripDigiSimLink( iDigi.channel(), iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex+1, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
278  } // end of loop over associationInfo
279  } // end of loop over the digis
280  } // end of check that iAssociationInfoByChannel is a valid iterator
281  outdigi.data = digis;
282  }
283 
284 
285  if(!zeroSuppression){
286  //if(noise){
287  // the constant pedestal offset is needed because
288  // negative adc counts are not allowed in case
289  // Pedestal and CMN subtraction is performed.
290  // The pedestal value read from the conditions
291  // is pedValue and after the pedestal subtraction
292  // the baseline is zero. The Common Mode Noise
293  // is not subtracted from the negative adc counts
294  // channels. Adding pedOffset the baseline is set
295  // to pedOffset after pedestal subtraction and CMN
296  // is subtracted to all the channels since none of
297  // them has negative adc value. The pedOffset is
298  // treated as a constant component in the CMN
299  // estimation and subtracted as CMN.
300 
301 
302  //calculating the charge deposited on each APV and subtracting the shift
303  //------------------------------------------------------
304  if(BaselineShift){
305  theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
306  }
307 
308  //Adding the strip noise
309  //------------------------------------------------------
310  if(noise){
311  std::vector<float> noiseRMSv;
312  noiseRMSv.clear();
313  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
314 
315  if(SingleStripNoise){
316  for(int strip=0; strip< numStrips; ++strip){
317  if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
318  }
319 
320  } else {
321  int RefStrip = 0; //int(numStrips/2.);
322  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
323  RefStrip++;
324  }
325  if(RefStrip<numStrips){
326  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
327  for(int strip=0; strip< numStrips; ++strip){
328  if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
329  }
330  }
331  }
332 
333  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
334  }
335 
336  //adding the CMN
337  //------------------------------------------------------
338  if(CommonModeNoise){
339  float cmnRMS = 0.;
340  DetId detId(detID);
341  uint32_t SubDet = detId.subdetId();
342  if(SubDet==3){
343  cmnRMS = cmnRMStib;
344  }else if(SubDet==4){
345  cmnRMS = cmnRMStid;
346  }else if(SubDet==5){
347  cmnRMS = cmnRMStob;
348  }else if(SubDet==6){
349  cmnRMS = cmnRMStec;
350  }
351  cmnRMS *= theElectronPerADC;
352  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels, engine);
353  }
354 
355 
356  //Adding the pedestals
357  //------------------------------------------------------
358 
359  std::vector<float> vPeds;
360  vPeds.clear();
361  vPeds.insert(vPeds.begin(),numStrips,0.);
362 
363  if(RealPedestals){
364  for(int strip=0; strip< numStrips; ++strip){
365  if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
366  }
367  } else {
368  for(int strip=0; strip< numStrips; ++strip){
369  if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
370  }
371  }
372 
373  theSiNoiseAdder->addPedestals(detAmpl, vPeds);
374 
375 
376  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
377  // edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
378  //}else{
379 
380  DigitalRawVecType rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
381 
382  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
383  // will be empty and the iterator will always equal associationInfoForDetId_.end().
384  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
385  // N.B. For the raw digis the channel is inferred from the position in the vector.
386  // I'VE NOT TESTED THIS YET!!!!!
387  // ToDo Test this properly.
388  for( size_t channel=0; channel<rawdigis.size(); ++channel ) {
389  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
390  const auto iAssociationInfo=associationInfoByChannel.find(channel);
391  if( iAssociationInfo==associationInfoByChannel.end() ) continue; // Skip if there is no sim information for this channel (i.e. it's all noise)
392  const std::vector<AssociationInfo>& associationInfo=iAssociationInfo->second;
393 
394  // Need to find the total from all sim hits, because this might not be the same as the total
395  // digitised due to noise or whatever.
396  float totalSimADC=0;
397  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
398  // Now I know that I can loop again and create the links
399  for( const auto& iAssociationInfo : associationInfo ) {
400  // Note simHitGlobalIndex has +1 because TrackerHitAssociator (the only place I can find this value being used)
401  // expects counting to start at 1, not 0.
402  outLink.push_back( StripDigiSimLink( channel, iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex+1, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
403  } // end of loop over associationInfo
404  } // end of loop over the digis
405  } // end of check that iAssociationInfoByChannel is a valid iterator
406 
407  outrawdigi.data = rawdigis;
408 
409  //}
410  }
411 
412  // Now that I've finished with this entry in the map of associations, I can remove it.
413  // Note that there might not be an association if the ADC reading is from noise in which
414  // case associationIsValid will be false.
415  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) associationInfoForDetId_.erase(iAssociationInfoByChannel);
416 }
#define LogDebug(id)
unsigned short range
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
const std::unique_ptr< SiPileUpSignals > theSiPileUpSignals
edm::ESHandle< SiStripLorentzAngle > lorentzAngleHandle
void push_back(const T &t)
Definition: DetSet.h:68
AssociationInfoForDetId associationInfoForDetId_
Structure that holds the information on the SimTrack contributions. Only filled if makeDigiSimLinks_ ...
SiDigitalConverter::DigitalRawVecType DigitalRawVecType
SiDigitalConverter::DigitalVecType DigitalVecType
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
std::vector< unsigned int >::const_iterator ContainerIterator
void initializeEvent(const edm::EventSetup &iSetup)
#define NULL
Definition: scimark2.h:8
std::pair< ContainerIterator, ContainerIterator > Range
void initializeDetUnit(StripGeomDetUnit const *det, const edm::EventSetup &iSetup)
void digitize(edm::DetSet< SiStripDigi > &outDigis, edm::DetSet< SiStripRawDigi > &outRawDigis, edm::DetSet< StripDigiSimLink > &outLink, const StripGeomDetUnit *stripdet, edm::ESHandle< SiStripGain > &, edm::ESHandle< SiStripThreshold > &, edm::ESHandle< SiStripNoises > &, edm::ESHandle< SiStripPedestals > &, CLHEP::HepRandomEngine *)
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:35
double charge(const std::vector< uint8_t > &Ampls)
const std::unique_ptr< const SiGaussianTailNoiseAdder > theSiNoiseAdder
const std::unique_ptr< SiHitDigitizer > theSiHitDigitizer
void getData(T &iHolder) const
Definition: EventSetup.h:67
std::map< unsigned int, size_t > firstChannelsWithSignal
std::map< int, Amplitude > SignalMapType
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, size_t inputBeginGlobalIndex, const StripGeomDetUnit *stripdet, const GlobalVector &bfield, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
const std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
const std::unique_ptr< SiStripFedZeroSuppression > theSiZeroSuppress
std::pair< ContainerIterator, ContainerIterator > Range
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:72
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::map< unsigned int, std::vector< bool > > allBadChannels
tuple conf
Definition: dbtoconf.py:185
void setParticleDataTable(const ParticleDataTable *pardt)
Definition: DetId.h:18
unsigned short firstStrip
std::map< int, std::vector< AssociationInfo > > AssociationInfoForChannel
const T & get() const
Definition: EventSetup.h:55
SiStripDigitizerAlgorithm(const edm::ParameterSet &conf)
collection_type data
Definition: DetSet.h:78
std::pair< ContainerIterator, ContainerIterator > Range
std::map< unsigned int, size_t > lastChannelsWithSignal
tuple cout
Definition: gather_cfg.py:121
volatile std::atomic< bool > shutdown_flag false
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:48
bool isValid() const
Definition: ESHandle.h:37
const ParticleDataTable * pdt