#include <HBHEIsolatedNoiseAlgos.h>
Definition at line 70 of file HBHEIsolatedNoiseAlgos.h.
ObjectValidator::ObjectValidator | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 42 of file HBHEIsolatedNoiseAlgos.cc.
References EBThreshold_, EcalAcceptSeverityLevel_, EEThreshold_, edm::ParameterSet::getParameter(), HBThreshold_, HcalAcceptSeverityLevel_, HEDThreshold_, HESThreshold_, MinValidTrackNHits_, MinValidTrackPt_, MinValidTrackPtBarrel_, theEBRecHitCollection_, theEcalChStatus_, theEcalSevLvlAlgo_, theEERecHitCollection_, theHcalChStatus_, theHcalSevLvlComputer_, UseEcalRecoveredHits_, and UseHcalRecoveredHits_.
{ HBThreshold_ = iConfig.getParameter<double>("HBThreshold"); HESThreshold_ = iConfig.getParameter<double>("HESThreshold"); HEDThreshold_ = iConfig.getParameter<double>("HEDThreshold"); EBThreshold_ = iConfig.getParameter<double>("EBThreshold"); EEThreshold_ = iConfig.getParameter<double>("EEThreshold"); HcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("HcalAcceptSeverityLevel"); EcalAcceptSeverityLevel_ = iConfig.getParameter<uint32_t>("EcalAcceptSeverityLevel"); UseHcalRecoveredHits_ = iConfig.getParameter<bool>("UseHcalRecoveredHits"); UseEcalRecoveredHits_ = iConfig.getParameter<bool>("UseEcalRecoveredHits"); MinValidTrackPt_ = iConfig.getParameter<double>("MinValidTrackPt"); MinValidTrackPtBarrel_ = iConfig.getParameter<double>("MinValidTrackPtBarrel"); MinValidTrackNHits_ = iConfig.getParameter<int>("MinValidTrackNHits"); theHcalChStatus_=0; theEcalChStatus_=0; theHcalSevLvlComputer_=0; theEcalSevLvlAlgo_=0; theEBRecHitCollection_=0; theEERecHitCollection_=0; return; }
ObjectValidator::ObjectValidator | ( | double | HBThreshold, |
double | HESThreshold, | ||
double | HEDThreshold, | ||
double | EBThreshold, | ||
double | EEThreshold, | ||
uint32_t | HcalAcceptSeverityLevel, | ||
uint32_t | EcalAcceptSeverityLevel, | ||
bool | UseHcalRecoveredHits, | ||
bool | UseEcalRecoveredHits, | ||
double | MinValidTrackPt, | ||
double | MinValidTrackPtBarrel, | ||
int | MinValidTrackNHits | ||
) | [inline] |
Definition at line 74 of file HBHEIsolatedNoiseAlgos.h.
: HBThreshold_(HBThreshold), HESThreshold_(HESThreshold), HEDThreshold_(HEDThreshold), EBThreshold_(EBThreshold), EEThreshold_(EEThreshold), HcalAcceptSeverityLevel_(HcalAcceptSeverityLevel), EcalAcceptSeverityLevel_(EcalAcceptSeverityLevel), UseHcalRecoveredHits_(UseHcalRecoveredHits), UseEcalRecoveredHits_(UseEcalRecoveredHits), MinValidTrackPt_(MinValidTrackPt), MinValidTrackPtBarrel_(MinValidTrackPtBarrel), MinValidTrackNHits_(MinValidTrackNHits), theHcalChStatus_(0), theEcalChStatus_(0), theHcalSevLvlComputer_(0), theEcalSevLvlAlgo_(0), theEBRecHitCollection_(0), theEERecHitCollection_(0) {}
ObjectValidator::~ObjectValidator | ( | ) | [virtual] |
Definition at line 69 of file HBHEIsolatedNoiseAlgos.cc.
{ }
void ObjectValidator::setEBRecHitCollection | ( | const EcalRecHitCollection * | q | ) | [inline] |
Definition at line 101 of file HBHEIsolatedNoiseAlgos.h.
References lumiQueryAPI::q, and theEBRecHitCollection_.
Referenced by HBHEIsolatedNoiseReflagger::produce().
{ theEBRecHitCollection_=q; }
void ObjectValidator::setEcalChannelStatus | ( | const EcalChannelStatus * | q | ) | [inline] |
Definition at line 98 of file HBHEIsolatedNoiseAlgos.h.
References lumiQueryAPI::q, and theEcalChStatus_.
Referenced by HBHEIsolatedNoiseReflagger::produce().
{ theEcalChStatus_=q; }
void ObjectValidator::setEcalSeverityLevelAlgo | ( | const EcalSeverityLevelAlgo * | q | ) | [inline] |
Definition at line 100 of file HBHEIsolatedNoiseAlgos.h.
References lumiQueryAPI::q, and theEcalSevLvlAlgo_.
Referenced by HBHEIsolatedNoiseReflagger::produce().
{ theEcalSevLvlAlgo_=q; }
void ObjectValidator::setEERecHitCollection | ( | const EcalRecHitCollection * | q | ) | [inline] |
Definition at line 102 of file HBHEIsolatedNoiseAlgos.h.
References lumiQueryAPI::q, and theEERecHitCollection_.
Referenced by HBHEIsolatedNoiseReflagger::produce().
{ theEERecHitCollection_=q; }
void ObjectValidator::setHcalChannelQuality | ( | const HcalChannelQuality * | q | ) | [inline] |
Definition at line 97 of file HBHEIsolatedNoiseAlgos.h.
References lumiQueryAPI::q, and theHcalChStatus_.
Referenced by HBHEIsolatedNoiseReflagger::produce().
{ theHcalChStatus_=q; }
void ObjectValidator::setHcalSeverityLevelComputer | ( | const HcalSeverityLevelComputer * | q | ) | [inline] |
Definition at line 99 of file HBHEIsolatedNoiseAlgos.h.
References lumiQueryAPI::q, and theHcalSevLvlComputer_.
Referenced by HBHEIsolatedNoiseReflagger::produce().
{ theHcalSevLvlComputer_=q; }
bool ObjectValidator::validHit | ( | const EcalRecHit & | hit | ) | const [virtual] |
Implements ObjectValidatorAbs.
Definition at line 95 of file HBHEIsolatedNoiseAlgos.cc.
References CaloRecHit::detid(), EBThreshold_, EcalAcceptSeverityLevel_, EcalBarrel, EcalEndcap, EEThreshold_, CaloRecHit::energy(), EcalSeverityLevelAlgo::kGood, EcalSeverityLevelAlgo::kRecovered, EcalSeverityLevelAlgo::kSwissCross, EcalSeverityLevelAlgo::severityLevel(), theEBRecHitCollection_, theEcalChStatus_, theEcalSevLvlAlgo_, theEERecHitCollection_, and UseEcalRecoveredHits_.
{ assert(theEcalSevLvlAlgo_!=0 && theEcalChStatus_!=0); // require the hit to pass a certain energy threshold const DetId id = hit.detid(); if(id.subdetId()==EcalBarrel && hit.energy()<EBThreshold_) return false; else if(id.subdetId()==EcalEndcap && hit.energy()<EEThreshold_) return false; // determine if the hit is good, bad, or recovered int severityLevel = 999; if (id.subdetId() == EcalBarrel && theEBRecHitCollection_!=0) severityLevel = theEcalSevLvlAlgo_->severityLevel(id, *theEBRecHitCollection_, *theEcalChStatus_, 5., EcalSeverityLevelAlgo::kSwissCross, 0.95, 2., 15., 0.999); else if(id.subdetId() == EcalEndcap && theEERecHitCollection_!=0) severityLevel = theEcalSevLvlAlgo_->severityLevel(id, *theEERecHitCollection_, *theEcalChStatus_, 5., EcalSeverityLevelAlgo::kSwissCross, 0.95, 2., 15., 0.999); else return false; if(severityLevel == EcalSeverityLevelAlgo::kGood) return true; if(severityLevel == EcalSeverityLevelAlgo::kRecovered) return UseEcalRecoveredHits_; if(severityLevel > static_cast<int>(EcalAcceptSeverityLevel_)) return false; else return true; }
bool ObjectValidator::validHit | ( | const HBHERecHit & | hit | ) | const [virtual] |
Implements ObjectValidatorAbs.
Definition at line 73 of file HBHEIsolatedNoiseAlgos.cc.
References CaloRecHit::detid(), CaloRecHit::energy(), CaloRecHit::flags(), HcalSeverityLevelComputer::getSeverityLevel(), HcalChannelStatus::getValue(), HcalCondObjectContainer< Item >::getValues(), HBThreshold_, HcalAcceptSeverityLevel_, HcalBarrel, HcalEndcap, HEDThreshold_, HESThreshold_, HBHERecHit::id(), HcalDetId::ietaAbs(), HcalSeverityLevelComputer::recoveredRecHit(), HcalDetId::subdet(), theHcalChStatus_, theHcalSevLvlComputer_, and UseHcalRecoveredHits_.
{ assert(theHcalSevLvlComputer_!=0 && theHcalChStatus_!=0); // require the hit to pass a certain energy threshold if(hit.id().subdet()==HcalBarrel && hit.energy()<HBThreshold_) return false; else if(hit.id().subdet()==HcalEndcap && hit.id().ietaAbs()<=20 && hit.energy()<HESThreshold_) return false; else if(hit.id().subdet()==HcalEndcap && hit.id().ietaAbs()>20 && hit.energy()<HEDThreshold_) return false; // determine if the hit is good, bad, or recovered const DetId id = hit.detid(); const uint32_t recHitFlag = hit.flags(); const uint32_t dbStatusFlag = theHcalChStatus_->getValues(id)->getValue(); int severityLevel = theHcalSevLvlComputer_->getSeverityLevel(id, recHitFlag, dbStatusFlag); bool isRecovered = theHcalSevLvlComputer_->recoveredRecHit(id, recHitFlag); if(severityLevel==0) return true; if(isRecovered) return UseHcalRecoveredHits_; if(severityLevel>static_cast<int>(HcalAcceptSeverityLevel_)) return false; else return true; }
bool ObjectValidator::validTrack | ( | const reco::Track & | trk | ) | const [virtual] |
Implements ObjectValidatorAbs.
Definition at line 116 of file HBHEIsolatedNoiseAlgos.cc.
References MinValidTrackNHits_, MinValidTrackPt_, MinValidTrackPtBarrel_, reco::TrackBase::momentum(), reco::TrackBase::numberOfValidHits(), and reco::TrackBase::pt().
{ if(trk.pt()<MinValidTrackPt_) return false; if(trk.pt()<MinValidTrackPtBarrel_ && std::fabs(trk.momentum().eta())<1.479) return false; if(trk.numberOfValidHits()<MinValidTrackNHits_) return false; return true; }
double ObjectValidator::EBThreshold_ [private] |
Definition at line 114 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
uint32_t ObjectValidator::EcalAcceptSeverityLevel_ [private] |
Definition at line 118 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
double ObjectValidator::EEThreshold_ [private] |
Definition at line 115 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
double ObjectValidator::HBThreshold_ [private] |
Definition at line 111 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
uint32_t ObjectValidator::HcalAcceptSeverityLevel_ [private] |
Definition at line 117 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
double ObjectValidator::HEDThreshold_ [private] |
Definition at line 113 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
double ObjectValidator::HESThreshold_ [private] |
Definition at line 112 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
int ObjectValidator::MinValidTrackNHits_ [private] |
Definition at line 124 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validTrack().
double ObjectValidator::MinValidTrackPt_ [private] |
Definition at line 122 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validTrack().
double ObjectValidator::MinValidTrackPtBarrel_ [private] |
Definition at line 123 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validTrack().
const EcalRecHitCollection* ObjectValidator::theEBRecHitCollection_ [private] |
Definition at line 135 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), setEBRecHitCollection(), and validHit().
const EcalChannelStatus* ObjectValidator::theEcalChStatus_ [private] |
Definition at line 128 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), setEcalChannelStatus(), and validHit().
const EcalSeverityLevelAlgo* ObjectValidator::theEcalSevLvlAlgo_ [private] |
Definition at line 132 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), setEcalSeverityLevelAlgo(), and validHit().
const EcalRecHitCollection* ObjectValidator::theEERecHitCollection_ [private] |
Definition at line 136 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), setEERecHitCollection(), and validHit().
const HcalChannelQuality* ObjectValidator::theHcalChStatus_ [private] |
Definition at line 127 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), setHcalChannelQuality(), and validHit().
const HcalSeverityLevelComputer* ObjectValidator::theHcalSevLvlComputer_ [private] |
Definition at line 131 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), setHcalSeverityLevelComputer(), and validHit().
bool ObjectValidator::UseEcalRecoveredHits_ [private] |
Definition at line 120 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().
bool ObjectValidator::UseHcalRecoveredHits_ [private] |
Definition at line 119 of file HBHEIsolatedNoiseAlgos.h.
Referenced by ObjectValidator(), and validHit().