Go to the documentation of this file.00001
00008 #include "HLTrigger/JetMET/interface/HLTPhi2METFilter.h"
00009
00010 #include "DataFormats/Common/interface/Handle.h"
00011
00012 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00013
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015
00016 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00017 #include "DataFormats/METReco/interface/CaloMET.h"
00018
00019 #include "FWCore/Framework/interface/ESHandle.h"
00020 #include "FWCore/Framework/interface/EventSetup.h"
00021
00022
00023
00024
00025
00026 HLTPhi2METFilter::HLTPhi2METFilter(const edm::ParameterSet& iConfig)
00027 {
00028 inputJetTag_ = iConfig.getParameter< edm::InputTag > ("inputJetTag");
00029 inputMETTag_ = iConfig.getParameter< edm::InputTag > ("inputMETTag");
00030 saveTags_ = iConfig.getUntrackedParameter<bool>("saveTags",false);
00031 minDPhi_ = iConfig.getParameter<double> ("minDeltaPhi");
00032 maxDPhi_ = iConfig.getParameter<double> ("maxDeltaPhi");
00033 minEtjet1_= iConfig.getParameter<double> ("minEtJet1");
00034 minEtjet2_= iConfig.getParameter<double> ("minEtJet2");
00035
00036
00037 produces<trigger::TriggerFilterObjectWithRefs>();
00038 }
00039
00040 HLTPhi2METFilter::~HLTPhi2METFilter(){}
00041
00042
00043
00044 bool
00045 HLTPhi2METFilter::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00046 {
00047 using namespace std;
00048 using namespace edm;
00049 using namespace reco;
00050 using namespace trigger;
00051
00052
00053 auto_ptr<trigger::TriggerFilterObjectWithRefs>
00054 filterobject (new trigger::TriggerFilterObjectWithRefs(path(),module()));
00055 if (saveTags_) {
00056 filterobject->addCollectionTag(inputJetTag_);
00057 filterobject->addCollectionTag(inputMETTag_);
00058 }
00059
00060 Handle<CaloJetCollection> recocalojets;
00061 iEvent.getByLabel(inputJetTag_,recocalojets);
00062 Handle<trigger::TriggerFilterObjectWithRefs> metcal;
00063 iEvent.getByLabel(inputMETTag_,metcal);
00064
00065
00066 int n(0);
00067
00068 VRcalomet vrefMET;
00069 metcal->getObjects(TriggerMET,vrefMET);
00070 CaloMETRef metRef=vrefMET.at(0);
00071 CaloJetRef ref1,ref2;
00072
00073 if(recocalojets->size() > 1){
00074
00075
00076 double etjet1=0.;
00077 double etjet2=0.;
00078 double phijet2=0.;
00079
00080 double phimiss = vrefMET.at(0)->phi();
00081 int countjets =0;
00082
00083 for (CaloJetCollection::const_iterator recocalojet = recocalojets->begin();
00084 recocalojet<=(recocalojets->begin()+1); recocalojet++) {
00085
00086 if(countjets==0) {
00087 etjet1 = recocalojet->et();
00088 ref1 = CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet));
00089 }
00090 if(countjets==1) {
00091 etjet2 = recocalojet->et();
00092 phijet2 = recocalojet->phi();
00093 ref2 = CaloJetRef(recocalojets,distance(recocalojets->begin(),recocalojet));
00094 }
00095 countjets++;
00096 }
00097 double Dphi= fabs(phimiss-phijet2);
00098 if (Dphi>M_PI) Dphi=2.0*M_PI-Dphi;
00099 if(etjet1>minEtjet1_ && etjet2>minEtjet2_ && Dphi>=minDPhi_ && Dphi<=maxDPhi_){
00100 filterobject->addObject(TriggerMET,metRef);
00101 filterobject->addObject(TriggerJet,ref1);
00102 filterobject->addObject(TriggerJet,ref2);
00103 n++;
00104 }
00105
00106 }
00107
00108
00109
00110
00111 bool accept(n>=1);
00112
00113
00114 iEvent.put(filterobject);
00115
00116 return accept;
00117 }