![]() |
![]() |
00001 00009 #include "HLTrigger/Egamma/interface/HLTEgammaGenericQuadraticFilter.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 #include "DataFormats/EgammaReco/interface/SuperCluster.h" 00019 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h" 00020 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidateIsolation.h" 00021 #include "DataFormats/Common/interface/AssociationMap.h" 00022 #include "DataFormats/RecoCandidate/interface/RecoCandidate.h" 00023 00024 // 00025 // constructors and destructor 00026 // 00027 HLTEgammaGenericQuadraticFilter::HLTEgammaGenericQuadraticFilter(const edm::ParameterSet& iConfig){ 00028 candTag_ = iConfig.getParameter< edm::InputTag > ("candTag"); 00029 isoTag_ = iConfig.getParameter< edm::InputTag > ("isoTag"); 00030 nonIsoTag_ = iConfig.getParameter< edm::InputTag > ("nonIsoTag"); 00031 00032 lessThan_ = iConfig.getParameter<bool> ("lessThan"); 00033 useEt_ = iConfig.getParameter<bool> ("useEt"); 00034 thrRegularEB_ = iConfig.getParameter<double> ("thrRegularEB"); 00035 thrRegularEE_ = iConfig.getParameter<double> ("thrRegularEE"); 00036 thrOverEEB_ = iConfig.getParameter<double> ("thrOverEEB"); 00037 thrOverEEE_ = iConfig.getParameter<double> ("thrOverEEE"); 00038 thrOverE2EB_ = iConfig.getParameter<double> ("thrOverE2EB"); 00039 thrOverE2EE_ = iConfig.getParameter<double> ("thrOverE2EE"); 00040 00041 ncandcut_ = iConfig.getParameter<int> ("ncandcut"); 00042 doIsolated_ = iConfig.getParameter<bool> ("doIsolated"); 00043 00044 store_ = iConfig.getUntrackedParameter<bool> ("SaveTag",false) ; 00045 L1IsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1IsoCand"); 00046 L1NonIsoCollTag_= iConfig.getParameter< edm::InputTag > ("L1NonIsoCand"); 00047 00048 //register your products 00049 produces<trigger::TriggerFilterObjectWithRefs>(); 00050 } 00051 00052 HLTEgammaGenericQuadraticFilter::~HLTEgammaGenericQuadraticFilter(){} 00053 00054 00055 // ------------ method called to produce the data ------------ 00056 bool 00057 HLTEgammaGenericQuadraticFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) 00058 { 00059 using namespace trigger; 00060 std::auto_ptr<trigger::TriggerFilterObjectWithRefs> filterproduct (new trigger::TriggerFilterObjectWithRefs(path(),module())); 00061 if( store_ ){filterproduct->addCollectionTag(L1IsoCollTag_);} 00062 if( store_ && !doIsolated_){filterproduct->addCollectionTag(L1NonIsoCollTag_);} 00063 00064 // Ref to Candidate object to be recorded in filter object 00065 edm::Ref<reco::RecoEcalCandidateCollection> ref; 00066 00067 // Set output format 00068 int trigger_type = trigger::TriggerCluster; 00069 if ( store_ ) trigger_type = trigger::TriggerPhoton; 00070 00071 edm::Handle<trigger::TriggerFilterObjectWithRefs> PrevFilterOutput; 00072 00073 iEvent.getByLabel (candTag_,PrevFilterOutput); 00074 00075 std::vector<edm::Ref<reco::RecoEcalCandidateCollection> > recoecalcands; 00076 PrevFilterOutput->getObjects(TriggerCluster, recoecalcands); 00077 00078 //get hold of isolated association map 00079 edm::Handle<reco::RecoEcalCandidateIsolationMap> depMap; 00080 iEvent.getByLabel (isoTag_,depMap); 00081 00082 //get hold of non-isolated association map 00083 edm::Handle<reco::RecoEcalCandidateIsolationMap> depNonIsoMap; 00084 if(!doIsolated_) iEvent.getByLabel (nonIsoTag_,depNonIsoMap); 00085 00086 // look at all photons, check cuts and add to filter object 00087 int n = 0; 00088 00089 for (unsigned int i=0; i<recoecalcands.size(); i++) { 00090 00091 ref = recoecalcands[i]; 00092 reco::RecoEcalCandidateIsolationMap::const_iterator mapi = (*depMap).find( ref ); 00093 if (mapi==(*depMap).end() && !doIsolated_) mapi = (*depNonIsoMap).find( ref ); 00094 00095 float vali = mapi->val; 00096 float energy = ref->superCluster()->energy(); 00097 float EtaSC = ref->eta(); 00098 if (useEt_) energy = energy * sin (2*atan(exp(-EtaSC))); 00099 if (energy < 0.) energy=0.; /* first and second order terms assume non-negative energies */ 00100 00101 if ( lessThan_ ) { 00102 if ((fabs(EtaSC) < 1.479 && vali <= thrRegularEB_ + energy*thrOverEEB_ + energy*energy*thrOverE2EB_) || (fabs(EtaSC) >= 1.479 && vali <= thrRegularEE_ + energy*thrOverEEE_ + energy*energy*thrOverE2EE_) ) { 00103 n++; 00104 filterproduct->addObject(trigger_type, ref); 00105 continue; 00106 } 00107 } else { 00108 if ((fabs(EtaSC) < 1.479 && vali >= thrRegularEB_ + energy*thrOverEEB_ + energy*energy*thrOverE2EB_) || (fabs(EtaSC) >= 1.479 && vali >= thrRegularEE_ + energy*thrOverEEE_ + energy*energy*thrOverE2EE_) ) { 00109 n++; 00110 filterproduct->addObject(trigger_type, ref); 00111 continue; 00112 } 00113 } 00114 } 00115 00116 // filter decision 00117 bool accept(n>=ncandcut_); 00118 00119 // put filter object into the Event 00120 iEvent.put(filterproduct); 00121 00122 return accept; 00123 } 00124