CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/RecoMET/METProducers/src/HcalNoiseInfoProducer.cc

Go to the documentation of this file.
00001 //
00002 // HcalNoiseInfoProducer.cc
00003 //
00004 //   description: Implementation of the producer for the HCAL noise information
00005 //
00006 //   author: J.P. Chou, Brown
00007 //
00008 //
00009 
00010 #include "RecoMET/METProducers/interface/HcalNoiseInfoProducer.h"
00011 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00012 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00013 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00014 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00015 #include "DataFormats/TrackReco/interface/Track.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017 #include "CalibFormats/HcalObjects/interface/HcalCoderDb.h"
00018 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00019 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00020 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00021 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
00022 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
00023 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
00024 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00025 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00026 
00027 using namespace reco;
00028 
00029 //
00030 // constructors and destructor
00031 //
00032 
00033 HcalNoiseInfoProducer::HcalNoiseInfoProducer(const edm::ParameterSet& iConfig) : algo_(iConfig)
00034 {
00035   // set the parameters
00036   fillDigis_         = iConfig.getParameter<bool>("fillDigis");
00037   fillRecHits_       = iConfig.getParameter<bool>("fillRecHits");
00038   fillCaloTowers_    = iConfig.getParameter<bool>("fillCaloTowers");
00039   fillTracks_        = iConfig.getParameter<bool>("fillTracks");
00040 
00041   maxProblemRBXs_    = iConfig.getParameter<int>("maxProblemRBXs");
00042 
00043   maxCaloTowerIEta_  = iConfig.getParameter<int>("maxCaloTowerIEta");
00044   maxTrackEta_       = iConfig.getParameter<double>("maxTrackEta");
00045   minTrackPt_        = iConfig.getParameter<double>("minTrackPt");
00046 
00047   digiCollName_      = iConfig.getParameter<std::string>("digiCollName");
00048   recHitCollName_    = iConfig.getParameter<std::string>("recHitCollName");
00049   caloTowerCollName_ = iConfig.getParameter<std::string>("caloTowerCollName");
00050   trackCollName_     = iConfig.getParameter<std::string>("trackCollName");
00051 
00052   minRecHitE_        = iConfig.getParameter<double>("minRecHitE");
00053   minLowHitE_        = iConfig.getParameter<double>("minLowHitE");
00054   minHighHitE_       = iConfig.getParameter<double>("minHighHitE");
00055 
00056   HcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("HcalAcceptSeverityLevel");
00057   if (iConfig.exists("HcalRecHitFlagsToBeExcluded"))
00058       HcalRecHitFlagsToBeExcluded_ = iConfig.getParameter<std::vector<int> >("HcalRecHitFlagsToBeExcluded");
00059   else{
00060     edm::LogWarning("MisConfiguration")<<"the module is missing the parameter HcalAcceptSeverityLevel. created empty.";
00061     HcalRecHitFlagsToBeExcluded_.resize(0);
00062   }
00063 
00064   TS4TS5EnergyThreshold_ = iConfig.getParameter<double>("TS4TS5EnergyThreshold");
00065 
00066   std::vector<double> TS4TS5UpperThresholdTemp = iConfig.getParameter<std::vector<double> >("TS4TS5UpperThreshold");
00067   std::vector<double> TS4TS5UpperCutTemp = iConfig.getParameter<std::vector<double> >("TS4TS5UpperCut");
00068   std::vector<double> TS4TS5LowerThresholdTemp = iConfig.getParameter<std::vector<double> >("TS4TS5LowerThreshold");
00069   std::vector<double> TS4TS5LowerCutTemp = iConfig.getParameter<std::vector<double> >("TS4TS5LowerCut");
00070 
00071   for(int i = 0; i < (int)TS4TS5UpperThresholdTemp.size() && i < (int)TS4TS5UpperCutTemp.size(); i++)
00072      TS4TS5UpperCut_.push_back(std::pair<double, double>(TS4TS5UpperThresholdTemp[i], TS4TS5UpperCutTemp[i]));
00073   sort(TS4TS5UpperCut_.begin(), TS4TS5UpperCut_.end());
00074 
00075   for(int i = 0; i < (int)TS4TS5LowerThresholdTemp.size() && i < (int)TS4TS5LowerCutTemp.size(); i++)
00076      TS4TS5LowerCut_.push_back(std::pair<double, double>(TS4TS5LowerThresholdTemp[i], TS4TS5LowerCutTemp[i]));
00077   sort(TS4TS5LowerCut_.begin(), TS4TS5LowerCut_.end());
00078 
00079   // if digis are filled, then rechits must also be filled
00080   if(fillDigis_ && !fillRecHits_) {
00081     fillRecHits_=true;
00082     edm::LogWarning("HCalNoiseInfoProducer") << " forcing fillRecHits to be true if fillDigis is true.\n";
00083   }
00084 
00085   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,
00086      13.5,15.,17.,19.,21.,23.,25.,27.,29.5,32.5,35.5,38.5,42.,46.,50.,54.5,59.5,
00087      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,
00088      124.5,129.5,137.,147.,157.,167.,177.,187.,197.,209.5,224.5,239.5,254.5,272.,
00089      292.,312.,334.5,359.5,384.5,359.5,384.5,409.5,434.5,459.5,484.5,509.5,534.5,
00090      559.5,584.5,609.5,634.5,659.5,684.5,709.5,747.,797.,847.,897.,947.,997.,
00091      1047.,1109.5,1184.5,1259.5,1334.5,1422.,1522.,1622.,1734.5,1859.5,1984.5,
00092      1859.5,1984.5,2109.5,2234.5,2359.5,2484.5,2609.5,2734.5,2859.5,2984.5,
00093      3109.5,3234.5,3359.5,3484.5,3609.5,3797.,4047.,4297.,4547.,4797.,5047.,
00094      5297.,5609.5,5984.5,6359.5,6734.5,7172.,7672.,8172.,8734.5,9359.5,9984.5};
00095   for(int i = 0; i < 128; i++)
00096      adc2fC[i] = adc2fCTemp[i];
00097 
00098   // we produce a vector of HcalNoiseRBXs
00099   produces<HcalNoiseRBXCollection>();
00100   // we also produce a noise summary
00101   produces<HcalNoiseSummary>();
00102 }
00103 
00104 
00105 HcalNoiseInfoProducer::~HcalNoiseInfoProducer()
00106 {
00107 }
00108 
00109 
00110 //
00111 // member functions
00112 //
00113 
00114 // ------------ method called to produce the data  ------------
00115 void
00116 HcalNoiseInfoProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00117 {
00118   // this is what we're going to actually write to the EDM
00119   std::auto_ptr<HcalNoiseRBXCollection> result1(new HcalNoiseRBXCollection);
00120   std::auto_ptr<HcalNoiseSummary> result2(new HcalNoiseSummary);
00121 
00122   // define an empty HcalNoiseRBXArray that we're going to fill
00123   HcalNoiseRBXArray rbxarray;
00124   HcalNoiseSummary &summary=*result2;
00125 
00126   // fill them with the various components
00127   // digi assumes that rechit information is available
00128   if(fillRecHits_)    fillrechits(iEvent, iSetup, rbxarray, summary);
00129   if(fillDigis_)      filldigis(iEvent, iSetup, rbxarray);
00130   if(fillCaloTowers_) fillcalotwrs(iEvent, iSetup, rbxarray, summary);
00131   if(fillTracks_)     filltracks(iEvent, iSetup, summary);
00132 
00133   if(fillDigis_)      summary.calibCharge_ = TotalCalibCharge;
00134 
00135   // select those RBXs which are interesting
00136   // also look for the highest energy RBX
00137   HcalNoiseRBXArray::iterator maxit=rbxarray.begin();
00138   double maxenergy=-999;
00139   bool maxwritten=false;
00140   for(HcalNoiseRBXArray::iterator rit = rbxarray.begin(); rit!=rbxarray.end(); ++rit) {
00141     HcalNoiseRBX &rbx=(*rit);
00142     CommonHcalNoiseRBXData data(rbx, minRecHitE_, minLowHitE_, minHighHitE_, TS4TS5EnergyThreshold_,
00143       TS4TS5UpperCut_, TS4TS5LowerCut_);
00144 
00145     // find the highest energy rbx
00146     if(data.energy()>maxenergy) {
00147       maxenergy=data.energy();
00148       maxit=rit;
00149       maxwritten=false;
00150     }
00151 
00152     // find out if the rbx is problematic/noisy/etc.
00153     bool writerbx = algo_.isProblematic(data) || !algo_.passLooseNoiseFilter(data) ||
00154       !algo_.passTightNoiseFilter(data) || !algo_.passHighLevelNoiseFilter(data);
00155 
00156     // fill variables in the summary object not filled elsewhere
00157     fillOtherSummaryVariables(summary, data);
00158 
00159     if(writerbx) {
00160       summary.nproblemRBXs_++;
00161       if(summary.nproblemRBXs_<=maxProblemRBXs_) {
00162         result1->push_back(rbx);
00163         if(maxit==rit) maxwritten=true;
00164       }
00165     }
00166   } // end loop over rbxs
00167 
00168   // if we still haven't written the maximum energy rbx, write it now
00169   if(!maxwritten) {
00170     HcalNoiseRBX &rbx=(*maxit);
00171 
00172     // add the RBX to the event
00173     result1->push_back(rbx);
00174   }
00175 
00176   // put the rbxcollection and summary into the EDM
00177   iEvent.put(result1);
00178   iEvent.put(result2);
00179 
00180   return;
00181 }
00182 
00183 // ------------ method called once each job just before starting event loop  ------------
00184 void
00185 HcalNoiseInfoProducer::beginJob()
00186 {
00187   return;
00188 }
00189 
00190 // ------------ method called once each job just after ending the event loop  ------------
00191 void
00192 HcalNoiseInfoProducer::endJob()
00193 {
00194   return;
00195 }
00196 
00197 
00198 // ------------ method called once each run just before starting event loop  ------------
00199 // ------------ fills the pedestals
00200 void
00201 HcalNoiseInfoProducer::beginRun(edm::Run&, const edm::EventSetup&)
00202 {
00203   return;
00204 }
00205 
00206 // ------------ method called once each job just after ending the event loop  ------------
00207 void
00208 HcalNoiseInfoProducer::endRun(edm::Run&, const edm::EventSetup&)
00209 {
00210   return;
00211 }
00212 
00213 // ------------ here we fill specific variables in the summary object not already accounted for earlier
00214 void
00215 HcalNoiseInfoProducer::fillOtherSummaryVariables(HcalNoiseSummary& summary, const CommonHcalNoiseRBXData& data) const
00216 {
00217   // charge ratio
00218   if(algo_.passRatioThreshold(data) && data.validRatio()) {
00219     if(data.ratio()<summary.minE2Over10TS()) {
00220       summary.mine2ts_ = data.e2ts();
00221       summary.mine10ts_ = data.e10ts();    }
00222     if(data.ratio()>summary.maxE2Over10TS()) {
00223       summary.maxe2ts_ = data.e2ts();
00224       summary.maxe10ts_ = data.e10ts();
00225     }
00226   }
00227 
00228   // ADC zeros
00229   if(algo_.passZerosThreshold(data)) {
00230     if(data.numZeros()>summary.maxZeros()) {
00231       summary.maxzeros_ = data.numZeros();
00232     }
00233   }
00234 
00235   // hits count
00236   if(data.numHPDHits() > summary.maxHPDHits()) {
00237     summary.maxhpdhits_ = data.numHPDHits();
00238   }
00239   if(data.numRBXHits() > summary.maxRBXHits()) {
00240     summary.maxrbxhits_ = data.numRBXHits();
00241   }
00242   if(data.numHPDNoOtherHits() > summary.maxHPDNoOtherHits()) {
00243     summary.maxhpdhitsnoother_ = data.numHPDNoOtherHits();
00244   }
00245 
00246   // TS4TS5
00247   if(data.PassTS4TS5() == false)
00248      summary.hasBadRBXTS4TS5_ = true;
00249 
00250   // hit timing
00251   if(data.minLowEHitTime()<summary.min10GeVHitTime()) {
00252     summary.min10_ = data.minLowEHitTime();
00253   }
00254   if(data.maxLowEHitTime()>summary.max10GeVHitTime()) {
00255     summary.max10_ = data.maxLowEHitTime();
00256   }
00257   summary.rms10_ += data.lowEHitTimeSqrd();
00258   summary.cnthit10_ += data.numLowEHits();
00259   if(data.minHighEHitTime()<summary.min25GeVHitTime()) {
00260     summary.min25_ = data.minHighEHitTime();
00261   }
00262   if(data.maxHighEHitTime()>summary.max25GeVHitTime()) {
00263     summary.max25_ = data.maxHighEHitTime();
00264   }
00265   summary.rms25_ += data.highEHitTimeSqrd();
00266   summary.cnthit25_ += data.numHighEHits();
00267 
00268   // EMF
00269   if(algo_.passEMFThreshold(data)) {
00270     if(summary.minHPDEMF() > data.HPDEMF()) {
00271       summary.minhpdemf_ = data.HPDEMF();
00272     }
00273     if(summary.minRBXEMF() > data.RBXEMF()) {
00274       summary.minrbxemf_ = data.RBXEMF();
00275     }
00276   }
00277 
00278   // summary flag
00279   if(!algo_.passLooseRatio(data))  summary.filterstatus_ |= 0x1;
00280   if(!algo_.passLooseHits(data))   summary.filterstatus_ |= 0x2;
00281   if(!algo_.passLooseZeros(data))  summary.filterstatus_ |= 0x4;
00282   if(!algo_.passLooseTiming(data)) summary.filterstatus_ |= 0x8;
00283 
00284   if(!algo_.passTightRatio(data))  summary.filterstatus_ |= 0x100;
00285   if(!algo_.passTightHits(data))   summary.filterstatus_ |= 0x200;
00286   if(!algo_.passTightZeros(data))  summary.filterstatus_ |= 0x400;
00287   if(!algo_.passTightTiming(data)) summary.filterstatus_ |= 0x800;
00288 
00289   if(!algo_.passHighLevelNoiseFilter(data)) summary.filterstatus_ |= 0x10000;
00290 
00291   // summary refvectors
00292   JoinCaloTowerRefVectorsWithoutDuplicates join;
00293   if(!algo_.passLooseNoiseFilter(data))
00294     join(summary.loosenoisetwrs_, data.rbxTowers());
00295   if(!algo_.passTightNoiseFilter(data))
00296     join(summary.tightnoisetwrs_, data.rbxTowers());
00297   if(!algo_.passHighLevelNoiseFilter(data))
00298     join(summary.hlnoisetwrs_, data.rbxTowers());
00299 
00300   return;
00301 }
00302 
00303 
00304 // ------------ fill the array with digi information
00305 void
00306 HcalNoiseInfoProducer::filldigis(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseRBXArray& array)
00307 {
00308   // Some initialization
00309   TotalCalibCharge = 0;
00310 
00311   // get the conditions and channel quality
00312   edm::ESHandle<HcalDbService> conditions;
00313   iSetup.get<HcalDbRecord>().get(conditions);
00314   const HcalQIEShape* shape = conditions->getHcalShape();
00315   edm::ESHandle<HcalChannelQuality> qualhandle;
00316   iSetup.get<HcalChannelQualityRcd>().get(qualhandle);
00317   const HcalChannelQuality* myqual = qualhandle.product();
00318   edm::ESHandle<HcalSeverityLevelComputer> mycomputer;
00319   iSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
00320   const HcalSeverityLevelComputer* mySeverity = mycomputer.product();
00321 
00322   // get the digis
00323   edm::Handle<HBHEDigiCollection> handle;
00324   iEvent.getByLabel(digiCollName_, handle);
00325   if(!handle.isValid()) {
00326     throw edm::Exception(edm::errors::ProductNotFound) << " could not find HBHEDigiCollection named " << digiCollName_ << "\n.";
00327     return;
00328   }
00329 
00330   // loop over all of the digi information
00331   for(HBHEDigiCollection::const_iterator it=handle->begin(); it!=handle->end(); ++it) {
00332     const HBHEDataFrame &digi=(*it);
00333     HcalDetId cell = digi.id();
00334     DetId detcell=(DetId)cell;
00335 
00336     // check on cells to be ignored and dropped
00337     const HcalChannelStatus* mydigistatus=myqual->getValues(detcell.rawId());
00338     if(mySeverity->dropChannel(mydigistatus->getValue())) continue;
00339     if(digi.zsMarkAndPass()) continue;
00340 
00341     // get the calibrations and coder
00342     const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
00343     const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00344     HcalCoderDb coder (*channelCoder, *shape);
00345 
00346     // match the digi to an rbx and hpd
00347     HcalNoiseRBX &rbx=(*array.findRBX(digi));
00348     HcalNoiseHPD &hpd=(*array.findHPD(digi));
00349 
00350     // determine if the digi is one the highest energy hits in the HPD
00351     // this works because the rechits are sorted by energy (see fillrechits() below)
00352     bool isBig=false, isBig5=false, isRBX=false;
00353     int counter=0;
00354     edm::RefVector<HBHERecHitCollection> &rechits=hpd.rechits_;
00355     for(edm::RefVector<HBHERecHitCollection>::const_iterator rit=rechits.begin();
00356         rit!=rechits.end(); ++rit, ++counter) {
00357       if((*rit)->id() == digi.id()) {
00358         if(counter==0) isBig=isBig5=true;  // digi is also the highest energy rechit
00359         if(counter<5) isBig5=true;         // digi is one of 5 highest energy rechits
00360         isRBX=true;
00361       }
00362     }
00363 
00364     // loop over each of the digi's time slices
00365     int totalzeros=0;
00366     CaloSamples tool;
00367     coder.adc2fC(digi,tool);
00368     for(int ts=0; ts<tool.size(); ++ts) {
00369 
00370       // count zero's
00371       if(digi[ts].adc()==0) {
00372         ++hpd.totalZeros_;
00373         ++totalzeros;
00374       }
00375 
00376       // get the fC's
00377       double corrfc = tool[ts]-calibrations.pedestal(digi[ts].capid());
00378 
00379       // fill the relevant digi arrays
00380       if(isBig)  hpd.bigCharge_[ts]+=corrfc;
00381       if(isBig5) hpd.big5Charge_[ts]+=corrfc;
00382       if(isRBX)  rbx.allCharge_[ts]+=corrfc;
00383     }
00384 
00385     // record the maximum number of zero's found
00386     if(totalzeros>hpd.maxZeros_)
00387       hpd.maxZeros_=totalzeros;
00388   }
00389 
00390   // get the calibration digis
00391   edm::Handle<HcalCalibDigiCollection> hCalib;
00392   iEvent.getByLabel("hcalDigis", hCalib);
00393 
00394   // get total charge in calibration channels
00395   if(hCalib.isValid() == true)
00396   {
00397      for(HcalCalibDigiCollection::const_iterator digi = hCalib->begin(); digi != hCalib->end(); digi++)
00398      {
00399         if(digi->id().hcalSubdet() == 0)
00400            continue;
00401 
00402         /*
00403         HcalCalibDetId cell = digi->id();
00404         DetId detcell = (DetId)cell;
00405         
00406         const HcalChannelStatus* mydigistatus=myqual->getValues(detcell.rawId());
00407 
00408         if(mySeverity->dropChannel(mydigistatus->getValue()))
00409            continue;
00410         if(digi->zsMarkAndPass())
00411            continue;
00412 
00413         const HcalQIECoder *channelCoder = conditions->getHcalCoder(cell);
00414         HcalCoderDb coder(*channelCoder, *shape);
00415 
00416         CaloSamples tool;
00417         coder.adc2fC(*digi, tool);
00418 
00419         for(int i = 0; i < (int)digi->size(); i++)
00420            TotalCalibCharge = TotalCalibCharge + tool[i];
00421         */
00422 
00423         for(int i = 0; i < (int)digi->size(); i++)
00424            TotalCalibCharge = TotalCalibCharge + adc2fC[digi->sample(i).adc()&0xff];
00425      }
00426   }
00427 
00428   return;
00429 }
00430 
00431 // ------------ fill the array with rec hit information
00432 void
00433 HcalNoiseInfoProducer::fillrechits(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseRBXArray& array, HcalNoiseSummary& summary) const
00434 {
00435   // get the HCAL channel status map
00436   edm::ESHandle<HcalChannelQuality> hcalChStatus;
00437   iSetup.get<HcalChannelQualityRcd>().get( hcalChStatus );
00438   const HcalChannelQuality* dbHcalChStatus = hcalChStatus.product();
00439 
00440   // get the severity level computer
00441   edm::ESHandle<HcalSeverityLevelComputer> hcalSevLvlComputerHndl;
00442   iSetup.get<HcalSeverityLevelComputerRcd>().get(hcalSevLvlComputerHndl);
00443   const HcalSeverityLevelComputer* hcalSevLvlComputer = hcalSevLvlComputerHndl.product();
00444 
00445   // get the calo geometry
00446   edm::ESHandle<CaloGeometry> pG;
00447   iSetup.get<CaloGeometryRecord>().get(pG);
00448   const CaloGeometry* geo = pG.product();
00449 
00450   // get the rechits
00451   edm::Handle<HBHERecHitCollection> handle;
00452   iEvent.getByLabel(recHitCollName_, handle);
00453   if(!handle.isValid()) {
00454     throw edm::Exception(edm::errors::ProductNotFound)
00455       << " could not find HBHERecHitCollection named " << recHitCollName_ << "\n.";
00456     return;
00457   }
00458 
00459   summary.rechitCount_ = 0;
00460   summary.rechitCount15_ = 0;
00461   summary.rechitEnergy_ = 0;
00462   summary.rechitEnergy15_ = 0;
00463 
00464   // loop over all of the rechit information
00465   for(HBHERecHitCollection::const_iterator it=handle->begin(); it!=handle->end(); ++it) {
00466     const HBHERecHit &rechit=(*it);
00467 
00468     // skip bad rechits (other than those flagged by the isolated noise, triangle, flat, and spike algorithms)
00469     const DetId id = rechit.detid();
00470     uint32_t recHitFlag = rechit.flags();
00471     uint32_t isolbitset = (1 << HcalCaloFlagLabels::HBHEIsolatedNoise);
00472     uint32_t flatbitset = (1 << HcalCaloFlagLabels::HBHEFlatNoise);
00473     uint32_t spikebitset = (1 << HcalCaloFlagLabels::HBHESpikeNoise);
00474     uint32_t trianglebitset = (1 << HcalCaloFlagLabels::HBHETriangleNoise);
00475     uint32_t ts4ts5bitset = (1 << HcalCaloFlagLabels::HBHETS4TS5Noise);
00476     for(unsigned int i=0; i<HcalRecHitFlagsToBeExcluded_.size(); i++) {
00477       uint32_t bitset = (1 << HcalRecHitFlagsToBeExcluded_[i]);
00478       recHitFlag = (recHitFlag & bitset) ? recHitFlag-bitset : recHitFlag;
00479     }
00480     const uint32_t dbStatusFlag = dbHcalChStatus->getValues(id)->getValue();
00481     int severityLevel = hcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
00482     bool isRecovered  = hcalSevLvlComputer->recoveredRecHit(id, recHitFlag);
00483     if(severityLevel!=0 && !isRecovered && severityLevel>static_cast<int>(HcalAcceptSeverityLevel_)) continue;
00484 
00485     // do some rechit counting and energies
00486     summary.rechitCount_ = summary.rechitCount_ + 1;
00487     summary.rechitEnergy_ = summary.rechitEnergy_ + rechit.energy();
00488     if(rechit.energy() > 1.5)
00489     {
00490       summary.rechitCount15_ = summary.rechitCount15_ + 1;
00491       summary.rechitEnergy15_ = summary.rechitEnergy15_ + rechit.energy();
00492     }
00493 
00494     // if it was ID'd as isolated noise, update the summary object
00495     if(rechit.flags() & isolbitset) {
00496       summary.nisolnoise_++;
00497       summary.isolnoisee_ += rechit.energy();
00498       GlobalPoint gp = geo->getPosition(rechit.id());
00499       double et = rechit.energy()*gp.perp()/gp.mag();
00500       summary.isolnoiseet_ += et;
00501     }
00502 
00503     if(rechit.flags() & flatbitset) {
00504       summary.nflatnoise_++;
00505       summary.flatnoisee_ += rechit.energy();
00506       GlobalPoint gp = geo->getPosition(rechit.id());
00507       double et = rechit.energy()*gp.perp()/gp.mag();
00508       summary.flatnoiseet_ += et;
00509     }
00510 
00511     if(rechit.flags() & spikebitset) {
00512       summary.nspikenoise_++;
00513       summary.spikenoisee_ += rechit.energy();
00514       GlobalPoint gp = geo->getPosition(rechit.id());
00515       double et = rechit.energy()*gp.perp()/gp.mag();
00516       summary.spikenoiseet_ += et;
00517     }
00518 
00519     if(rechit.flags() & trianglebitset) {
00520       summary.ntrianglenoise_++;
00521       summary.trianglenoisee_ += rechit.energy();
00522       GlobalPoint gp = geo->getPosition(rechit.id());
00523       double et = rechit.energy()*gp.perp()/gp.mag();
00524       summary.trianglenoiseet_ += et;
00525     }
00526 
00527     if(rechit.flags() & ts4ts5bitset) {
00528       summary.nts4ts5noise_++;
00529       summary.ts4ts5noisee_ += rechit.energy();
00530       GlobalPoint gp = geo->getPosition(rechit.id());
00531       double et = rechit.energy()*gp.perp()/gp.mag();
00532       summary.ts4ts5noiseet_ += et;
00533     }
00534 
00535     // find the hpd that the rechit is in
00536     HcalNoiseHPD& hpd=(*array.findHPD(rechit));
00537 
00538     // create a persistent reference to the rechit
00539     edm::Ref<HBHERecHitCollection> myRef(handle, it-handle->begin());
00540 
00541     // store it in a place so that it remains sorted by energy
00542     hpd.refrechitset_.insert(myRef);
00543 
00544   } // end loop over rechits
00545 
00546   // loop over all HPDs and transfer the information from refrechitset_ to rechits_;
00547   for(HcalNoiseRBXArray::iterator rbxit=array.begin(); rbxit!=array.end(); ++rbxit) {
00548     for(std::vector<HcalNoiseHPD>::iterator hpdit=rbxit->hpds_.begin(); hpdit!=rbxit->hpds_.end(); ++hpdit) {
00549 
00550       // loop over all of the entries in the set and add them to rechits_
00551       for(std::set<edm::Ref<HBHERecHitCollection>, RefHBHERecHitEnergyComparison>::const_iterator
00552             it=hpdit->refrechitset_.begin(); it!=hpdit->refrechitset_.end(); ++it) {
00553         hpdit->rechits_.push_back(*it);
00554       }
00555     }
00556   }
00557   // now the rechits in all the HPDs are sorted by energy!
00558 
00559   return;
00560 }
00561 
00562 // ------------ fill the array with calo tower information
00563 void
00564 HcalNoiseInfoProducer::fillcalotwrs(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseRBXArray& array, HcalNoiseSummary& summary) const
00565 {
00566   // get the calotowers
00567   edm::Handle<CaloTowerCollection> handle;
00568   iEvent.getByLabel(caloTowerCollName_, handle);
00569   if(!handle.isValid()) {
00570     throw edm::Exception(edm::errors::ProductNotFound)
00571       << " could not find CaloTowerCollection named " << caloTowerCollName_ << "\n.";
00572     return;
00573   }
00574 
00575   summary.emenergy_ = summary.hadenergy_ = 0.0;
00576 
00577   // loop over all of the calotower information
00578   for(CaloTowerCollection::const_iterator it = handle->begin(); it!=handle->end(); ++it) {
00579     const CaloTower& twr=(*it);
00580 
00581     // create a persistent reference to the tower
00582     edm::Ref<CaloTowerCollection> myRef(handle, it-handle->begin());
00583 
00584     // get all of the hpd's that are pointed to by the calotower
00585     std::vector<std::vector<HcalNoiseHPD>::iterator> hpditervec;
00586     array.findHPD(twr, hpditervec);
00587 
00588     // loop over the hpd's and add the reference to the RefVectors
00589     for(std::vector<std::vector<HcalNoiseHPD>::iterator>::iterator it=hpditervec.begin();
00590         it!=hpditervec.end(); ++it)
00591       (*it)->calotowers_.push_back(myRef);
00592 
00593     // skip over anything with |ieta|>maxCaloTowerIEta
00594     if(twr.ietaAbs()>maxCaloTowerIEta_) {
00595       summary.emenergy_ += twr.emEnergy();
00596       summary.hadenergy_ += twr.hadEnergy();
00597     }
00598   }
00599 
00600   return;
00601 }
00602 
00603 // ------------ fill the summary with track information
00604 void
00605 HcalNoiseInfoProducer::filltracks(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseSummary& summary) const
00606 {
00607   edm::Handle<reco::TrackCollection> handle;
00608   iEvent.getByLabel(trackCollName_, handle);
00609 
00610   // don't throw exception, just return quietly
00611   if(!handle.isValid()) {
00612     //    throw edm::Exception(edm::errors::ProductNotFound)
00613     //      << " could not find trackCollection named " << trackCollName_ << "\n.";
00614     return;
00615   }
00616 
00617   summary.trackenergy_=0.0;
00618   for(reco::TrackCollection::const_iterator iTrack = handle->begin(); iTrack!=handle->end(); ++iTrack) {
00619     reco::Track trk=*iTrack;
00620     if(trk.pt()<minTrackPt_ || fabs(trk.eta())>maxTrackEta_) continue;
00621 
00622     summary.trackenergy_ += trk.p();
00623   }
00624 
00625   return;
00626 }
00627 
00628 //define this as a plug-in
00629 DEFINE_FWK_MODULE(HcalNoiseInfoProducer);