CMS 3D CMS Logo

DataMixingSiStripMCDigiWorker.cc
Go to the documentation of this file.
1 // File: DataMixingSiStripMCDigiWorker.cc
2 // Description: see DataMixingSiStripMCDigiWorker.h
3 // Author: Mike Hildreth, University of Notre Dame
4 //
5 //--------------------------------------------
6 
7 #include <map>
8 #include <memory>
19 //
20 #include "CLHEP/Random/RandFlat.h"
21 //
23 
24 using namespace std;
25 
26 namespace edm
27 {
28 
29  // Virtual constructor
30 
31  DataMixingSiStripMCDigiWorker::DataMixingSiStripMCDigiWorker() { }
32 
33  // Constructor
34  DataMixingSiStripMCDigiWorker::DataMixingSiStripMCDigiWorker(const edm::ParameterSet& ps,
35  edm::ConsumesCollector && iC) :
36  label_(ps.getParameter<std::string>("Label")),
37  gainLabel(ps.getParameter<std::string>("Gain")),
38  SingleStripNoise(ps.getParameter<bool>("SingleStripNoise")),
39  peakMode(ps.getParameter<bool>("APVpeakmode")),
40  theThreshold(ps.getParameter<double>("NoiseSigmaThreshold")),
41  theElectronPerADC(ps.getParameter<double>( peakMode ? "electronPerAdcPeak" : "electronPerAdcDec" )),
42  APVSaturationFromHIP_(ps.getParameter<bool>("APVSaturationFromHIP")),
43  theFedAlgo(ps.getParameter<int>("FedAlgorithm_PM")),
44  geometryType(ps.getParameter<std::string>("GeometryType")),
45  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo)),
46  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC, false)) // no premixing
47 
48  {
49 
50  // get the subdetector names
51  // this->getSubdetectorNames(); //something like this may be useful to check what we are supposed to do...
52 
53  // declare the products to produce
54 
55  SistripLabelSig_ = ps.getParameter<edm::InputTag>("SistripLabelSig");
56  SiStripPileInputTag_ = ps.getParameter<edm::InputTag>("SiStripPileInputTag");
57 
58  SiStripDigiCollectionDM_ = ps.getParameter<std::string>("SiStripDigiCollectionDM");
59  SistripAPVListDM_= ps.getParameter<std::string>("SiStripAPVListDM");
60 
61 
63  SistripAPVLabelSig_ = ps.getParameter<edm::InputTag>("SistripAPVLabelSig");
64  SiStripAPVPileInputTag_ = ps.getParameter<edm::InputTag>("SistripAPVPileInputTag");
65  iC.consumes< std::vector<std::pair<int,std::bitset<6>> > >(SistripAPVLabelSig_);
66  }
68  // clear local storage for this event
69  SiHitStorage_.clear();
70 
72  if ( ! rng.isAvailable()) {
73  throw cms::Exception("Psiguration")
74  << "SiStripDigitizer requires the RandomNumberGeneratorService\n"
75  "which is not present in the psiguration file. You must add the service\n"
76  "in the configuration file or remove the modules that require it.";
77  }
78 
80 
81  // theSiZeroSuppress = new SiStripFedZeroSuppression(theFedAlgo);
82  //theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC));
83 
84  }
85 
86 
87  // Virtual destructor needed.
89  }
90 
92  // initialize individual detectors so we can copy real digitization code:
93 
95 
96  for(auto iu = pDD->detUnits().begin(); iu != pDD->detUnits().end(); ++iu) {
97  unsigned int detId = (*iu)->geographicalId().rawId();
98  DetId idet=DetId(detId);
99  unsigned int isub=idet.subdetId();
100  if((isub == StripSubdetector::TIB) ||
101  (isub == StripSubdetector::TID) ||
102  (isub == StripSubdetector::TOB) ||
103  (isub == StripSubdetector::TEC)) {
104 
105  auto stripdet = dynamic_cast<StripGeomDetUnit const*>((*iu));
106  assert(stripdet != 0);
107  DMinitializeDetUnit(stripdet, iSetup);
108  }
109  }
110  }
111 
112 
114 
115  edm::ESHandle<SiStripBadStrip> deadChannelHandle;
116  iSetup.get<SiStripBadChannelRcd>().get(deadChannelHandle);
117 
118  unsigned int detId = det->geographicalId().rawId();
119  int numStrips = (det->specificTopology()).nstrips();
120 
121  SiStripBadStrip::Range detBadStripRange = deadChannelHandle->getRange(detId);
122  //storing the bad strip of the the module. the module is not removed but just signal put to 0
123  std::vector<bool>& badChannels = allBadChannels[detId];
124  std::vector<bool>& hipChannels = allHIPChannels[detId];
125  badChannels.clear();
126  badChannels.insert(badChannels.begin(), numStrips, false);
127  hipChannels.clear();
128  hipChannels.insert(hipChannels.begin(), numStrips, false);
129 
130  for(SiStripBadStrip::ContainerIterator it = detBadStripRange.first; it != detBadStripRange.second; ++it) {
131  SiStripBadStrip::data fs = deadChannelHandle->decode(*it);
132  for(int strip = fs.firstStrip; strip < fs.firstStrip + fs.range; ++strip) badChannels[strip] = true;
133  }
134  firstChannelsWithSignal[detId] = numStrips;
135  lastChannelsWithSignal[detId]= 0;
136  }
137 
138 
140  // fill in maps of hits
141 
143 
144  if( e.getByLabel(SistripLabelSig_,input) ) {
145  OneDetectorMap LocalMap;
146 
147  //loop on all detsets (detectorIDs) inside the input collection
148  edm::DetSetVector<SiStripDigi>::const_iterator DSViter=input->begin();
149  for (; DSViter!=input->end();DSViter++){
150 
151 #ifdef DEBUG
152  LogDebug("DataMixingSiStripMCDigiWorker") << "Processing DetID " << DSViter->id;
153 #endif
154 
155  LocalMap.clear();
156  LocalMap.reserve((DSViter->data).size());
157  LocalMap.insert(LocalMap.end(),(DSViter->data).begin(),(DSViter->data).end());
158 
159  SiHitStorage_.insert( SiGlobalIndex::value_type( DSViter->id, LocalMap ) );
160  }
161 
162  }
163 
164  // keep here for future reference. In current implementation, HIP killing is done once in PU file
165  /* if(APVSaturationFromHIP_) {
166  Handle<std::vector<std::pair<int,std::bitset<6>> > > APVinput;
167 
168  if( e.getByLabel(SistripAPVLabelSig_,APVinput) ) {
169 
170  std::vector<std::pair<int,std::bitset<6>> >::const_iterator entry = APVinput->begin();
171  for( ; entry != APVinput->end(); entry++) {
172  theAffectedAPVmap_.insert(APVMap::value_type(entry->first, entry->second));
173  }
174  }
175  } */
176 
177  } // end of addSiStripSignals
178 
179 
180 
181  void DataMixingSiStripMCDigiWorker::addSiStripPileups(const int bcr, const EventPrincipal *ep, unsigned int eventNr,
182  ModuleCallingContext const* mcc) {
183  LogDebug("DataMixingSiStripMCDigiWorker") <<"\n===============> adding pileups from event "<<ep->id()<<" for bunchcrossing "<<bcr;
184 
185  // fill in maps of hits; same code as addSignals, except now applied to the pileup events
186 
187  std::shared_ptr<Wrapper<edm::DetSetVector<SiStripDigi> > const> inputPTR =
188  getProductByTag<edm::DetSetVector<SiStripDigi> >(*ep, SiStripPileInputTag_, mcc);
189 
190  if(inputPTR ) {
191 
192  const edm::DetSetVector<SiStripDigi> *input = const_cast< edm::DetSetVector<SiStripDigi> * >(inputPTR->product());
193 
194  // Handle< edm::DetSetVector<SiStripDigi> > input;
195 
196  // if( e->getByLabel(Sistripdigi_collectionPile_.label(),SistripLabelPile_.label(),input) ) {
197 
198  OneDetectorMap LocalMap;
199 
200  //loop on all detsets (detectorIDs) inside the input collection
202  for (; DSViter!=input->end();DSViter++){
203 
204 #ifdef DEBUG
205  LogDebug("DataMixingSiStripMCDigiWorker") << "Pileups: Processing DetID " << DSViter->id;
206 #endif
207 
208  // find correct local map (or new one) for this detector ID
209 
210  SiGlobalIndex::const_iterator itest;
211 
212  itest = SiHitStorage_.find(DSViter->id);
213 
214  if(itest!=SiHitStorage_.end()) { // this detID already has hits, add to existing map
215 
216  LocalMap = itest->second;
217 
218  // fill in local map with extra channels
219  LocalMap.insert(LocalMap.end(),(DSViter->data).begin(),(DSViter->data).end());
220  std::stable_sort(LocalMap.begin(),LocalMap.end(),DataMixingSiStripMCDigiWorker::StrictWeakOrdering());
221  SiHitStorage_[DSViter->id]=LocalMap;
222 
223  }
224  else{ // fill local storage with this information, put in global collection
225 
226  LocalMap.clear();
227  LocalMap.reserve((DSViter->data).size());
228  LocalMap.insert(LocalMap.end(),(DSViter->data).begin(),(DSViter->data).end());
229 
230  SiHitStorage_.insert( SiGlobalIndex::value_type( DSViter->id, LocalMap ) );
231  }
232  }
233 
235  std::shared_ptr<Wrapper<std::vector<std::pair<int,std::bitset<6>> > > const> inputAPVPTR =
236  getProductByTag< std::vector<std::pair<int,std::bitset<6>> > >(*ep, SiStripAPVPileInputTag_, mcc);
237 
238  if(inputAPVPTR) {
239 
240  const std::vector<std::pair<int,std::bitset<6>> > *APVinput = const_cast< std::vector<std::pair<int,std::bitset<6>> > * >(inputAPVPTR->product());
241 
242  std::vector<std::pair<int,std::bitset<6>> >::const_iterator entry = APVinput->begin();
243  for( ; entry != APVinput->end(); entry++) {
244  theAffectedAPVmap_.insert( APVMap::value_type(entry->first, entry->second ));
245  }
246  }
247  }
248  }
249  }
250 
251 
252 
254 
255  // set up machinery to do proper noise adding:
256  edm::ESHandle<SiStripGain> gainHandle;
257  edm::ESHandle<SiStripNoises> noiseHandle;
258  edm::ESHandle<SiStripThreshold> thresholdHandle;
259  edm::ESHandle<SiStripPedestals> pedestalHandle;
260  edm::ESHandle<SiStripBadStrip> deadChannelHandle;
261  iSetup.get<SiStripGainSimRcd>().get(gainLabel,gainHandle);
262  iSetup.get<SiStripNoisesRcd>().get(noiseHandle);
263  iSetup.get<SiStripThresholdRcd>().get(thresholdHandle);
264  iSetup.get<SiStripPedestalsRcd>().get(pedestalHandle);
265 
267  CLHEP::HepRandomEngine* engine = &rng->getEngine(e.streamID());
268 
269  std::map< int,std::bitset<6>> DeadAPVList;
270  DeadAPVList.clear();
271 
272 
273  // First, have to convert all ADC counts to raw pulse heights so that values can be added properly
274  // In PreMixing, pulse heights are saved with ADC = sqrt(9.0*PulseHeight) - have to undo.
275 
276  // This is done here because it's the only place we have access to EventSetup
277  // Simultaneously, merge lists of hit channels in each DetId.
278  // Signal Digis are in the list first, have to merge lists of hit strips on the fly,
279  // add signals on duplicates later
280 
281  OneDetectorRawMap LocalRawMap;
282 
283  // Now, loop over hits and add them to the map in the proper sorted order
284  // Note: We are assuming that the hits from the Signal events have been created in
285  // "PreMix" mode, rather than in the standard ADC conversion routines. If not, this
286  // doesn't work at all.
287 
288  // At the moment, both Signal and Reconstituted PU hits have the same compression algorithm.
289  // If this were different, and one needed gains, the conversion back to pulse height can only
290  // be done in this routine. So, yes, there is an extra loop over hits here in the current code,
291  // because, in principle, one could convert to pulse height during the read/store phase.
292 
293  for(SiGlobalIndex::const_iterator IDet = SiHitStorage_.begin();
294  IDet != SiHitStorage_.end(); IDet++) {
295 
296  uint32_t detID = IDet->first;
297 
298  OneDetectorMap LocalMap = IDet->second;
299 
300  //loop over hit strips for this DetId, do conversion to pulse height, store.
301 
302  LocalRawMap.clear();
303 
304  OneDetectorMap::const_iterator iLocal = LocalMap.begin();
305  for(;iLocal != LocalMap.end(); ++iLocal) {
306 
307  uint16_t currentStrip = iLocal->strip();
308  float signal = float(iLocal->adc());
309  if(iLocal->adc() == 1022) signal = 1500.; // average values for overflows
310  if(iLocal->adc() == 1023) signal = 3000.;
311 
312  //convert signals back to raw counts
313 
314  float ReSignal = signal*signal/9.0; // The PreMixing conversion is adc = sqrt(9.0*pulseHeight)
315 
316  RawDigi NewRawDigi = std::make_pair(currentStrip,ReSignal);
317 
318  LocalRawMap.push_back(NewRawDigi);
319 
320  }
321 
322  // save information for this detiD into global map
323  SiRawDigis_.insert( SiGlobalRawIndex::value_type( detID, LocalRawMap ) );
324  }
325 
326  // If we are killing APVs, merge list of dead ones before we digitize
327 
328  int NumberOfBxBetweenHIPandEvent=1e3;
329 
330  if(APVSaturationFromHIP_) {
331 
332  // calculate affected BX parameter
333 
334  bool HasAtleastOneAffectedAPV=false;
335  while(!HasAtleastOneAffectedAPV){
336  for(int bx=floor(300.0/25.0);bx>0;bx--){ //Reminder: make these numbers not hard coded!!
337  float temp=CLHEP::RandFlat::shoot(engine)<0.5?1:0;
338  if(temp==1 && bx<NumberOfBxBetweenHIPandEvent){
339  NumberOfBxBetweenHIPandEvent=bx;
340  HasAtleastOneAffectedAPV=true;
341  }
342  }
343  }
344 
345  APVMap::const_iterator iAPVchk;
346  uint32_t formerID = 0;
347  uint32_t currentID;
348  std::bitset<6> NewAPVBits;
349 
350  for(APVMap::const_iterator iAPV = theAffectedAPVmap_.begin();
351  iAPV != theAffectedAPVmap_.end(); ++iAPV) {
352 
353  currentID = iAPV->first;
354 
355  if (currentID == formerID) { // we have to OR these
356  for( int ibit=0; ibit<6; ++ibit){
357  NewAPVBits[ibit] = NewAPVBits[ibit]||(iAPV->second)[ibit];
358  }
359  }
360  else {
361  DeadAPVList[currentID]=NewAPVBits;
362  //save pointers for next iteration
363  formerID = currentID;
364  NewAPVBits = iAPV->second;
365  }
366 
367  iAPVchk = iAPV;
368  if((++iAPVchk) == theAffectedAPVmap_.end()) { //make sure not to lose the last one
369  DeadAPVList[currentID]=NewAPVBits;
370  }
371  }
372 
373  }
374  //
375 
376  // Ok, done with merging raw signals and APVs - now add signals on duplicate strips
377 
378  // collection of Digis to put in the event
379  std::vector< edm::DetSet<SiStripDigi> > vSiStripDigi;
380 
381  // loop through our collection of detectors, merging hits and making a new list of "signal" digis
382 
383  // clear some temporary storage for later digitization:
384 
385  signals_.clear();
386 
387  // big loop over Detector IDs:
388  for(SiGlobalRawIndex::const_iterator IDet = SiRawDigis_.begin();
389  IDet != SiRawDigis_.end(); IDet++) {
390 
391  uint32_t detID = IDet->first;
392 
393  SignalMapType Signals;
394  Signals.clear();
395 
396  OneDetectorRawMap LocalMap = IDet->second;
397 
398  //counter variables
399  int formerStrip = -1;
400  int currentStrip;
401  float ADCSum = 0;
402 
403  //loop over hit strips for this DetId, add duplicates
404 
405  OneDetectorRawMap::const_iterator iLocalchk;
406  OneDetectorRawMap::const_iterator iLocal = LocalMap.begin();
407  for(;iLocal != LocalMap.end(); ++iLocal) {
408 
409  currentStrip = iLocal->first; // strip is first element
410 
411  if (currentStrip == formerStrip) { // we have to add these digis together
412 
413  ADCSum+=iLocal->second ; // raw pulse height is second element.
414  }
415  else{
416  if(formerStrip!=-1){
417  Signals.insert( std::make_pair(formerStrip, ADCSum));
418  }
419  // save pointers for next iteration
420  formerStrip = currentStrip;
421  ADCSum = iLocal->second; // lone ADC
422  }
423 
424  iLocalchk = iLocal;
425  if((++iLocalchk) == LocalMap.end()) { //make sure not to lose the last one
426  Signals.insert( std::make_pair(formerStrip, ADCSum));
427  }
428  }
429  // save merged map:
430  signals_.insert( std::make_pair( detID, Signals));
431  }
432 
433  //Now, do noise, zero suppression, take into account bad channels, etc.
434  // This section stolen from SiStripDigitizerAlgorithm
435  // must loop over all detIds in the tracker to get all of the noise added properly.
436  for(TrackingGeometry::DetUnitContainer::const_iterator iu = pDD->detUnits().begin(); iu != pDD->detUnits().end(); iu ++){
437 
438  const StripGeomDetUnit* sgd = dynamic_cast<const StripGeomDetUnit*>((*iu));
439  if (sgd != 0){
440 
441  uint32_t detID = sgd->geographicalId().rawId();
442 
443  edm::DetSet<SiStripDigi> SSD(detID); // Make empty collection with this detector ID
444 
445  int numStrips = (sgd->specificTopology()).nstrips();
446 
447  // see if there is some signal on this detector
448 
449  const SignalMapType* theSignal(getSignal(detID));
450 
451  std::vector<float> detAmpl(numStrips, 0.);
452  if(theSignal) {
453  for(const auto& amp : *theSignal) {
454  detAmpl[amp.first] = amp.second;
455  }
456  }
457 
458  //removing signal from the dead (and HIP effected) strips
459  std::vector<bool>& badChannels = allBadChannels[detID];
460 
461  for(int strip =0; strip < numStrips; ++strip) {
462  if(badChannels[strip]) detAmpl[strip] = 0.;
463  }
464 
466  std::bitset<6> & bs=DeadAPVList[detID];
467 
468  if(bs.any()){
469  // Here below is the scaling function which describes the evolution of the baseline (i.e. how the charge is suppressed).
470  // This must be replaced as soon as we have a proper modeling of the baseline evolution from VR runs
471  float Shift=1-NumberOfBxBetweenHIPandEvent/floor(300.0/25.0); //Reminder: make these numbers not hardcoded!!
472  float randomX=CLHEP::RandFlat::shoot(engine);
473  float scalingValue=(randomX-Shift)*10.0/7.0-3.0/7.0;
474 
475  for(int strip =0; strip < numStrips; ++strip) {
476  if(!badChannels[strip] && bs[strip/128]==1){
477  detAmpl[strip] *=scalingValue>0?scalingValue:0.0;
478  }
479  }
480  }
481  }
482 
483  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
484  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
485 
486  // Gain conversion is already done during signal adding
487  //convert our signals back to raw counts so that we can add noise properly:
488 
489  /*
490  if(theSignal) {
491  for(unsigned int iv = 0; iv!=detAmpl.size(); iv++) {
492  float signal = detAmpl[iv];
493  if(signal > 0) {
494  float gainValue = gainHandle->getStripGain(iv, detGainRange);
495  signal *= theElectronPerADC/gainValue;
496  detAmpl[iv] = signal;
497  }
498  }
499  }
500  */
501 
502  //SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
503 
504  // -----------------------------------------------------------
505 
506  size_t firstChannelWithSignal = 0;
507  size_t lastChannelWithSignal = numStrips;
508 
509  if(SingleStripNoise){
510  // std::cout<<"In SSN, detId="<<detID<<std::endl;
511  std::vector<float> noiseRMSv;
512  noiseRMSv.clear();
513  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
514  for(int strip=0; strip< numStrips; ++strip){
515  if(!badChannels[strip]){
516  float gainValue = gainHandle->getStripGain(strip, detGainRange);
517  noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC/gainValue;
518  //std::cout<<"<SiStripDigitizerAlgorithm::digitize>: gainValue: "<<gainValue<<"\tnoiseRMSv["<<strip<<"]: "<<noiseRMSv[strip]<<std::endl;
519  }
520  }
521  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
522  } else {
523  int RefStrip = int(numStrips/2.);
524  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
525  RefStrip++;
526  }
527  if(RefStrip<numStrips){
528  float RefgainValue = gainHandle->getStripGain(RefStrip, detGainRange);
529  float RefnoiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC/RefgainValue;
530 
531  theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,RefnoiseRMS, engine);
532  //std::cout<<"<SiStripDigitizerAlgorithm::digitize>: RefgainValue: "<<RefgainValue<<"\tRefnoiseRMS: "<<RefnoiseRMS<<std::endl;
533  }
534  }
535 
536  DigitalVecType digis;
537  theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
538 
539 
540  SSD.data = digis;
541  // if(digis.size() > 0) {
542  // std::cout << " Real SiS Mixed Digi: " << detID << " ADC values ";
543  // for(const auto& iDigi : digis) { std::cout << iDigi.adc() << " " ;}
544  // std::cout << std::endl;
545  //}
546 
547  // stick this into the global vector of detector info
548  vSiStripDigi.push_back(SSD);
549 
550  } // end of loop over one detector
551 
552  } // end of big loop over all detector IDs
553 
554  // put the collection of digis in the event
555  LogInfo("DataMixingSiStripMCDigiWorker") << "total # Merged strips: " << vSiStripDigi.size() ;
556 
557  // make new digi collection
558 
559  std::unique_ptr< edm::DetSetVector<SiStripDigi> > MySiStripDigis(new edm::DetSetVector<SiStripDigi>(vSiStripDigi) );
560 
561  // put collection
562 
563  e.put(std::move(MySiStripDigis), SiStripDigiCollectionDM_ );
564 
565  // clear local storage for this event
566  SiHitStorage_.clear();
567  SiRawDigis_.clear();
568  signals_.clear();
569  }
570 
571 } //edm
#define LogDebug(id)
unsigned short range
T getParameter(std::string const &) const
std::map< unsigned int, size_t > firstChannelsWithSignal
std::map< unsigned int, std::vector< bool > > allHIPChannels
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
void reserve(size_t s)
Definition: DetSetVector.h:150
std::vector< unsigned int >::const_iterator ContainerIterator
std::pair< uint16_t, Amplitude > RawDigi
EventID const & id() const
virtual void initializeEvent(const edm::Event &e, edm::EventSetup const &iSetup)
virtual const StripTopology & specificTopology() const
Returns a reference to the strip proxy topology.
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
static std::string const input
Definition: EdmProvDump.cc:44
SiDigitalConverter::DigitalVecType DigitalVecType
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
static float getNoise(uint16_t strip, const Range &range)
Definition: SiStripNoises.h:72
static float getStripGain(const uint16_t &strip, const SiStripApvGain::Range &range)
Definition: SiStripGain.h:68
void addSiStripPileups(const int bcr, const edm::EventPrincipal *, unsigned int EventId, ModuleCallingContext const *)
bool isAvailable() const
Definition: Service.h:46
std::pair< ContainerIterator, ContainerIterator > Range
void DMinitializeDetUnit(StripGeomDetUnit const *det, const edm::EventSetup &iSetup)
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
#define end
Definition: vmac.h:37
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
std::unique_ptr< SiGaussianTailNoiseAdder > theSiNoiseAdder
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:413
iterator end()
Return the off-the-end iterator.
Definition: DetSetVector.h:361
Definition: DetId.h:18
unsigned short firstStrip
const T & get() const
Definition: EventSetup.h:56
const SignalMapType * getSignal(uint32_t detID) const
void insert(detset const &s)
Insert the given DetSet.
Definition: DetSetVector.h:235
collection_type data
Definition: DetSet.h:78
std::map< unsigned int, std::vector< bool > > allBadChannels
#define begin
Definition: vmac.h:30
const Range getRange(const uint32_t detID) const
HLT enums.
const Range getRange(const uint32_t detID) const
std::pair< ContainerIterator, ContainerIterator > Range
std::unique_ptr< SiTrivialDigitalConverter > theSiDigitalConverter
StreamID streamID() const
Definition: Event.h:81
std::unique_ptr< SiStripFedZeroSuppression > theSiZeroSuppress
iterator begin()
Return an iterator to the first DetSet.
Definition: DetSetVector.h:346
std::map< unsigned int, size_t > lastChannelsWithSignal
std::pair< ContainerIterator, ContainerIterator > Range
Definition: SiStripNoises.h:48
void putSiStrip(edm::Event &e, edm::EventSetup const &iSetup)
collection_type::const_iterator const_iterator
Definition: DetSetVector.h:104
def move(src, dest)
Definition: eostools.py:510
data decode(const unsigned int &value) const
const SiStripApvGain::Range getRange(uint32_t detID) const
Definition: SiStripGain.h:66
const DetUnitContainer & detUnits() const
Returm a vector of all GeomDetUnit.