![]() |
![]() |
#include <HLTEgammaDoubleEtDeltaPhiFilter.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HLTEgammaDoubleEtDeltaPhiFilter (const edm::ParameterSet &) | |
~HLTEgammaDoubleEtDeltaPhiFilter () | |
Private Attributes | |
double | etcut_ |
edm::InputTag | inputTag_ |
edm::InputTag | L1IsoCollTag_ |
edm::InputTag | L1NonIsoCollTag_ |
double | minDeltaPhi_ |
bool | relaxed_ |
bool | store_ |
This filter will require only two HLT photons and DeltaPhi between the two photons larger than 2.5
Definition at line 16 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
HLTEgammaDoubleEtDeltaPhiFilter::HLTEgammaDoubleEtDeltaPhiFilter | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 22 of file HLTEgammaDoubleEtDeltaPhiFilter.cc.
References etcut_, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), inputTag_, L1IsoCollTag_, L1NonIsoCollTag_, minDeltaPhi_, relaxed_, and store_.
{ inputTag_ = iConfig.getParameter< edm::InputTag > ("inputTag"); etcut_ = iConfig.getParameter<double> ("etcut"); minDeltaPhi_ = iConfig.getParameter<double> ("minDeltaPhi"); // ncandcut_ = iConfig.getParameter<int> ("ncandcut",2); store_ = iConfig.getUntrackedParameter<bool> ("saveTag",false) ; relaxed_ = iConfig.getUntrackedParameter<bool> ("relaxed",true) ; L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand"); L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand"); //register your products produces<trigger::TriggerFilterObjectWithRefs>(); }
HLTEgammaDoubleEtDeltaPhiFilter::~HLTEgammaDoubleEtDeltaPhiFilter | ( | ) |
Definition at line 37 of file HLTEgammaDoubleEtDeltaPhiFilter.cc.
{}
bool HLTEgammaDoubleEtDeltaPhiFilter::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements HLTFilter.
Definition at line 41 of file HLTEgammaDoubleEtDeltaPhiFilter.cc.
References accept(), Geom::deltaPhi(), etcut_, edm::Event::getByLabel(), i, inputTag_, L1IsoCollTag_, L1NonIsoCollTag_, M_PI, minDeltaPhi_, module(), n, path(), edm::Event::put(), relaxed_, store_, and trigger::TriggerCluster.
{ using namespace trigger; // The filter object std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module())); if( store_ ){filterproduct->addCollectionTag(L1IsoCollTag_);} if( store_ && relaxed_){filterproduct->addCollectionTag(L1NonIsoCollTag_);} // get hold of filtered candidates edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput; iEvent.getByLabel (inputTag_,PrevFilterOutput); std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands; PrevFilterOutput->getObjects(TriggerCluster, recoecalcands); // Refs to the two Candidate objects used to calculate deltaPhi edm::Ref<reco::RecoEcalCandidateCollection> ref1, ref2; // look at all candidates, check cuts int n(0); for(unsigned int i=0; i<recoecalcands.size(); i++) { const edm::Ref<reco::RecoEcalCandidateCollection> & ref = recoecalcands[i]; if( ref->et() >= etcut_) { ++n; if(n==1) ref1 = ref; if(n==2) ref2 = ref; } } // if there are only two Candidate objects, calculate deltaPhi double deltaPhi(0.0); if(n==2) { deltaPhi = fabs(ref1->phi()-ref2->phi()); if(deltaPhi>M_PI) deltaPhi = 2*M_PI - deltaPhi; filterproduct->addObject(TriggerCluster, ref1); filterproduct->addObject(TriggerCluster, ref2); } // filter decision bool accept(n==2 && deltaPhi>minDeltaPhi_); // put filter object into the Event iEvent.put(filterproduct); return accept; }
double HLTEgammaDoubleEtDeltaPhiFilter::etcut_ [private] |
Definition at line 25 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
Referenced by filter(), and HLTEgammaDoubleEtDeltaPhiFilter().
Definition at line 24 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
Referenced by filter(), and HLTEgammaDoubleEtDeltaPhiFilter().
Definition at line 30 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
Referenced by filter(), and HLTEgammaDoubleEtDeltaPhiFilter().
Definition at line 31 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
Referenced by filter(), and HLTEgammaDoubleEtDeltaPhiFilter().
double HLTEgammaDoubleEtDeltaPhiFilter::minDeltaPhi_ [private] |
Definition at line 26 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
Referenced by filter(), and HLTEgammaDoubleEtDeltaPhiFilter().
bool HLTEgammaDoubleEtDeltaPhiFilter::relaxed_ [private] |
Definition at line 29 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
Referenced by filter(), and HLTEgammaDoubleEtDeltaPhiFilter().
bool HLTEgammaDoubleEtDeltaPhiFilter::store_ [private] |
Definition at line 28 of file HLTEgammaDoubleEtDeltaPhiFilter.h.
Referenced by filter(), and HLTEgammaDoubleEtDeltaPhiFilter().