CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalNoiseInfoProducer.cc
Go to the documentation of this file.
1 //
2 // HcalNoiseInfoProducer.cc
3 //
4 // description: Implementation of the producer for the HCAL noise information
5 //
6 // author: J.P. Chou, Brown
7 //
8 //
9 
26 
27 using namespace reco;
28 
29 //
30 // constructors and destructor
31 //
32 
34 {
35  // set the parameters
36  fillDigis_ = iConfig.getParameter<bool>("fillDigis");
37  fillRecHits_ = iConfig.getParameter<bool>("fillRecHits");
38  fillCaloTowers_ = iConfig.getParameter<bool>("fillCaloTowers");
39  fillTracks_ = iConfig.getParameter<bool>("fillTracks");
40 
41  maxProblemRBXs_ = iConfig.getParameter<int>("maxProblemRBXs");
42 
43  maxCaloTowerIEta_ = iConfig.getParameter<int>("maxCaloTowerIEta");
44  maxTrackEta_ = iConfig.getParameter<double>("maxTrackEta");
45  minTrackPt_ = iConfig.getParameter<double>("minTrackPt");
46 
47  digiCollName_ = iConfig.getParameter<std::string>("digiCollName");
48  recHitCollName_ = iConfig.getParameter<std::string>("recHitCollName");
49  caloTowerCollName_ = iConfig.getParameter<std::string>("caloTowerCollName");
50  trackCollName_ = iConfig.getParameter<std::string>("trackCollName");
51 
52  minRecHitE_ = iConfig.getParameter<double>("minRecHitE");
53  minLowHitE_ = iConfig.getParameter<double>("minLowHitE");
54  minHighHitE_ = iConfig.getParameter<double>("minHighHitE");
55 
56  HcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("HcalAcceptSeverityLevel");
57  if (iConfig.exists("HcalRecHitFlagsToBeExcluded"))
58  HcalRecHitFlagsToBeExcluded_ = iConfig.getParameter<std::vector<int> >("HcalRecHitFlagsToBeExcluded");
59  else{
60  edm::LogWarning("MisConfiguration")<<"the module is missing the parameter HcalAcceptSeverityLevel. created empty.";
62  }
63 
64  // Digi threshold and time slices to use for HBHE and HF calibration digis
65  useCalibDigi_ = true;
66  if(iConfig.existsAs<double>("calibdigiHBHEthreshold") == false) useCalibDigi_ = false;
67  if(iConfig.existsAs<double>("calibdigiHFthreshold") == false) useCalibDigi_ = false;
68  if(iConfig.existsAs<std::vector<int> >("calibdigiHBHEtimeslices") == false) useCalibDigi_ = false;
69  if(iConfig.existsAs<std::vector<int> >("calibdigiHFtimeslices") == false) useCalibDigi_ = false;
70 
71  if(useCalibDigi_ == true)
72  {
73  calibdigiHBHEthreshold_ = iConfig.getParameter<double>("calibdigiHBHEthreshold");
74  calibdigiHBHEtimeslices_ = iConfig.getParameter<std::vector<int> >("calibdigiHBHEtimeslices");
75  calibdigiHFthreshold_ = iConfig.getParameter<double>("calibdigiHFthreshold");
76  calibdigiHFtimeslices_ = iConfig.getParameter<std::vector<int> >("calibdigiHFtimeslices");
77  }
78  else
79  {
81  calibdigiHBHEtimeslices_ = std::vector<int>();
83  calibdigiHFtimeslices_ = std::vector<int>();
84  }
85 
86  TS4TS5EnergyThreshold_ = iConfig.getParameter<double>("TS4TS5EnergyThreshold");
87 
88  std::vector<double> TS4TS5UpperThresholdTemp = iConfig.getParameter<std::vector<double> >("TS4TS5UpperThreshold");
89  std::vector<double> TS4TS5UpperCutTemp = iConfig.getParameter<std::vector<double> >("TS4TS5UpperCut");
90  std::vector<double> TS4TS5LowerThresholdTemp = iConfig.getParameter<std::vector<double> >("TS4TS5LowerThreshold");
91  std::vector<double> TS4TS5LowerCutTemp = iConfig.getParameter<std::vector<double> >("TS4TS5LowerCut");
92 
93  for(int i = 0; i < (int)TS4TS5UpperThresholdTemp.size() && i < (int)TS4TS5UpperCutTemp.size(); i++)
94  TS4TS5UpperCut_.push_back(std::pair<double, double>(TS4TS5UpperThresholdTemp[i], TS4TS5UpperCutTemp[i]));
95  sort(TS4TS5UpperCut_.begin(), TS4TS5UpperCut_.end());
96 
97  for(int i = 0; i < (int)TS4TS5LowerThresholdTemp.size() && i < (int)TS4TS5LowerCutTemp.size(); i++)
98  TS4TS5LowerCut_.push_back(std::pair<double, double>(TS4TS5LowerThresholdTemp[i], TS4TS5LowerCutTemp[i]));
99  sort(TS4TS5LowerCut_.begin(), TS4TS5LowerCut_.end());
100 
101  // if digis are filled, then rechits must also be filled
102  if(fillDigis_ && !fillRecHits_) {
103  fillRecHits_=true;
104  edm::LogWarning("HCalNoiseInfoProducer") << " forcing fillRecHits to be true if fillDigis is true.\n";
105  }
106 
107  const float adc2fCTemp[128]={-0.5,0.5,1.5,2.5,3.5,4.5,5.5,6.5,7.5,8.5,9.5,10.5,11.5,12.5,
108  13.5,15.,17.,19.,21.,23.,25.,27.,29.5,32.5,35.5,38.5,42.,46.,50.,54.5,59.5,
109  64.5,59.5,64.5,69.5,74.5,79.5,84.5,89.5,94.5,99.5,104.5,109.5,114.5,119.5,
110  124.5,129.5,137.,147.,157.,167.,177.,187.,197.,209.5,224.5,239.5,254.5,272.,
111  292.,312.,334.5,359.5,384.5,359.5,384.5,409.5,434.5,459.5,484.5,509.5,534.5,
112  559.5,584.5,609.5,634.5,659.5,684.5,709.5,747.,797.,847.,897.,947.,997.,
113  1047.,1109.5,1184.5,1259.5,1334.5,1422.,1522.,1622.,1734.5,1859.5,1984.5,
114  1859.5,1984.5,2109.5,2234.5,2359.5,2484.5,2609.5,2734.5,2859.5,2984.5,
115  3109.5,3234.5,3359.5,3484.5,3609.5,3797.,4047.,4297.,4547.,4797.,5047.,
116  5297.,5609.5,5984.5,6359.5,6734.5,7172.,7672.,8172.,8734.5,9359.5,9984.5};
117  for(int i = 0; i < 128; i++)
118  adc2fC[i] = adc2fCTemp[i];
119 
120  // we produce a vector of HcalNoiseRBXs
121  produces<HcalNoiseRBXCollection>();
122  // we also produce a noise summary
123  produces<HcalNoiseSummary>();
124 }
125 
126 
128 {
129 }
130 
131 
132 //
133 // member functions
134 //
135 
136 // ------------ method called to produce the data ------------
137 void
139 {
140  // this is what we're going to actually write to the EDM
141  std::auto_ptr<HcalNoiseRBXCollection> result1(new HcalNoiseRBXCollection);
142  std::auto_ptr<HcalNoiseSummary> result2(new HcalNoiseSummary);
143 
144  // define an empty HcalNoiseRBXArray that we're going to fill
145  HcalNoiseRBXArray rbxarray;
146  HcalNoiseSummary &summary=*result2;
147 
148  // fill them with the various components
149  // digi assumes that rechit information is available
150  if(fillRecHits_) fillrechits(iEvent, iSetup, rbxarray, summary);
151  if(fillDigis_) filldigis(iEvent, iSetup, rbxarray, summary);
152  if(fillCaloTowers_) fillcalotwrs(iEvent, iSetup, rbxarray, summary);
153  if(fillTracks_) filltracks(iEvent, iSetup, summary);
154 
155  // Why is this here? Shouldn't it have been in the filldigis method? Any reason for TotalCalibCharge to be defined outside filldigis(...) ?-- Jeff, 7/2/12
156  //if(fillDigis_) summary.calibCharge_ = TotalCalibCharge;
157 
158  // select those RBXs which are interesting
159  // also look for the highest energy RBX
160  HcalNoiseRBXArray::iterator maxit=rbxarray.begin();
161  double maxenergy=-999;
162  bool maxwritten=false;
163  for(HcalNoiseRBXArray::iterator rit = rbxarray.begin(); rit!=rbxarray.end(); ++rit) {
164  HcalNoiseRBX &rbx=(*rit);
167 
168  // find the highest energy rbx
169  if(data.energy()>maxenergy) {
170  maxenergy=data.energy();
171  maxit=rit;
172  maxwritten=false;
173  }
174 
175  // find out if the rbx is problematic/noisy/etc.
176  bool writerbx = algo_.isProblematic(data) || !algo_.passLooseNoiseFilter(data) ||
178 
179  // fill variables in the summary object not filled elsewhere
180  fillOtherSummaryVariables(summary, data);
181 
182  if(writerbx) {
183  summary.nproblemRBXs_++;
184  if(summary.nproblemRBXs_<=maxProblemRBXs_) {
185  result1->push_back(rbx);
186  if(maxit==rit) maxwritten=true;
187  }
188  }
189  } // end loop over rbxs
190 
191  // if we still haven't written the maximum energy rbx, write it now
192  if(!maxwritten) {
193  HcalNoiseRBX &rbx=(*maxit);
194 
195  // add the RBX to the event
196  result1->push_back(rbx);
197  }
198 
199  // put the rbxcollection and summary into the EDM
200  iEvent.put(result1);
201  iEvent.put(result2);
202 
203  return;
204 }
205 
206 
207 // ------------ here we fill specific variables in the summary object not already accounted for earlier
208 void
210 {
211  // charge ratio
212  if(algo_.passRatioThreshold(data) && data.validRatio()) {
213  if(data.ratio()<summary.minE2Over10TS()) {
214  summary.mine2ts_ = data.e2ts();
215  summary.mine10ts_ = data.e10ts(); }
216  if(data.ratio()>summary.maxE2Over10TS()) {
217  summary.maxe2ts_ = data.e2ts();
218  summary.maxe10ts_ = data.e10ts();
219  }
220  }
221 
222  // ADC zeros
223  if(algo_.passZerosThreshold(data)) {
224  if(data.numZeros()>summary.maxZeros()) {
225  summary.maxzeros_ = data.numZeros();
226  }
227  }
228 
229  // hits count
230  if(data.numHPDHits() > summary.maxHPDHits()) {
231  summary.maxhpdhits_ = data.numHPDHits();
232  }
233  if(data.numRBXHits() > summary.maxRBXHits()) {
234  summary.maxrbxhits_ = data.numRBXHits();
235  }
236  if(data.numHPDNoOtherHits() > summary.maxHPDNoOtherHits()) {
237  summary.maxhpdhitsnoother_ = data.numHPDNoOtherHits();
238  }
239 
240  // TS4TS5
241  if(data.PassTS4TS5() == false)
242  summary.hasBadRBXTS4TS5_ = true;
243 
244  // hit timing
245  if(data.minLowEHitTime()<summary.min10GeVHitTime()) {
246  summary.min10_ = data.minLowEHitTime();
247  }
248  if(data.maxLowEHitTime()>summary.max10GeVHitTime()) {
249  summary.max10_ = data.maxLowEHitTime();
250  }
251  summary.rms10_ += data.lowEHitTimeSqrd();
252  summary.cnthit10_ += data.numLowEHits();
253  if(data.minHighEHitTime()<summary.min25GeVHitTime()) {
254  summary.min25_ = data.minHighEHitTime();
255  }
256  if(data.maxHighEHitTime()>summary.max25GeVHitTime()) {
257  summary.max25_ = data.maxHighEHitTime();
258  }
259  summary.rms25_ += data.highEHitTimeSqrd();
260  summary.cnthit25_ += data.numHighEHits();
261 
262  // EMF
263  if(algo_.passEMFThreshold(data)) {
264  if(summary.minHPDEMF() > data.HPDEMF()) {
265  summary.minhpdemf_ = data.HPDEMF();
266  }
267  if(summary.minRBXEMF() > data.RBXEMF()) {
268  summary.minrbxemf_ = data.RBXEMF();
269  }
270  }
271 
272  // summary flag
273  if(!algo_.passLooseRatio(data)) summary.filterstatus_ |= 0x1;
274  if(!algo_.passLooseHits(data)) summary.filterstatus_ |= 0x2;
275  if(!algo_.passLooseZeros(data)) summary.filterstatus_ |= 0x4;
276  if(!algo_.passLooseTiming(data)) summary.filterstatus_ |= 0x8;
277 
278  if(!algo_.passTightRatio(data)) summary.filterstatus_ |= 0x100;
279  if(!algo_.passTightHits(data)) summary.filterstatus_ |= 0x200;
280  if(!algo_.passTightZeros(data)) summary.filterstatus_ |= 0x400;
281  if(!algo_.passTightTiming(data)) summary.filterstatus_ |= 0x800;
282 
283  if(!algo_.passHighLevelNoiseFilter(data)) summary.filterstatus_ |= 0x10000;
284 
285  // summary refvectors
287  if(!algo_.passLooseNoiseFilter(data))
288  join(summary.loosenoisetwrs_, data.rbxTowers());
289  if(!algo_.passTightNoiseFilter(data))
290  join(summary.tightnoisetwrs_, data.rbxTowers());
292  join(summary.hlnoisetwrs_, data.rbxTowers());
293 
294  return;
295 }
296 
297 
298 // ------------ fill the array with digi information
299 void
301 {
302  // Some initialization
303  TotalCalibCharge = 0;
304 
305  // Starting with this version (updated by Jeff Temple on Dec. 6, 2012), the "TS45" names in the variables are mis-nomers. The actual time slices used are determined from the digiTimeSlices_ variable, which may not be limited to only time slices 4 and 5. For now, "TS45" name kept, because that is what is used in HcalNoiseSummary object (in GetCalibCountTS45, etc.). Likewise, the charge value in 'gt15' is now configurable, though the name remains the same. For HBHE, we track both the number of calibration channels (NcalibTS45) and the number of calibration channels above threshold (NcalibTS45gt15). For HF, we track only the number of channels above the given threshold in the given time window (NcalibHFgtX). Default for HF in 2012 is to use the full time sample with effectively no threshold (threshold=-999)
306  int NcalibTS45=0;
307  int NcalibTS45gt15=0;
308  int NcalibHFgtX=0;
309 
310  double chargecalibTS45=0;
311  double chargecalibgt15TS45=0;
312 
313  // get the conditions and channel quality
314  edm::ESHandle<HcalDbService> conditions;
315  iSetup.get<HcalDbRecord>().get(conditions);
317  iSetup.get<HcalChannelQualityRcd>().get(qualhandle);
318  const HcalChannelQuality* myqual = qualhandle.product();
320  iSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
321  const HcalSeverityLevelComputer* mySeverity = mycomputer.product();
322 
323  // get the digis
325  iEvent.getByLabel(digiCollName_, handle);
326  if(!handle.isValid()) {
327  throw edm::Exception(edm::errors::ProductNotFound) << " could not find HBHEDigiCollection named " << digiCollName_ << "\n.";
328  return;
329  }
330 
331  // loop over all of the digi information
332  for(HBHEDigiCollection::const_iterator it=handle->begin(); it!=handle->end(); ++it) {
333  const HBHEDataFrame &digi=(*it);
334  HcalDetId cell = digi.id();
335  DetId detcell=(DetId)cell;
336 
337  // check on cells to be ignored and dropped
338  const HcalChannelStatus* mydigistatus=myqual->getValues(detcell.rawId());
339  if(mySeverity->dropChannel(mydigistatus->getValue())) continue;
340  if(digi.zsMarkAndPass()) continue;
341  // Drop if exclude bit set
342  if ((mydigistatus->getValue() & (1 <<HcalChannelStatus::HcalCellExcludeFromHBHENoiseSummary))==1) continue;
343 
344  // get the calibrations and coder
345  const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
346  const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
347  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
348  HcalCoderDb coder (*channelCoder, *shape);
349 
350  // match the digi to an rbx and hpd
351  HcalNoiseRBX &rbx=(*array.findRBX(digi));
352  HcalNoiseHPD &hpd=(*array.findHPD(digi));
353 
354  // determine if the digi is one the highest energy hits in the HPD
355  // this works because the rechits are sorted by energy (see fillrechits() below)
356  bool isBig=false, isBig5=false, isRBX=false;
357  int counter=0;
360  rit!=rechits.end(); ++rit, ++counter) {
361  if((*rit)->id() == digi.id()) {
362  if(counter==0) isBig=isBig5=true; // digi is also the highest energy rechit
363  if(counter<5) isBig5=true; // digi is one of 5 highest energy rechits
364  isRBX=true;
365  }
366  }
367 
368  // loop over each of the digi's time slices
369  int totalzeros=0;
370  CaloSamples tool;
371  coder.adc2fC(digi,tool);
372  for(int ts=0; ts<tool.size(); ++ts) {
373 
374  // count zero's
375  if(digi[ts].adc()==0) {
376  ++hpd.totalZeros_;
377  ++totalzeros;
378  }
379 
380  // get the fC's
381  double corrfc = tool[ts]-calibrations.pedestal(digi[ts].capid());
382 
383  // fill the relevant digi arrays
384  if(isBig) hpd.bigCharge_[ts]+=corrfc;
385  if(isBig5) hpd.big5Charge_[ts]+=corrfc;
386  if(isRBX) rbx.allCharge_[ts]+=corrfc;
387  }
388 
389  // record the maximum number of zero's found
390  if(totalzeros>hpd.maxZeros_)
391  hpd.maxZeros_=totalzeros;
392  }
393 
394  // get the calibration digis
396  iEvent.getByLabel("hcalDigis", hCalib);
397 
398  // get total charge in calibration channels
399  if(hCalib.isValid() == true)
400  {
401  for(HcalCalibDigiCollection::const_iterator digi = hCalib->begin(); digi != hCalib->end(); digi++)
402  {
403  if(digi->id().hcalSubdet() == 0)
404  continue;
405 
406  /*
407  HcalCalibDetId cell = digi->id();
408  DetId detcell = (DetId)cell;
409 
410  const HcalChannelStatus* mydigistatus=myqual->getValues(detcell.rawId());
411 
412  if(mySeverity->dropChannel(mydigistatus->getValue()))
413  continue;
414  if(digi->zsMarkAndPass())
415  continue;
416 
417  const HcalQIECoder *channelCoder = conditions->getHcalCoder(cell);
418  const HcalQIEShape* shape = conditions->getHcalShape(channelCoder);
419  HcalCoderDb coder(*channelCoder, *shape);
420 
421  CaloSamples tool;
422  coder.adc2fC(*digi, tool);
423 
424  for(int i = 0; i < (int)digi->size(); i++)
425  TotalCalibCharge = TotalCalibCharge + tool[i];
426  */
427 
428  // Original code computes total calib charge over all digis. While I think it would be more useful to skip
429  // zs mark-and-pass channels, I keep this computation as is. Individual HBHE and HF variables do skip
430  // the m-p channels. -- Jeff Temple, 6 December 2012
431 
432  for(int i = 0; i < (int)digi->size(); i++)
433  TotalCalibCharge = TotalCalibCharge + adc2fC[digi->sample(i).adc()&0xff];
434 
435 
436  HcalCalibDetId myid=(HcalCalibDetId)digi->id();
438  continue; // ignore HOCrosstalk channels
439  if(digi->zsMarkAndPass()) continue; // skip "mark-and-pass" channels when computing charge in calib channels
440 
441 
442  if (digi->id().hcalSubdet()==HcalForward) // check HF
443  {
444  double sumChargeHF=0;
445  for (unsigned int i=0;i<calibdigiHFtimeslices_.size();++i)
446  {
447  // skip unphysical time slices
448  if (calibdigiHFtimeslices_[i]<0 || calibdigiHFtimeslices_[i]>digi->size())
449  continue;
450  sumChargeHF+=adc2fC[digi->sample(calibdigiHFtimeslices_[i]).adc()&0xff];
451  }
452  if (sumChargeHF>calibdigiHFthreshold_) ++NcalibHFgtX;
453  } // end of HF check
454  else if (digi->id().hcalSubdet()==HcalBarrel || digi->id().hcalSubdet()==HcalEndcap) // now check HBHE
455  {
456  double sumChargeHBHE=0;
457  for (unsigned int i=0;i<calibdigiHBHEtimeslices_.size();++i)
458  {
459  // skip unphysical time slices
460  if (calibdigiHBHEtimeslices_[i]<0 || calibdigiHBHEtimeslices_[i]>digi->size())
461  continue;
462  sumChargeHBHE+=adc2fC[digi->sample(calibdigiHBHEtimeslices_[i]).adc()&0xff];
463  }
464  ++NcalibTS45;
465  chargecalibTS45+=sumChargeHBHE;
466  if (sumChargeHBHE>calibdigiHBHEthreshold_)
467  {
468  ++NcalibTS45gt15;
469  chargecalibgt15TS45+=sumChargeHBHE;
470  }
471  } // end of HBHE check
472  } // loop on HcalCalibDigiCollection
473  } // if (hCalib.isValid()==true)
474 
475  summary.calibCharge_ = TotalCalibCharge;
476  summary.calibCountTS45_=NcalibTS45;
477  summary.calibCountgt15TS45_=NcalibTS45gt15;
478  summary.calibChargeTS45_=chargecalibTS45;
479  summary.calibChargegt15TS45_=chargecalibgt15TS45;
480  summary.calibCountHF_=NcalibHFgtX;
481 
482  return;
483 }
484 
485 // ------------ fill the array with rec hit information
486 void
488 {
489  // get the HCAL channel status map
491  iSetup.get<HcalChannelQualityRcd>().get( hcalChStatus );
492  const HcalChannelQuality* dbHcalChStatus = hcalChStatus.product();
493 
494  // get the severity level computer
495  edm::ESHandle<HcalSeverityLevelComputer> hcalSevLvlComputerHndl;
496  iSetup.get<HcalSeverityLevelComputerRcd>().get(hcalSevLvlComputerHndl);
497  const HcalSeverityLevelComputer* hcalSevLvlComputer = hcalSevLvlComputerHndl.product();
498 
499  // get the calo geometry
501  iSetup.get<CaloGeometryRecord>().get(pG);
502  const CaloGeometry* geo = pG.product();
503 
504  // get the rechits
506  iEvent.getByLabel(recHitCollName_, handle);
507  if(!handle.isValid()) {
509  << " could not find HBHERecHitCollection named " << recHitCollName_ << "\n.";
510  return;
511  }
512 
513  summary.rechitCount_ = 0;
514  summary.rechitCount15_ = 0;
515  summary.rechitEnergy_ = 0;
516  summary.rechitEnergy15_ = 0;
517 
518  summary.hitsInLaserRegion_=0;
519  summary.hitsInNonLaserRegion_=0;
520  summary.energyInLaserRegion_=0;
521  summary.energyInNonLaserRegion_=0;
522 
523 
524 
525  // loop over all of the rechit information
526  for(HBHERecHitCollection::const_iterator it=handle->begin(); it!=handle->end(); ++it) {
527  const HBHERecHit &rechit=(*it);
528 
529  // skip bad rechits (other than those flagged by the isolated noise, triangle, flat, and spike algorithms)
530  const DetId id = rechit.detid();
531  uint32_t recHitFlag = rechit.flags();
532  uint32_t isolbitset = (1 << HcalCaloFlagLabels::HBHEIsolatedNoise);
533  uint32_t flatbitset = (1 << HcalCaloFlagLabels::HBHEFlatNoise);
534  uint32_t spikebitset = (1 << HcalCaloFlagLabels::HBHESpikeNoise);
535  uint32_t trianglebitset = (1 << HcalCaloFlagLabels::HBHETriangleNoise);
536  uint32_t ts4ts5bitset = (1 << HcalCaloFlagLabels::HBHETS4TS5Noise);
537  for(unsigned int i=0; i<HcalRecHitFlagsToBeExcluded_.size(); i++) {
538  uint32_t bitset = (1 << HcalRecHitFlagsToBeExcluded_[i]);
539  recHitFlag = (recHitFlag & bitset) ? recHitFlag-bitset : recHitFlag;
540  }
541  const uint32_t dbStatusFlag = dbHcalChStatus->getValues(id)->getValue();
542 
543  // Ignore rechit if exclude bit set, regardless of severity of other bits
544  if ((dbStatusFlag & (1 <<HcalChannelStatus::HcalCellExcludeFromHBHENoiseSummary))==1) continue;
545 
546  int severityLevel = hcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
547  bool isRecovered = hcalSevLvlComputer->recoveredRecHit(id, recHitFlag);
548  if(severityLevel!=0 && !isRecovered && severityLevel>static_cast<int>(HcalAcceptSeverityLevel_)) continue;
549 
550  // do some rechit counting and energies
551  summary.rechitCount_ = summary.rechitCount_ + 1;
552  summary.rechitEnergy_ = summary.rechitEnergy_ + rechit.energy();
553  if ((dbStatusFlag & (1 <<HcalChannelStatus::HcalBadLaserSignal))==1) // hit comes from a region where no laser calibration pulse is normally seen
554  {
555  ++summary.hitsInNonLaserRegion_;
556  summary.energyInNonLaserRegion_+=rechit.energy();
557  }
558  else // hit comes from region where laser calibration pulse is seen
559  {
560  ++summary.hitsInLaserRegion_;
561  summary.energyInLaserRegion_+=rechit.energy();
562  }
563 
564  if(rechit.energy() > 1.5)
565  {
566  summary.rechitCount15_ = summary.rechitCount15_ + 1;
567  summary.rechitEnergy15_ = summary.rechitEnergy15_ + rechit.energy();
568  }
569 
570  // if it was ID'd as isolated noise, update the summary object
571  if(rechit.flags() & isolbitset) {
572  summary.nisolnoise_++;
573  summary.isolnoisee_ += rechit.energy();
574  GlobalPoint gp = geo->getPosition(rechit.id());
575  double et = rechit.energy()*gp.perp()/gp.mag();
576  summary.isolnoiseet_ += et;
577  }
578 
579  if(rechit.flags() & flatbitset) {
580  summary.nflatnoise_++;
581  summary.flatnoisee_ += rechit.energy();
582  GlobalPoint gp = geo->getPosition(rechit.id());
583  double et = rechit.energy()*gp.perp()/gp.mag();
584  summary.flatnoiseet_ += et;
585  }
586 
587  if(rechit.flags() & spikebitset) {
588  summary.nspikenoise_++;
589  summary.spikenoisee_ += rechit.energy();
590  GlobalPoint gp = geo->getPosition(rechit.id());
591  double et = rechit.energy()*gp.perp()/gp.mag();
592  summary.spikenoiseet_ += et;
593  }
594 
595  if(rechit.flags() & trianglebitset) {
596  summary.ntrianglenoise_++;
597  summary.trianglenoisee_ += rechit.energy();
598  GlobalPoint gp = geo->getPosition(rechit.id());
599  double et = rechit.energy()*gp.perp()/gp.mag();
600  summary.trianglenoiseet_ += et;
601  }
602 
603  if(rechit.flags() & ts4ts5bitset) {
604  if ((dbStatusFlag & (1 <<HcalChannelStatus::HcalCellExcludeFromHBHENoiseSummaryR45))==0) // only add to TS4TS5 if the bit is not marked as "HcalCellExcludeFromHBHENoiseSummaryR45"
605  {
606  summary.nts4ts5noise_++;
607  summary.ts4ts5noisee_ += rechit.energy();
608  GlobalPoint gp = geo->getPosition(rechit.id());
609  double et = rechit.energy()*gp.perp()/gp.mag();
610  summary.ts4ts5noiseet_ += et;
611  }
612  }
613 
614  // find the hpd that the rechit is in
615  HcalNoiseHPD& hpd=(*array.findHPD(rechit));
616 
617  // create a persistent reference to the rechit
618  edm::Ref<HBHERecHitCollection> myRef(handle, it-handle->begin());
619 
620  // store it in a place so that it remains sorted by energy
621  hpd.refrechitset_.insert(myRef);
622 
623  } // end loop over rechits
624 
625  // loop over all HPDs and transfer the information from refrechitset_ to rechits_;
626  for(HcalNoiseRBXArray::iterator rbxit=array.begin(); rbxit!=array.end(); ++rbxit) {
627  for(std::vector<HcalNoiseHPD>::iterator hpdit=rbxit->hpds_.begin(); hpdit!=rbxit->hpds_.end(); ++hpdit) {
628 
629  // loop over all of the entries in the set and add them to rechits_
631  it=hpdit->refrechitset_.begin(); it!=hpdit->refrechitset_.end(); ++it) {
632  hpdit->rechits_.push_back(*it);
633  }
634  }
635  }
636  // now the rechits in all the HPDs are sorted by energy!
637 
638  return;
639 }
640 
641 // ------------ fill the array with calo tower information
642 void
644 {
645  // get the calotowers
647  iEvent.getByLabel(caloTowerCollName_, handle);
648  if(!handle.isValid()) {
650  << " could not find CaloTowerCollection named " << caloTowerCollName_ << "\n.";
651  return;
652  }
653 
654  summary.emenergy_ = summary.hadenergy_ = 0.0;
655 
656  // loop over all of the calotower information
657  for(CaloTowerCollection::const_iterator it = handle->begin(); it!=handle->end(); ++it) {
658  const CaloTower& twr=(*it);
659 
660  // create a persistent reference to the tower
661  edm::Ref<CaloTowerCollection> myRef(handle, it-handle->begin());
662 
663  // get all of the hpd's that are pointed to by the calotower
664  std::vector<std::vector<HcalNoiseHPD>::iterator> hpditervec;
665  array.findHPD(twr, hpditervec);
666 
667  // loop over the hpd's and add the reference to the RefVectors
668  for(std::vector<std::vector<HcalNoiseHPD>::iterator>::iterator it=hpditervec.begin();
669  it!=hpditervec.end(); ++it)
670  (*it)->calotowers_.push_back(myRef);
671 
672  // skip over anything with |ieta|>maxCaloTowerIEta
673  if(twr.ietaAbs()>maxCaloTowerIEta_) {
674  summary.emenergy_ += twr.emEnergy();
675  summary.hadenergy_ += twr.hadEnergy();
676  }
677  }
678 
679  return;
680 }
681 
682 // ------------ fill the summary with track information
683 void
685 {
687  iEvent.getByLabel(trackCollName_, handle);
688 
689  // don't throw exception, just return quietly
690  if(!handle.isValid()) {
691  // throw edm::Exception(edm::errors::ProductNotFound)
692  // << " could not find trackCollection named " << trackCollName_ << "\n.";
693  return;
694  }
695 
696  summary.trackenergy_=0.0;
697  for(reco::TrackCollection::const_iterator iTrack = handle->begin(); iTrack!=handle->end(); ++iTrack) {
698  reco::Track trk=*iTrack;
699  if(trk.pt()<minTrackPt_ || fabs(trk.eta())>maxTrackEta_) continue;
700 
701  summary.trackenergy_ += trk.p();
702  }
703 
704  return;
705 }
706 
707 //define this as a plug-in
int adc(sample_type sample)
get the ADC sample (12 bits)
double e10ts(void) const
Definition: HcalNoiseAlgo.h:23
double HPDEMF(void) const
Definition: HcalNoiseAlgo.h:38
double e2ts(void) const
Definition: HcalNoiseAlgo.h:22
double p() const
momentum vector magnitude
Definition: TrackBase.h:129
T getParameter(std::string const &) const
edm::RefVector< CaloTowerCollection > hlnoisetwrs_
int i
Definition: DBlmapReader.cc:9
double lowEHitTimeSqrd(void) const
Definition: HcalNoiseAlgo.h:31
std::vector< int > calibdigiHFtimeslices_
std::vector< HcalNoiseRBX > HcalNoiseRBXCollection
Definition: HcalNoiseRBX.h:26
std::vector< int > HcalRecHitFlagsToBeExcluded_
bool isProblematic(const CommonHcalNoiseRBXData &) const
bool passTightRatio(const CommonHcalNoiseRBXData &) const
bool zsMarkAndPass() const
was ZS MarkAndPass?
Definition: HBHEDataFrame.h:30
double minLowEHitTime(void) const
Definition: HcalNoiseAlgo.h:29
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:187
bool passLooseRatio(const CommonHcalNoiseRBXData &) const
CalibDetType calibFlavor() const
get the flavor of this calibration detid
void fillrechits(edm::Event &, const edm::EventSetup &, HcalNoiseRBXArray &, HcalNoiseSummary &) const
const DetId & detid() const
Definition: CaloRecHit.h:22
double maxHighEHitTime(void) const
Definition: HcalNoiseAlgo.h:34
bool passTightTiming(const CommonHcalNoiseRBXData &) const
double maxLowEHitTime(void) const
Definition: HcalNoiseAlgo.h:30
int numLowEHits(void) const
Definition: HcalNoiseAlgo.h:32
virtual void produce(edm::Event &, const edm::EventSetup &) override
double ratio(void) const
Definition: HcalNoiseAlgo.h:21
std::vector< std::pair< double, double > > TS4TS5LowerCut_
HcalDetId id() const
Definition: HBHERecHit.h:25
std::set< edm::Ref< HBHERecHitCollection >, RefHBHERecHitEnergyComparison > refrechitset_
Definition: HcalNoiseHPD.h:142
std::vector< T >::const_iterator const_iterator
std::vector< float > bigCharge_
Definition: HcalNoiseHPD.h:134
float maxE2Over10TS(void) const
edm::RefVector< CaloTowerCollection > loosenoisetwrs_
bool exists(std::string const &parameterName) const
checks if a parameter exists
void fillcalotwrs(edm::Event &, const edm::EventSetup &, HcalNoiseRBXArray &, HcalNoiseSummary &) const
std::vector< float > big5Charge_
Definition: HcalNoiseHPD.h:135
double pedestal(int fCapId) const
get pedestal for capid=0..3
const Item * getValues(DetId fId, bool throwOnFail=true) const
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
float min25GeVHitTime(void) const
bool passHighLevelNoiseFilter(const CommonHcalNoiseRBXData &) const
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
double RBXEMF(void) const
Definition: HcalNoiseAlgo.h:37
edm::RefVector< HBHERecHitCollection > rechits_
Definition: HcalNoiseHPD.h:138
int numRBXHits(void) const
Definition: HcalNoiseAlgo.h:26
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
float min10GeVHitTime(void) const
HcalNoiseRBXArray::iterator findRBX(int rbxindex)
double highEHitTimeSqrd(void) const
Definition: HcalNoiseAlgo.h:35
int iEvent
Definition: GenABIO.cc:243
float minE2Over10TS(void) const
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:141
bool passLooseZeros(const CommonHcalNoiseRBXData &) const
double emEnergy() const
Definition: CaloTower.h:79
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
int maxRBXHits(void) const
float energy() const
Definition: CaloRecHit.h:19
virtual void adc2fC(const HBHEDataFrame &df, CaloSamples &lf) const
Definition: HcalCoderDb.cc:46
bool passTightZeros(const CommonHcalNoiseRBXData &) const
double pt() const
track transverse momentum
Definition: TrackBase.h:131
bool recoveredRecHit(const DetId &myid, const uint32_t &myflag) const
uint32_t flags() const
Definition: CaloRecHit.h:23
int numHighEHits(void) const
Definition: HcalNoiseAlgo.h:36
tuple handle
Definition: patZpeak.py:22
bool passRatioThreshold(const CommonHcalNoiseRBXData &) const
std::vector< HcalNoiseHPD >::iterator findHPD(int hpdindex)
bool dropChannel(const uint32_t &mystatus) const
const GlobalPoint & getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:68
int numHPDNoOtherHits(void) const
Definition: HcalNoiseAlgo.h:27
int numZeros(void) const
Definition: HcalNoiseAlgo.h:28
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
bool passZerosThreshold(const CommonHcalNoiseRBXData &) const
std::vector< std::pair< double, double > > TS4TS5UpperCut_
float minHPDEMF(void) const
bool passLooseNoiseFilter(const CommonHcalNoiseRBXData &) const
edm::RefVector< CaloTowerCollection > rbxTowers(void) const
Definition: HcalNoiseAlgo.h:40
double hadEnergy() const
Definition: CaloTower.h:80
std::vector< float > allCharge_
Definition: HcalNoiseRBX.h:105
DEFINE_FWK_MODULE(CosmicTrackingParticleSelector)
Definition: DetId.h:20
bool passLooseHits(const CommonHcalNoiseRBXData &) const
static std::string join(char **cmd)
Definition: RemoteFile.cc:18
int size() const
get the size
Definition: CaloSamples.h:26
void filldigis(edm::Event &, const edm::EventSetup &, HcalNoiseRBXArray &, HcalNoiseSummary &)
int maxZeros(void) const
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
edm::RefVector< CaloTowerCollection > tightnoisetwrs_
int getSeverityLevel(const DetId &myid, const uint32_t &myflag, const uint32_t &mystatus) const
double energy(void) const
Definition: HcalNoiseAlgo.h:20
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:18
int ietaAbs() const
Definition: CaloTower.h:155
bool passLooseTiming(const CommonHcalNoiseRBXData &) const
bool passTightNoiseFilter(const CommonHcalNoiseRBXData &) const
int maxHPDNoOtherHits(void) const
char data[epos_bytes_allocation]
Definition: EPOS_Wrapper.h:82
double minHighEHitTime(void) const
Definition: HcalNoiseAlgo.h:33
double energyInNonLaserRegion_
float minRBXEMF(void) const
std::vector< int > calibdigiHBHEtimeslices_
const HcalDetId & id() const
Definition: HBHEDataFrame.h:22
bool PassTS4TS5(void) const
Definition: HcalNoiseAlgo.h:39
uint32_t getValue() const
bool passEMFThreshold(const CommonHcalNoiseRBXData &) const
void filltracks(edm::Event &, const edm::EventSetup &, HcalNoiseSummary &) const
HcalNoiseInfoProducer(const edm::ParameterSet &)
void fillOtherSummaryVariables(HcalNoiseSummary &summary, const CommonHcalNoiseRBXData &data) const
int maxHPDHits(void) const
int numHPDHits(void) const
Definition: HcalNoiseAlgo.h:25
float max25GeVHitTime(void) const
float max10GeVHitTime(void) const
bool validRatio(void) const
Definition: HcalNoiseAlgo.h:24
bool passTightHits(const CommonHcalNoiseRBXData &) const