![]() |
![]() |
#include <Work/HLTHcalSimpleRecHitFilter/src/HLTHcalSimpleRecHitFilter.cc>
Public Member Functions | |
HLTHcalSimpleRecHitFilter (const edm::ParameterSet &) | |
~HLTHcalSimpleRecHitFilter () | |
Private Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
Private Attributes | |
bool | doCoincidence_ |
edm::InputTag | HcalRecHitCollection_ |
std::vector< int > | maskedList_ |
int | minNHitsNeg_ |
int | minNHitsPos_ |
double | threshold_ |
Description: <one line="" class="" summary>="">
Implementation: <Notes on="" implementation>="">
Definition at line 40 of file HLTHcalSimpleRecHitFilter.cc.
HLTHcalSimpleRecHitFilter::HLTHcalSimpleRecHitFilter | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 61 of file HLTHcalSimpleRecHitFilter.cc.
References doCoincidence_, edm::ParameterSet::getParameter(), HcalRecHitCollection_, maskedList_, minNHitsNeg_, minNHitsPos_, and threshold_.
{ //now do what ever initialization is needed threshold_ = iConfig.getParameter<double>("threshold"); minNHitsNeg_ = iConfig.getParameter<int>("minNHitsNeg"); minNHitsPos_ = iConfig.getParameter<int>("minNHitsPos"); doCoincidence_ = iConfig.getParameter<bool>("doCoincidence"); maskedList_ = iConfig.getParameter<std::vector<int> >("maskedChannels"); //this is using the hashed index HcalRecHitCollection_ = iConfig.getParameter<edm::InputTag>("HFRecHitCollection"); }
HLTHcalSimpleRecHitFilter::~HLTHcalSimpleRecHitFilter | ( | ) |
Definition at line 74 of file HLTHcalSimpleRecHitFilter.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
bool HLTHcalSimpleRecHitFilter::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements HLTFilter.
Definition at line 89 of file HLTHcalSimpleRecHitFilter.cc.
References accept(), doCoincidence_, CaloRecHit::energy(), exception, funct::false, spr::find(), edm::Event::getByLabel(), HcalRecHitCollection_, HFRecHit::id(), maskedList_, minNHitsNeg_, minNHitsPos_, query::result, threshold_, and HcalDetId::zside().
{ // using namespace edm; // getting very basic uncalRH edm::Handle<HFRecHitCollection> crudeHits; try { iEvent.getByLabel(HcalRecHitCollection_, crudeHits); } catch ( std::exception& ex) { edm::LogWarning("HLTHcalSimpleRecHitFilter") << HcalRecHitCollection_ << " not available"; } bool accept = false ; int nHitsNeg=0, nHitsPos=0; for ( HFRecHitCollection::const_iterator hitItr = crudeHits->begin(); hitItr != crudeHits->end(); ++hitItr ) { HFRecHit hit = (*hitItr); // masking noisy channels std::vector<int>::iterator result; result = find( maskedList_.begin(), maskedList_.end(), HcalDetId(hit.id()).hashed_index() ); if (result != maskedList_.end()) continue; // only count tower above threshold if ( hit.energy() < threshold_ ) continue; // count if (hit.id().zside()<0) ++nHitsNeg; else ++nHitsPos; } // Logic if (!doCoincidence_) accept = (nHitsNeg>=minNHitsNeg_) || (nHitsPos>=minNHitsPos_); else accept = (nHitsNeg>=minNHitsNeg_) && (nHitsPos>=minNHitsPos_); // edm::LogInfo("HcalFilter") << "at evet: " << iEvent.id().event() // << " and run: " << iEvent.id().run() // << " Total HF hits: " << crudeHits->size() << " Above Threshold - nNeg: " << nHitsNeg << " nPos " << nHitsPos // << " doCoinc: " << doCoincidence_ << " accept: " << accept << std::endl; // result return accept; }
bool HLTHcalSimpleRecHitFilter::doCoincidence_ [private] |
Definition at line 53 of file HLTHcalSimpleRecHitFilter.cc.
Referenced by filter(), and HLTHcalSimpleRecHitFilter().
Definition at line 49 of file HLTHcalSimpleRecHitFilter.cc.
Referenced by filter(), and HLTHcalSimpleRecHitFilter().
std::vector<int> HLTHcalSimpleRecHitFilter::maskedList_ [private] |
Definition at line 54 of file HLTHcalSimpleRecHitFilter.cc.
Referenced by filter(), and HLTHcalSimpleRecHitFilter().
int HLTHcalSimpleRecHitFilter::minNHitsNeg_ [private] |
Definition at line 51 of file HLTHcalSimpleRecHitFilter.cc.
Referenced by filter(), and HLTHcalSimpleRecHitFilter().
int HLTHcalSimpleRecHitFilter::minNHitsPos_ [private] |
Definition at line 52 of file HLTHcalSimpleRecHitFilter.cc.
Referenced by filter(), and HLTHcalSimpleRecHitFilter().
double HLTHcalSimpleRecHitFilter::threshold_ [private] |
Definition at line 50 of file HLTHcalSimpleRecHitFilter.cc.
Referenced by filter(), and HLTHcalSimpleRecHitFilter().