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  PreMixing_(conf.getParameter<bool>("PreMixingMode")),
55  theSiHitDigitizer(new SiHitDigitizer(conf)),
56  theSiPileUpSignals(new SiPileUpSignals()),
57  theSiNoiseAdder(new SiGaussianTailNoiseAdder(theThreshold)),
58  theSiDigitalConverter(new SiTrivialDigitalConverter(theElectronPerADC, PreMixing_)),
59  theSiZeroSuppress(new SiStripFedZeroSuppression(theFedAlgo)) {
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  unsigned int tofBin,
113  const StripGeomDetUnit* det,
114  const GlobalVector& bfield,
115  const TrackerTopology *tTopo,
116  CLHEP::HepRandomEngine* engine) {
117  // produce SignalPoints for all SimHits in detector
118  unsigned int detID = det->geographicalId().rawId();
119  int numStrips = (det->specificTopology()).nstrips();
120 
121  std::vector<bool>& badChannels = allBadChannels[detID];
122  size_t thisFirstChannelWithSignal = numStrips;
123  size_t thisLastChannelWithSignal = 0;
124 
125  float langle = (lorentzAngleHandle.isValid()) ? lorentzAngleHandle->getLorentzAngle(detID) : 0.;
126 
127  std::vector<float> locAmpl(numStrips, 0.);
128 
129  // Loop over hits
130 
131  uint32_t detId = det->geographicalId().rawId();
132  // First: loop on the SimHits
133  if(CLHEP::RandFlat::shoot(engine) > inefficiency) {
134  AssociationInfoForChannel* pDetIDAssociationInfo; // I only need this if makeDigiSimLinks_ is true...
135  if( makeDigiSimLinks_ ) pDetIDAssociationInfo=&(associationInfoForDetId_[detId]); // ...so only search the map if that is the case
136  std::vector<float> previousLocalAmplitude; // Only used if makeDigiSimLinks_ is true. Needed to work out the change in amplitude.
137 
138  size_t simHitGlobalIndex=inputBeginGlobalIndex; // This needs to stored to create the digi-sim link later
139  for (std::vector<PSimHit>::const_iterator simHitIter = inputBegin; simHitIter != inputEnd; ++simHitIter, ++simHitGlobalIndex ) {
140  // skip hits not in this detector.
141  if((*simHitIter).detUnitId() != detId) {
142  continue;
143  }
144  // check TOF
145  if (std::fabs(simHitIter->tof() - cosmicShift - det->surface().toGlobal(simHitIter->localPosition()).mag()/30.) < tofCut && simHitIter->energyLoss()>0) {
146  if( makeDigiSimLinks_ ) previousLocalAmplitude=locAmpl; // Not needed except to make the sim link association.
147  size_t localFirstChannel = numStrips;
148  size_t localLastChannel = 0;
149  // process the hit
150  theSiHitDigitizer->processHit(&*simHitIter, *det, bfield, langle, locAmpl, localFirstChannel, localLastChannel, tTopo, engine);
151 
152  //APV Killer to simulate HIP effect
153  //------------------------------------------------------
154 
156  int pdg_id = simHitIter->particleType();
157  particle = pdt->particle(pdg_id);
158  if(particle != NULL){
159  float charge = particle->charge();
160  bool isHadron = particle->isHadron();
161  if(charge!=0 && isHadron){
162  if(CLHEP::RandFlat::shoot(engine) < APVSaturationProb){
163  int FirstAPV = localFirstChannel/128;
164  int LastAPV = localLastChannel/128;
165  std::cout << "-------------------HIP--------------" << std::endl;
166  std::cout << "Killing APVs " << FirstAPV << " - " <<LastAPV << " " << detID <<std::endl;
167  for(int strip = FirstAPV*128; strip < LastAPV*128 +128; ++strip) {
168  badChannels[strip] = true;
169  }
170  //doing like that I remove the signal information only after the
171  //stip that got the HIP but it remains the signal of the previous
172  //one. I'll make a further loop to remove all signal
173  }
174  }
175  }
176  }
177 
178 
179  if(thisFirstChannelWithSignal > localFirstChannel) thisFirstChannelWithSignal = localFirstChannel;
180  if(thisLastChannelWithSignal < localLastChannel) thisLastChannelWithSignal = localLastChannel;
181 
182  if( makeDigiSimLinks_ ) { // No need to do any of this if truth association was turned off in the configuration
183  for( size_t stripIndex=0; stripIndex<locAmpl.size(); ++stripIndex ) {
184  // Work out the amplitude from this SimHit from the difference of what it was before and what it is now
185  float signalFromThisSimHit=locAmpl[stripIndex]-previousLocalAmplitude[stripIndex];
186  if( signalFromThisSimHit!=0 ) { // If this SimHit had any contribution I need to record it.
187  auto& associationVector=(*pDetIDAssociationInfo)[stripIndex];
188  bool addNewEntry=true;
189  // 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
190  // it's something to do with the stereo strips.
191  for( auto& associationInfo : associationVector ) {
192  if( associationInfo.trackID==simHitIter->trackId() && associationInfo.eventID==simHitIter->eventId() ) {
193  // The hit is already in, so add this second contribution and move on
194  associationInfo.contributionToADC+=signalFromThisSimHit;
195  addNewEntry=false;
196  break;
197  }
198  } // end of loop over associationVector
199  // If the hit wasn't already in create a new association info structure.
200  if( addNewEntry ) associationVector.push_back( AssociationInfo{ simHitIter->trackId(), simHitIter->eventId(), signalFromThisSimHit, simHitGlobalIndex, tofBin } );
201  } // end of "if( signalFromThisSimHit!=0 )"
202  } // end of loop over locAmpl strips
203  } // end of "if( makeDigiSimLinks_ )"
204  } // end of TOF check
205  } // end for
206  }
207  theSiPileUpSignals->add(detID, locAmpl, thisFirstChannelWithSignal, thisLastChannelWithSignal);
208 
209  if(firstChannelsWithSignal[detID] > thisFirstChannelWithSignal) firstChannelsWithSignal[detID] = thisFirstChannelWithSignal;
210  if(lastChannelsWithSignal[detID] < thisLastChannelWithSignal) lastChannelsWithSignal[detID] = thisLastChannelWithSignal;
211 }
212 
213 void
215  edm::DetSet<SiStripDigi>& outdigi,
216  edm::DetSet<SiStripRawDigi>& outrawdigi,
218  const StripGeomDetUnit *det,
219  edm::ESHandle<SiStripGain> & gainHandle,
220  edm::ESHandle<SiStripThreshold> & thresholdHandle,
221  edm::ESHandle<SiStripNoises> & noiseHandle,
222  edm::ESHandle<SiStripPedestals> & pedestalHandle,
223  CLHEP::HepRandomEngine* engine) {
224  unsigned int detID = det->geographicalId().rawId();
225  int numStrips = (det->specificTopology()).nstrips();
226 
227  const SiPileUpSignals::SignalMapType* theSignal(theSiPileUpSignals->getSignal(detID));
228 
229  std::vector<float> detAmpl(numStrips, 0.);
230  if(theSignal) {
231  for(const auto& amp : *theSignal) {
232  detAmpl[amp.first] = amp.second;
233  }
234  }
235 
236  //removing signal from the dead (and HIP effected) strips
237  std::vector<bool>& badChannels = allBadChannels[detID];
238  for(int strip =0; strip < numStrips; ++strip) if(badChannels[strip]) detAmpl[strip] = 0.;
239 
240  SiStripNoises::Range detNoiseRange = noiseHandle->getRange(detID);
241  SiStripApvGain::Range detGainRange = gainHandle->getRange(detID);
242  SiStripPedestals::Range detPedestalRange = pedestalHandle->getRange(detID);
243 
244 // -----------------------------------------------------------
245 
246  auto& firstChannelWithSignal = firstChannelsWithSignal[detID];
247  auto& lastChannelWithSignal = lastChannelsWithSignal[detID];
248  auto iAssociationInfoByChannel=associationInfoForDetId_.find(detID); // Use an iterator so that I can easily remove it once finished
249 
250  if(zeroSuppression){
251  if(noise){
252  int RefStrip = int(numStrips/2.);
253  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
254  RefStrip++;
255  }
256  if(RefStrip<numStrips){
257  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange);
258  float gainValue = gainHandle->getStripGain(RefStrip, detGainRange);
259  theSiNoiseAdder->addNoise(detAmpl,firstChannelWithSignal,lastChannelWithSignal,numStrips,noiseRMS*theElectronPerADC/gainValue, engine);
260  }
261  }
262  DigitalVecType digis;
263  theSiZeroSuppress->suppress(theSiDigitalConverter->convert(detAmpl, gainHandle, detID), digis, detID,noiseHandle,thresholdHandle);
264  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
265  // will be empty and the iterator will always equal associationInfoForDetId_.end().
266  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
267  for( const auto& iDigi : digis ) {
268  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
269  const std::vector<AssociationInfo>& associationInfo=associationInfoByChannel[iDigi.channel()];
270 
271  // Need to find the total from all sim hits, because this might not be the same as the total
272  // digitised due to noise or whatever.
273  float totalSimADC=0;
274  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
275  // Now I know that I can loop again and create the links
276  for( const auto& iAssociationInfo : associationInfo ) {
277  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
278  // expected counting to start at 1, not 0. Now changed.
279  outLink.push_back( StripDigiSimLink( iDigi.channel(), iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex, iAssociationInfo.tofBin, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
280  } // end of loop over associationInfo
281  } // end of loop over the digis
282  } // end of check that iAssociationInfoByChannel is a valid iterator
283  outdigi.data = digis;
284  }
285 
286 
287  if(!zeroSuppression){
288  //if(noise){
289  // the constant pedestal offset is needed because
290  // negative adc counts are not allowed in case
291  // Pedestal and CMN subtraction is performed.
292  // The pedestal value read from the conditions
293  // is pedValue and after the pedestal subtraction
294  // the baseline is zero. The Common Mode Noise
295  // is not subtracted from the negative adc counts
296  // channels. Adding pedOffset the baseline is set
297  // to pedOffset after pedestal subtraction and CMN
298  // is subtracted to all the channels since none of
299  // them has negative adc value. The pedOffset is
300  // treated as a constant component in the CMN
301  // estimation and subtracted as CMN.
302 
303 
304  //calculating the charge deposited on each APV and subtracting the shift
305  //------------------------------------------------------
306  if(BaselineShift){
307  theSiNoiseAdder->addBaselineShift(detAmpl, badChannels);
308  }
309 
310  //Adding the strip noise
311  //------------------------------------------------------
312  if(noise){
313  std::vector<float> noiseRMSv;
314  noiseRMSv.clear();
315  noiseRMSv.insert(noiseRMSv.begin(),numStrips,0.);
316 
317  if(SingleStripNoise){
318  for(int strip=0; strip< numStrips; ++strip){
319  if(!badChannels[strip]) noiseRMSv[strip] = (noiseHandle->getNoise(strip,detNoiseRange))* theElectronPerADC;
320  }
321 
322  } else {
323  int RefStrip = 0; //int(numStrips/2.);
324  while(RefStrip<numStrips&&badChannels[RefStrip]){ //if the refstrip is bad, I move up to when I don't find it
325  RefStrip++;
326  }
327  if(RefStrip<numStrips){
328  float noiseRMS = noiseHandle->getNoise(RefStrip,detNoiseRange) *theElectronPerADC;
329  for(int strip=0; strip< numStrips; ++strip){
330  if(!badChannels[strip]) noiseRMSv[strip] = noiseRMS;
331  }
332  }
333  }
334 
335  theSiNoiseAdder->addNoiseVR(detAmpl, noiseRMSv, engine);
336  }
337 
338  //adding the CMN
339  //------------------------------------------------------
340  if(CommonModeNoise){
341  float cmnRMS = 0.;
342  DetId detId(detID);
343  uint32_t SubDet = detId.subdetId();
344  if(SubDet==3){
345  cmnRMS = cmnRMStib;
346  }else if(SubDet==4){
347  cmnRMS = cmnRMStid;
348  }else if(SubDet==5){
349  cmnRMS = cmnRMStob;
350  }else if(SubDet==6){
351  cmnRMS = cmnRMStec;
352  }
353  cmnRMS *= theElectronPerADC;
354  theSiNoiseAdder->addCMNoise(detAmpl, cmnRMS, badChannels, engine);
355  }
356 
357 
358  //Adding the pedestals
359  //------------------------------------------------------
360 
361  std::vector<float> vPeds;
362  vPeds.clear();
363  vPeds.insert(vPeds.begin(),numStrips,0.);
364 
365  if(RealPedestals){
366  for(int strip=0; strip< numStrips; ++strip){
367  if(!badChannels[strip]) vPeds[strip] = (pedestalHandle->getPed(strip,detPedestalRange)+pedOffset)* theElectronPerADC;
368  }
369  } else {
370  for(int strip=0; strip< numStrips; ++strip){
371  if(!badChannels[strip]) vPeds[strip] = pedOffset* theElectronPerADC;
372  }
373  }
374 
375  theSiNoiseAdder->addPedestals(detAmpl, vPeds);
376 
377 
378  //if(!RealPedestals&&!CommonModeNoise&&!noise&&!BaselineShift&&!APVSaturationFromHIP){
379  // edm::LogWarning("SiStripDigitizer")<<"You are running the digitizer without Noise generation and without applying Zero Suppression. ARE YOU SURE???";
380  //}else{
381 
382  DigitalRawVecType rawdigis = theSiDigitalConverter->convertRaw(detAmpl, gainHandle, detID);
383 
384  // Now do the association to truth. Note that if truth association was turned off in the configuration this map
385  // will be empty and the iterator will always equal associationInfoForDetId_.end().
386  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) { // make sure the readings for this DetID aren't completely from noise
387  // N.B. For the raw digis the channel is inferred from the position in the vector.
388  // I'VE NOT TESTED THIS YET!!!!!
389  // ToDo Test this properly.
390  for( size_t channel=0; channel<rawdigis.size(); ++channel ) {
391  auto& associationInfoByChannel=iAssociationInfoByChannel->second;
392  const auto iAssociationInfo=associationInfoByChannel.find(channel);
393  if( iAssociationInfo==associationInfoByChannel.end() ) continue; // Skip if there is no sim information for this channel (i.e. it's all noise)
394  const std::vector<AssociationInfo>& associationInfo=iAssociationInfo->second;
395 
396  // Need to find the total from all sim hits, because this might not be the same as the total
397  // digitised due to noise or whatever.
398  float totalSimADC=0;
399  for( const auto& iAssociationInfo : associationInfo ) totalSimADC+=iAssociationInfo.contributionToADC;
400  // Now I know that I can loop again and create the links
401  for( const auto& iAssociationInfo : associationInfo ) {
402  // Note simHitGlobalIndex used to have +1 because TrackerHitAssociator (the only place I can find this value being used)
403  // expected counting to start at 1, not 0. Now changed.
404  outLink.push_back( StripDigiSimLink( channel, iAssociationInfo.trackID, iAssociationInfo.simHitGlobalIndex, iAssociationInfo.tofBin, iAssociationInfo.eventID, iAssociationInfo.contributionToADC/totalSimADC ) );
405  } // end of loop over associationInfo
406  } // end of loop over the digis
407  } // end of check that iAssociationInfoByChannel is a valid iterator
408 
409  outrawdigi.data = rawdigis;
410 
411  //}
412  }
413 
414  // Now that I've finished with this entry in the map of associations, I can remove it.
415  // Note that there might not be an association if the ADC reading is from noise in which
416  // case associationIsValid will be false.
417  if( iAssociationInfoByChannel!=associationInfoForDetId_.end() ) associationInfoForDetId_.erase(iAssociationInfoByChannel);
418 }
#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:40
const std::unique_ptr< const SiGaussianTailNoiseAdder > theSiNoiseAdder
const std::unique_ptr< SiHitDigitizer > theSiHitDigitizer
void getData(T &iHolder) const
Definition: EventSetup.h:78
std::map< unsigned int, size_t > firstChannelsWithSignal
std::map< int, Amplitude > SignalMapType
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
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:77
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
void accumulateSimHits(const std::vector< PSimHit >::const_iterator inputBegin, const std::vector< PSimHit >::const_iterator inputEnd, size_t inputBeginGlobalIndex, unsigned int tofBin, const StripGeomDetUnit *stripdet, const GlobalVector &bfield, const TrackerTopology *tTopo, CLHEP::HepRandomEngine *)
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:47
const ParticleDataTable * pdt