CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/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 
00058   // if digis are filled, then rechits must also be filled
00059   if(fillDigis_ && !fillRecHits_) {
00060     fillRecHits_=true;
00061     edm::LogWarning("HCalNoiseInfoProducer") << " forcing fillRecHits to be true if fillDigis is true.\n";
00062   }
00063 
00064   // we produce a vector of HcalNoiseRBXs
00065   produces<HcalNoiseRBXCollection>();
00066   // we also produce a noise summary
00067   produces<HcalNoiseSummary>();
00068 }
00069 
00070 
00071 HcalNoiseInfoProducer::~HcalNoiseInfoProducer()
00072 {
00073 }
00074 
00075 
00076 //
00077 // member functions
00078 //
00079 
00080 // ------------ method called to produce the data  ------------
00081 void
00082 HcalNoiseInfoProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00083 {
00084   // this is what we're going to actually write to the EDM
00085   std::auto_ptr<HcalNoiseRBXCollection> result1(new HcalNoiseRBXCollection);
00086   std::auto_ptr<HcalNoiseSummary> result2(new HcalNoiseSummary);
00087   
00088   // define an empty HcalNoiseRBXArray that we're going to fill
00089   HcalNoiseRBXArray rbxarray;
00090   HcalNoiseSummary &summary=*result2;
00091 
00092   // fill them with the various components
00093   // digi assumes that rechit information is available
00094   if(fillRecHits_)    fillrechits(iEvent, iSetup, rbxarray, summary);
00095   if(fillDigis_)      filldigis(iEvent, iSetup, rbxarray);
00096   if(fillCaloTowers_) fillcalotwrs(iEvent, iSetup, rbxarray, summary);
00097   if(fillTracks_)     filltracks(iEvent, iSetup, summary);
00098 
00099   // select those RBXs which are interesting
00100   // also look for the highest energy RBX
00101   HcalNoiseRBXArray::iterator maxit=rbxarray.begin();
00102   double maxenergy=-999;
00103   bool maxwritten=false;
00104   for(HcalNoiseRBXArray::iterator rit = rbxarray.begin(); rit!=rbxarray.end(); ++rit) {
00105     HcalNoiseRBX &rbx=(*rit);
00106     CommonHcalNoiseRBXData data(rbx, minRecHitE_, minLowHitE_, minHighHitE_);
00107 
00108     // find the highest energy rbx
00109     if(data.energy()>maxenergy) {
00110       maxenergy=data.energy();
00111       maxit=rit;
00112       maxwritten=false;
00113     }
00114     
00115     // find out if the rbx is problematic/noisy/etc.
00116     bool writerbx = algo_.isProblematic(data) || !algo_.passLooseNoiseFilter(data) ||
00117       !algo_.passTightNoiseFilter(data) || !algo_.passHighLevelNoiseFilter(data);
00118     
00119     // fill variables in the summary object not filled elsewhere
00120     fillOtherSummaryVariables(summary, data);
00121 
00122     if(writerbx) {
00123       summary.nproblemRBXs_++;
00124       if(summary.nproblemRBXs_<=maxProblemRBXs_) {
00125         result1->push_back(rbx);
00126         if(maxit==rit) maxwritten=true;
00127       }
00128     }
00129   } // end loop over rbxs
00130 
00131   // if we still haven't written the maximum energy rbx, write it now
00132   if(!maxwritten) {
00133     HcalNoiseRBX &rbx=(*maxit);
00134     
00135     // add the RBX to the event
00136     result1->push_back(rbx);
00137   }
00138   
00139   // put the rbxcollection and summary into the EDM
00140   iEvent.put(result1);
00141   iEvent.put(result2);
00142 
00143   return;
00144 }
00145 
00146 // ------------ method called once each job just before starting event loop  ------------
00147 void 
00148 HcalNoiseInfoProducer::beginJob()
00149 {
00150   return;
00151 }
00152 
00153 // ------------ method called once each job just after ending the event loop  ------------
00154 void 
00155 HcalNoiseInfoProducer::endJob()
00156 {
00157   return;
00158 }
00159 
00160 
00161 // ------------ method called once each run just before starting event loop  ------------
00162 // ------------ fills the pedestals
00163 void 
00164 HcalNoiseInfoProducer::beginRun(edm::Run&, const edm::EventSetup&)
00165 {
00166   return;
00167 }
00168 
00169 // ------------ method called once each job just after ending the event loop  ------------
00170 void 
00171 HcalNoiseInfoProducer::endRun(edm::Run&, const edm::EventSetup&)
00172 {
00173   return;
00174 }
00175 
00176 // ------------ here we fill specific variables in the summary object not already accounted for earlier
00177 void 
00178 HcalNoiseInfoProducer::fillOtherSummaryVariables(HcalNoiseSummary& summary, const CommonHcalNoiseRBXData& data) const
00179 {
00180   // charge ratio
00181   if(algo_.passRatioThreshold(data) && data.validRatio()) {
00182     if(data.ratio()<summary.minE2Over10TS()) {
00183       summary.mine2ts_ = data.e2ts();
00184       summary.mine10ts_ = data.e10ts();    }
00185     if(data.ratio()>summary.maxE2Over10TS()) {
00186       summary.maxe2ts_ = data.e2ts();
00187       summary.maxe10ts_ = data.e10ts();
00188     }
00189   }
00190 
00191   // ADC zeros
00192   if(algo_.passZerosThreshold(data)) {
00193     if(data.numZeros()>summary.maxZeros()) {
00194       summary.maxzeros_ = data.numZeros();
00195     }
00196   }
00197 
00198   // hits count
00199   if(data.numHPDHits() > summary.maxHPDHits()) {
00200     summary.maxhpdhits_ = data.numHPDHits();
00201   }
00202   if(data.numRBXHits() > summary.maxRBXHits()) {
00203     summary.maxrbxhits_ = data.numRBXHits();
00204   }
00205   if(data.numHPDNoOtherHits() > summary.maxHPDNoOtherHits()) {
00206     summary.maxhpdhitsnoother_ = data.numHPDNoOtherHits();
00207   }
00208 
00209   // hit timing
00210   if(data.minLowEHitTime()<summary.min10GeVHitTime()) {
00211     summary.min10_ = data.minLowEHitTime();
00212   }
00213   if(data.maxLowEHitTime()>summary.max10GeVHitTime()) {
00214     summary.max10_ = data.maxLowEHitTime();
00215   }
00216   summary.rms10_ += data.lowEHitTimeSqrd();
00217   summary.cnthit10_ += data.numLowEHits();
00218   if(data.minHighEHitTime()<summary.min25GeVHitTime()) {
00219     summary.min25_ = data.minHighEHitTime();
00220   }
00221   if(data.maxHighEHitTime()>summary.max25GeVHitTime()) {
00222     summary.max25_ = data.maxHighEHitTime();
00223   }
00224   summary.rms25_ += data.highEHitTimeSqrd();
00225   summary.cnthit25_ += data.numHighEHits();
00226 
00227   // EMF
00228   if(algo_.passEMFThreshold(data)) {
00229     if(summary.minHPDEMF() > data.HPDEMF()) {
00230       summary.minhpdemf_ = data.HPDEMF();
00231     }
00232     if(summary.minRBXEMF() > data.RBXEMF()) {
00233       summary.minrbxemf_ = data.RBXEMF();
00234     }
00235   }
00236 
00237   // summary flag
00238   if(!algo_.passLooseRatio(data))  summary.filterstatus_ |= 0x1;
00239   if(!algo_.passLooseHits(data))   summary.filterstatus_ |= 0x2;
00240   if(!algo_.passLooseZeros(data))  summary.filterstatus_ |= 0x4;
00241   if(!algo_.passLooseTiming(data)) summary.filterstatus_ |= 0x8;
00242 
00243   if(!algo_.passTightRatio(data))  summary.filterstatus_ |= 0x100;
00244   if(!algo_.passTightHits(data))   summary.filterstatus_ |= 0x200;
00245   if(!algo_.passTightZeros(data))  summary.filterstatus_ |= 0x400;
00246   if(!algo_.passTightTiming(data)) summary.filterstatus_ |= 0x800;
00247 
00248   if(!algo_.passHighLevelNoiseFilter(data)) summary.filterstatus_ |= 0x10000;
00249 
00250   // summary refvectors
00251   JoinCaloTowerRefVectorsWithoutDuplicates join;
00252   if(!algo_.passLooseNoiseFilter(data))
00253     join(summary.loosenoisetwrs_, data.rbxTowers());
00254   if(!algo_.passTightNoiseFilter(data))
00255     join(summary.tightnoisetwrs_, data.rbxTowers());
00256   if(!algo_.passHighLevelNoiseFilter(data))
00257     join(summary.hlnoisetwrs_, data.rbxTowers());
00258 
00259   return;
00260 }
00261 
00262 
00263 // ------------ fill the array with digi information
00264 void
00265 HcalNoiseInfoProducer::filldigis(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseRBXArray& array) const
00266 {
00267 
00268   // get the conditions and channel quality
00269   edm::ESHandle<HcalDbService> conditions;
00270   iSetup.get<HcalDbRecord>().get(conditions);
00271   const HcalQIEShape* shape = conditions->getHcalShape();
00272   edm::ESHandle<HcalChannelQuality> qualhandle;
00273   iSetup.get<HcalChannelQualityRcd>().get(qualhandle);
00274   const HcalChannelQuality* myqual = qualhandle.product();
00275   edm::ESHandle<HcalSeverityLevelComputer> mycomputer;
00276   iSetup.get<HcalSeverityLevelComputerRcd>().get(mycomputer);
00277   const HcalSeverityLevelComputer* mySeverity = mycomputer.product();
00278 
00279   // get the digis
00280   edm::Handle<HBHEDigiCollection> handle;
00281   iEvent.getByLabel(digiCollName_, handle);
00282   if(!handle.isValid()) {
00283     throw edm::Exception(edm::errors::ProductNotFound) << " could not find HBHEDigiCollection named " << digiCollName_ << "\n.";
00284     return;
00285   }
00286 
00287   // loop over all of the digi information
00288   for(HBHEDigiCollection::const_iterator it=handle->begin(); it!=handle->end(); ++it) {
00289     const HBHEDataFrame &digi=(*it);
00290     HcalDetId cell = digi.id();
00291     DetId detcell=(DetId)cell;
00292 
00293     // check on cells to be ignored and dropped
00294     const HcalChannelStatus* mydigistatus=myqual->getValues(detcell.rawId());
00295     if(mySeverity->dropChannel(mydigistatus->getValue())) continue;
00296     if(digi.zsMarkAndPass()) continue;
00297 
00298     // get the calibrations and coder
00299     const HcalCalibrations& calibrations=conditions->getHcalCalibrations(cell);
00300     const HcalQIECoder* channelCoder = conditions->getHcalCoder (cell);
00301     HcalCoderDb coder (*channelCoder, *shape);
00302 
00303     // match the digi to an rbx and hpd
00304     HcalNoiseRBX &rbx=(*array.findRBX(digi));
00305     HcalNoiseHPD &hpd=(*array.findHPD(digi));
00306 
00307     // determine if the digi is one the highest energy hits in the HPD
00308     // this works because the rechits are sorted by energy (see fillrechits() below)
00309     bool isBig=false, isBig5=false, isRBX=false;
00310     int counter=0;
00311     edm::RefVector<HBHERecHitCollection> &rechits=hpd.rechits_;
00312     for(edm::RefVector<HBHERecHitCollection>::const_iterator rit=rechits.begin();
00313         rit!=rechits.end(); ++rit, ++counter) {
00314       if((*rit)->id() == digi.id()) {
00315         if(counter==0) isBig=isBig5=true;  // digi is also the highest energy rechit
00316         if(counter<5) isBig5=true;         // digi is one of 5 highest energy rechits
00317         isRBX=true;
00318       }
00319     }
00320 
00321     // loop over each of the digi's time slices
00322     int totalzeros=0;
00323     CaloSamples tool;
00324     coder.adc2fC(digi,tool);
00325     for(int ts=0; ts<tool.size(); ++ts) {
00326 
00327       // count zero's
00328       if(digi[ts].adc()==0) {
00329         ++hpd.totalZeros_;
00330         ++totalzeros;
00331       }
00332 
00333       // get the fC's
00334       double corrfc = tool[ts]-calibrations.pedestal(digi[ts].capid());
00335 
00336       // fill the relevant digi arrays
00337       if(isBig)  hpd.bigCharge_[ts]+=corrfc;
00338       if(isBig5) hpd.big5Charge_[ts]+=corrfc;
00339       if(isRBX)  rbx.allCharge_[ts]+=corrfc;
00340     }
00341 
00342     // record the maximum number of zero's found
00343     if(totalzeros>hpd.maxZeros_)
00344       hpd.maxZeros_=totalzeros;
00345   }
00346 
00347   return;
00348 }
00349 
00350 // ------------ fill the array with rec hit information
00351 void
00352 HcalNoiseInfoProducer::fillrechits(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseRBXArray& array, HcalNoiseSummary& summary) const
00353 {
00354   // get the HCAL channel status map
00355   edm::ESHandle<HcalChannelQuality> hcalChStatus;    
00356   iSetup.get<HcalChannelQualityRcd>().get( hcalChStatus );
00357   const HcalChannelQuality* dbHcalChStatus = hcalChStatus.product();
00358 
00359   // get the severity level computer
00360   edm::ESHandle<HcalSeverityLevelComputer> hcalSevLvlComputerHndl;
00361   iSetup.get<HcalSeverityLevelComputerRcd>().get(hcalSevLvlComputerHndl);
00362   const HcalSeverityLevelComputer* hcalSevLvlComputer = hcalSevLvlComputerHndl.product();
00363 
00364   // get the calo geometry
00365   edm::ESHandle<CaloGeometry> pG;
00366   iSetup.get<CaloGeometryRecord>().get(pG);
00367   const CaloGeometry* geo = pG.product();
00368 
00369   // get the rechits
00370   edm::Handle<HBHERecHitCollection> handle;
00371   iEvent.getByLabel(recHitCollName_, handle);
00372   if(!handle.isValid()) {
00373     throw edm::Exception(edm::errors::ProductNotFound)
00374       << " could not find HBHERecHitCollection named " << recHitCollName_ << "\n.";
00375     return;
00376   }
00377 
00378   // loop over all of the rechit information
00379   for(HBHERecHitCollection::const_iterator it=handle->begin(); it!=handle->end(); ++it) {
00380     const HBHERecHit &rechit=(*it);
00381 
00382     // skip bad rechits (other than those flagged by the isolated noise algorithm)
00383     const DetId id = rechit.detid();
00384     uint32_t recHitFlag = rechit.flags();    
00385     uint32_t noisebitset = (1 << HcalCaloFlagLabels::HBHEIsolatedNoise);
00386     recHitFlag = (recHitFlag & noisebitset) ? recHitFlag-noisebitset : recHitFlag;
00387     const uint32_t dbStatusFlag = dbHcalChStatus->getValues(id)->getValue();
00388     int severityLevel = hcalSevLvlComputer->getSeverityLevel(id, recHitFlag, dbStatusFlag);
00389     bool isRecovered  = hcalSevLvlComputer->recoveredRecHit(id, recHitFlag);
00390     if(severityLevel!=0 && !isRecovered && severityLevel>static_cast<int>(HcalAcceptSeverityLevel_)) continue;
00391 
00392     // if it was ID'd as isolated noise, update the summary object
00393     if(rechit.flags() & noisebitset) {
00394       ++summary.nisolnoise_;
00395       summary.isolnoisee_ += rechit.energy();
00396       GlobalPoint gp = geo->getPosition(rechit.id());
00397       double et = rechit.energy()*gp.perp()/gp.mag();
00398       summary.isolnoiseet_ += et;
00399     }
00400 
00401     // find the hpd that the rechit is in
00402     HcalNoiseHPD& hpd=(*array.findHPD(rechit));
00403 
00404     // create a persistent reference to the rechit
00405     edm::Ref<HBHERecHitCollection> myRef(handle, it-handle->begin());
00406     
00407     // store it in a place so that it remains sorted by energy
00408     hpd.refrechitset_.insert(myRef);
00409 
00410   } // end loop over rechits
00411 
00412   // loop over all HPDs and transfer the information from refrechitset_ to rechits_;
00413   for(HcalNoiseRBXArray::iterator rbxit=array.begin(); rbxit!=array.end(); ++rbxit) {
00414     for(std::vector<HcalNoiseHPD>::iterator hpdit=rbxit->hpds_.begin(); hpdit!=rbxit->hpds_.end(); ++hpdit) {
00415       
00416       // loop over all of the entries in the set and add them to rechits_
00417       for(std::set<edm::Ref<HBHERecHitCollection>, RefHBHERecHitEnergyComparison>::const_iterator
00418             it=hpdit->refrechitset_.begin(); it!=hpdit->refrechitset_.end(); ++it) {
00419         hpdit->rechits_.push_back(*it);
00420       }
00421     }
00422   }
00423   // now the rechits in all the HPDs are sorted by energy!
00424 
00425   return;
00426 }
00427 
00428 // ------------ fill the array with calo tower information
00429 void
00430 HcalNoiseInfoProducer::fillcalotwrs(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseRBXArray& array, HcalNoiseSummary& summary) const
00431 {
00432   // get the calotowers
00433   edm::Handle<CaloTowerCollection> handle;
00434   iEvent.getByLabel(caloTowerCollName_, handle);
00435   if(!handle.isValid()) {
00436     throw edm::Exception(edm::errors::ProductNotFound)
00437       << " could not find CaloTowerCollection named " << caloTowerCollName_ << "\n.";
00438     return;
00439   }
00440 
00441   summary.emenergy_ = summary.hadenergy_ = 0.0;
00442 
00443   // loop over all of the calotower information
00444   for(CaloTowerCollection::const_iterator it = handle->begin(); it!=handle->end(); ++it) {
00445     const CaloTower& twr=(*it);
00446 
00447     // create a persistent reference to the tower
00448     edm::Ref<CaloTowerCollection> myRef(handle, it-handle->begin());
00449 
00450     // get all of the hpd's that are pointed to by the calotower
00451     std::vector<std::vector<HcalNoiseHPD>::iterator> hpditervec;
00452     array.findHPD(twr, hpditervec);
00453 
00454     // loop over the hpd's and add the reference to the RefVectors
00455     for(std::vector<std::vector<HcalNoiseHPD>::iterator>::iterator it=hpditervec.begin();
00456         it!=hpditervec.end(); ++it)
00457       (*it)->calotowers_.push_back(myRef);
00458 
00459     // skip over anything with |ieta|>maxCaloTowerIEta
00460     if(twr.ietaAbs()>maxCaloTowerIEta_) {
00461       summary.emenergy_ += twr.emEnergy();
00462       summary.hadenergy_ += twr.hadEnergy();
00463     }
00464   }
00465 
00466   return;
00467 }
00468 
00469 // ------------ fill the summary with track information
00470 void
00471 HcalNoiseInfoProducer::filltracks(edm::Event& iEvent, const edm::EventSetup& iSetup, HcalNoiseSummary& summary) const
00472 {
00473   edm::Handle<reco::TrackCollection> handle;
00474   iEvent.getByLabel(trackCollName_, handle);
00475 
00476   // don't throw exception, just return quietly
00477   if(!handle.isValid()) {
00478     //    throw edm::Exception(edm::errors::ProductNotFound)
00479     //      << " could not find trackCollection named " << trackCollName_ << "\n.";
00480     return;
00481   }
00482 
00483   summary.trackenergy_=0.0;
00484   for(reco::TrackCollection::const_iterator iTrack = handle->begin(); iTrack!=handle->end(); ++iTrack) {
00485     reco::Track trk=*iTrack;
00486     if(trk.pt()<minTrackPt_ || fabs(trk.eta())>maxTrackEta_) continue;
00487 
00488     summary.trackenergy_ += trk.p();
00489   }
00490 
00491   return;
00492 }
00493 
00494 //define this as a plug-in
00495 DEFINE_FWK_MODULE(HcalNoiseInfoProducer);