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, eng)),
55  theSiPileUpSignals(new SiPileUpSignals()),
56  theSiNoiseAdder(new SiGaussianTailNoiseAdder(theThreshold, eng)),
57  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC)),
58  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo)),
59  theFlatDistribution(new CLHEP::RandFlat(eng, 0., 1.)) {
60 
61  if (peakMode) {
62  LogDebug("StripDigiInfo")<<"APVs running in peak mode (poor time resolution)";
63  } else {
64  LogDebug("StripDigiInfo")<<"APVs running in deconvolution mode (good time resolution)";
65  };
66 }
67 
69 }
70 
71 void
73  edm::ESHandle<SiStripBadStrip> deadChannelHandle;
74  iSetup.get<SiStripBadChannelRcd>().get(deadChannelHandle);
75 
76  unsigned int detId = det->geographicalId().rawId();
77  int numStrips = (det->specificTopology()).nstrips();
78 
79  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detId);
80  //storing the bad strip of the the module. the module is not removed but just signal put to 0
81  std::vector<bool>& badChannels = allBadChannels[detId];
82  badChannels.clear();
83  badChannels.insert(badChannels.begin(), numStrips, false);
84  for(SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
85  SiStripBadStrip::data fs = deadChannelHandle->decode(*it);
86  for(int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip) badChannels[strip] = true;
87  }
88  firstChannelsWithSignal[detId] = numStrips;
89  lastChannelsWithSignal[detId]= 0;
90 }
91 
92 void
94  theSiPileUpSignals->reset();
95  // This should be clear by after all calls to digitize(), but I might as well make sure
97 
98  //get gain noise pedestal lorentzAngle from ES handle
100  iSetup.getData(pdt);
101  setParticleDataTable(&*pdt);
103 }
104 
105 // Run the algorithm for a given module
106 // ------------------------------------
107 
108 void
109 SiStripDigitizerAlgorithm::accumulateSimHits(std::vector<PSimHit>::const_iterator inputBegin,
110  std::vector<PSimHit>::const_iterator inputEnd,
111  size_t inputBeginGlobalIndex,
112  const StripGeomDetUnit* det,
113  const GlobalVector& bfield,
114  const TrackerTopology *tTopo) {
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(theFlatDistribution->fire()>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);
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){
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  unsigned int detID = det->geographicalId().rawId();
222  int numStrips = (det->specificTopology()).nstrips();
223 
224  const SiPileUpSignals::SignalMapType* theSignal(theSiPileUpSignals->getSignal(detID));
225 
226  std::vector<float> detAmpl(numStrips, 0.);
227  if(theSignal) {
228  for(const auto& amp : *theSignal) {
229  detAmpl[amp.first] = amp.second;
230  }
231  }
232 
233  //removing signal from the dead (and HIP effected) strips
234  std::vector<bool>& badChannels = allBadChannels[detID];
235  for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] = 0.;
236 
237  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
238  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
239  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
240 
241 // -----------------------------------------------------------
242 
243  auto& firstChannelWithSignal = firstChannelsWithSignal[detID];
244  auto& lastChannelWithSignal = lastChannelsWithSignal[detID];
245  auto iAssociationInfoByChannel=associationInfoForDetId_.find(detID); // Use an iterator so that I can easily remove it once finished
246 
247  if(zeroSuppression){
248  if(noise){
249  int RefStrip = int(numStrips/2.);
250  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
251  RefStrip++;
252  }
253  if(RefStrip<numStrips){
254  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange);
255  float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
256  theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,noiseRMS*theElectronPerADC/gainValue);
257  }
258  }
259  DigitalVecType digis;
260  theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
261  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
262  // will be empty and the iterator will always equal associationInfoForDetId_.end().
263  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
264  for( const auto& iDigi : digis ) {
265  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
266  const std::vector<AssociationInfo>& associationInfo=associationInfoByChannel[iDigi.channel()];
267 
268  // Need to find the total from all sim hits, because this might not be the same as the total
269  // digitised due to noise or whatever.
270  float totalSimADC=0;
271  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
272  // Now I know that I can loop again and create the links
273  for( const auto& iAssociationInfo : associationInfo ) {
274  // Note simHitGlobalIndex has +1 because TrackerHitAssociator (the only place I can find this value being used)
275  // expects counting to start at 1, not 0.
276  outLink.push_back( StripDigiSimLink( iDigi.channel(), iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex+1, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
277  } // end of loop over associationInfo
278  } // end of loop over the digis
279  } // end of check that iAssociationInfoByChannel is a valid iterator
280  outdigi.data = digis;
281  }
282 
283 
284  if(!zeroSuppression){
285  //if(noise){
286  // the constant pedestal offset is needed because
287  // negative adc counts are not allowed in case
288  // Pedestal and CMN subtraction is performed.
289  // The pedestal value read from the conditions
290  // is pedValue and after the pedestal subtraction
291  // the baseline is zero. The Common Mode Noise
292  // is not subtracted from the negative adc counts
293  // channels. Adding pedOffset the baseline is set
294  // to pedOffset after pedestal subtraction and CMN
295  // is subtracted to all the channels since none of
296  // them has negative adc value. The pedOffset is
297  // treated as a constant component in the CMN
298  // estimation and subtracted as CMN.
299 
300 
301  //calculating the charge deposited on each APV and subtracting the shift
302  //------------------------------------------------------
303  if(BaselineShift){
304  theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
305  }
306 
307  //Adding the strip noise
308  //------------------------------------------------------
309  if(noise){
310  std::vector<float> noiseRMSv;
311  noiseRMSv.clear();
312  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
313 
314  if(SingleStripNoise){
315  for(int strip=0; strip< numStrips; ++strip){
316  if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
317  }
318 
319  } else {
320  int RefStrip = 0; //int(numStrips/2.);
321  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
322  RefStrip++;
323  }
324  if(RefStrip<numStrips){
325  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
326  for(int strip=0; strip< numStrips; ++strip){
327  if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
328  }
329  }
330  }
331 
332  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv);
333  }
334 
335  //adding the CMN
336  //------------------------------------------------------
337  if(CommonModeNoise){
338  float cmnRMS = 0.;
339  DetId detId(detID);
340  uint32_t SubDet = detId.subdetId();
341  if(SubDet==3){
342  cmnRMS = cmnRMStib;
343  }else if(SubDet==4){
344  cmnRMS = cmnRMStid;
345  }else if(SubDet==5){
346  cmnRMS = cmnRMStob;
347  }else if(SubDet==6){
348  cmnRMS = cmnRMStec;
349  }
350  cmnRMS *= theElectronPerADC;
351  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels);
352  }
353 
354 
355  //Adding the pedestals
356  //------------------------------------------------------
357 
358  std::vector<float> vPeds;
359  vPeds.clear();
360  vPeds.insert(vPeds.begin(),numStrips,0.);
361 
362  if(RealPedestals){
363  for(int strip=0; strip< numStrips; ++strip){
364  if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
365  }
366  } else {
367  for(int strip=0; strip< numStrips; ++strip){
368  if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
369  }
370  }
371 
372  theSiNoiseAdder->addPedestals(detAmpl, vPeds);
373 
374 
375  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
376  // edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
377  //}else{
378 
379  DigitalRawVecType rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
380 
381  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
382  // will be empty and the iterator will always equal associationInfoForDetId_.end().
383  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
384  // N.B. For the raw digis the channel is inferred from the position in the vector.
385  // I'VE NOT TESTED THIS YET!!!!!
386  // ToDo Test this properly.
387  for( size_t channel=0; channel<rawdigis.size(); ++channel ) {
388  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
389  const auto iAssociationInfo=associationInfoByChannel.find(channel);
390  if( iAssociationInfo==associationInfoByChannel.end() ) continue; // Skip if there is no sim information for this channel (i.e. it's all noise)
391  const std::vector<AssociationInfo>& associationInfo=iAssociationInfo->second;
392 
393  // Need to find the total from all sim hits, because this might not be the same as the total
394  // digitised due to noise or whatever.
395  float totalSimADC=0;
396  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
397  // Now I know that I can loop again and create the links
398  for( const auto& iAssociationInfo : associationInfo ) {
399  // Note simHitGlobalIndex has +1 because TrackerHitAssociator (the only place I can find this value being used)
400  // expects counting to start at 1, not 0.
401  outLink.push_back( StripDigiSimLink( channel, iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex+1, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
402  } // end of loop over associationInfo
403  } // end of loop over the digis
404  } // end of check that iAssociationInfoByChannel is a valid iterator
405 
406  outrawdigi.data = rawdigis;
407 
408  //}
409  }
410 
411  // Now that I've finished with this entry in the map of associations, I can remove it.
412  // Note that there might not be an association if the ADC reading is from noise in which
413  // case associationIsValid will be false.
414  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) associationInfoForDetId_.erase(iAssociationInfoByChannel);
415 }
#define LogDebug(id)
unsigned short range
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
const std::unique_ptr< CLHEP::RandFlat > theFlatDistribution
const std::unique_ptr< SiPileUpSignals > theSiPileUpSignals
edm::ESHandle< SiStripLorentzAngle > lorentzAngleHandle
void push_back(const T &t)
Definition: DetSet.h:68
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 > &)
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
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)
const std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
SiStripDigitizerAlgorithm(const edm::ParameterSet &conf, CLHEP::HepRandomEngine &)
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
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:44
bool isValid() const
Definition: ESHandle.h:37
void initializeDetUnit(StripGeomDetUnit *det, const edm::EventSetup &iSetup)
const ParticleDataTable * pdt