![]() |
![]() |
00001 00009 #include "HLTrigger/Egamma/interface/HLTEgammaDoubleEtDeltaPhiFilter.h" 00010 00011 #include "DataFormats/Common/interface/Handle.h" 00012 00013 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h" 00014 00015 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00016 00017 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h" 00018 00019 // 00020 // constructors and destructor 00021 // 00022 HLTEgammaDoubleEtDeltaPhiFilter::HLTEgammaDoubleEtDeltaPhiFilter(const edm::ParameterSet& iConfig) 00023 { 00024 inputTag_ = iConfig.getParameter< edm::InputTag > ("inputTag"); 00025 etcut_ = iConfig.getParameter<double> ("etcut"); 00026 minDeltaPhi_ = iConfig.getParameter<double> ("minDeltaPhi"); 00027 // ncandcut_ = iConfig.getParameter<int> ("ncandcut",2); 00028 store_ = iConfig.getUntrackedParameter<bool> ("saveTag",false) ; 00029 relaxed_ = iConfig.getUntrackedParameter<bool> ("relaxed",true) ; 00030 L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand"); 00031 L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand"); 00032 00033 //register your products 00034 produces<trigger::TriggerFilterObjectWithRefs>(); 00035 } 00036 00037 HLTEgammaDoubleEtDeltaPhiFilter::~HLTEgammaDoubleEtDeltaPhiFilter(){} 00038 00039 // ------------ method called to produce the data ------------ 00040 bool 00041 HLTEgammaDoubleEtDeltaPhiFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) 00042 { 00043 using namespace trigger; 00044 // The filter object 00045 std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module())); 00046 if( store_ ){filterproduct->addCollectionTag(L1IsoCollTag_);} 00047 if( store_ && relaxed_){filterproduct->addCollectionTag(L1NonIsoCollTag_);} 00048 00049 // get hold of filtered candidates 00050 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput; 00051 iEvent.getByLabel (inputTag_,PrevFilterOutput); 00052 00053 std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands; 00054 PrevFilterOutput->getObjects(TriggerCluster, recoecalcands); 00055 00056 // Refs to the two Candidate objects used to calculate deltaPhi 00057 edm::Ref<reco::RecoEcalCandidateCollection> ref1, ref2; 00058 00059 // look at all candidates, check cuts 00060 int n(0); 00061 for(unsigned int i=0; i<recoecalcands.size(); i++) { 00062 const edm::Ref<reco::RecoEcalCandidateCollection> & ref = recoecalcands[i]; 00063 if( ref->et() >= etcut_) { 00064 ++n; 00065 if(n==1) ref1 = ref; 00066 if(n==2) ref2 = ref; 00067 } 00068 } 00069 00070 // if there are only two Candidate objects, calculate deltaPhi 00071 double deltaPhi(0.0); 00072 if(n==2) { 00073 deltaPhi = fabs(ref1->phi()-ref2->phi()); 00074 if(deltaPhi>M_PI) deltaPhi = 2*M_PI - deltaPhi; 00075 00076 filterproduct->addObject(TriggerCluster, ref1); 00077 filterproduct->addObject(TriggerCluster, ref2); 00078 } 00079 00080 // filter decision 00081 bool accept(n==2 && deltaPhi>minDeltaPhi_); 00082 00083 // put filter object into the Event 00084 iEvent.put(filterproduct); 00085 00086 return accept; 00087 } 00088 00089 // define as a framework module 00090 // #include "FWCore/Framework/interface/MakerMacros.h" 00091 // DEFINE_FWK_MODULE(HLTEgammaDoubleEtDeltaPhiFilter); 00092