#include <HLTRHemisphere.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HLTRHemisphere (const edm::ParameterSet &) | |
~HLTRHemisphere () | |
Static Public Member Functions | |
static void | fillDescriptions (edm::ConfigurationDescriptions &descriptions) |
Private Member Functions | |
void | ComputeHemispheres (std::auto_ptr< std::vector< math::XYZTLorentzVector > > &hlist, const std::vector< math::XYZTLorentzVector > &JETS, std::vector< math::XYZTLorentzVector > *extraJets=0) |
Private Attributes | |
bool | accNJJets_ |
bool | doMuonCorrection_ |
edm::InputTag | inputTag_ |
double | max_Eta_ |
int | max_NJ_ |
double | min_Jet_Pt_ |
double | muonEta_ |
edm::InputTag | muonTag_ |
Definition at line 18 of file HLTRHemisphere.h.
HLTRHemisphere::HLTRHemisphere | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 37 of file HLTRHemisphere.cc.
References accNJJets_, edm::InputTag::encode(), inputTag_, LogDebug, max_Eta_, max_NJ_, and min_Jet_Pt_.
: inputTag_ (iConfig.getParameter<edm::InputTag>("inputTag")), muonTag_ (iConfig.getParameter<edm::InputTag>("muonTag")), doMuonCorrection_(iConfig.getParameter<bool> ("doMuonCorrection" )), muonEta_ (iConfig.getParameter<double> ("maxMuonEta" )), min_Jet_Pt_ (iConfig.getParameter<double> ("minJetPt" )), max_Eta_ (iConfig.getParameter<double> ("maxEta" )), max_NJ_ (iConfig.getParameter<int> ("maxNJ" )), accNJJets_ (iConfig.getParameter<bool> ("acceptNJ" )) { LogDebug("") << "Input/minJetPt/maxEta/maxNJ/acceptNJ : " << inputTag_.encode() << " " << min_Jet_Pt_ << "/" << max_Eta_ << "/" << max_NJ_ << "/" << accNJJets_ << "."; //register your products produces<std::vector<math::XYZTLorentzVector> >(); }
HLTRHemisphere::~HLTRHemisphere | ( | ) |
Definition at line 58 of file HLTRHemisphere.cc.
{ }
void HLTRHemisphere::ComputeHemispheres | ( | std::auto_ptr< std::vector< math::XYZTLorentzVector > > & | hlist, |
const std::vector< math::XYZTLorentzVector > & | JETS, | ||
std::vector< math::XYZTLorentzVector > * | extraJets = 0 |
||
) | [private] |
Definition at line 165 of file HLTRHemisphere.cc.
References prof2calltree::count, i, funct::pow(), and dt_dqm_sourceclient_common_cff::reco.
Referenced by filter().
{ using namespace math; using namespace reco; XYZTLorentzVector j1R(0.1, 0., 0., 0.1); XYZTLorentzVector j2R(0.1, 0., 0., 0.1); int nJets = JETS.size(); if(extraJets) nJets+=extraJets->size(); if(nJets<2){ // put empty hemispheres if not enough jets hlist->push_back(j1R); hlist->push_back(j2R); return; } unsigned int N_comb = pow(2,nJets); // compute the number of combinations of jets possible //Make the hemispheres double M_minR = 9999999999.0; unsigned int j_count; for (unsigned int i = 0; i < N_comb; i++) { XYZTLorentzVector j_temp1, j_temp2; unsigned int itemp = i; j_count = N_comb/2; unsigned int count = 0; while (j_count > 0) { if (itemp/j_count == 1){ if(count<JETS.size()) j_temp1 += JETS.at(count); else j_temp1 +=extraJets->at(count-JETS.size()); } else { if(count<JETS.size()) j_temp2 += JETS.at(count); else j_temp2 +=extraJets->at(count-JETS.size()); } itemp -= j_count * (itemp/j_count); j_count /= 2; count++; } double M_temp = j_temp1.M2() + j_temp2.M2(); if (M_temp < M_minR) { M_minR = M_temp; j1R = j_temp1; j2R = j_temp2; } } hlist->push_back(j1R); hlist->push_back(j2R); return; }
void HLTRHemisphere::fillDescriptions | ( | edm::ConfigurationDescriptions & | descriptions | ) | [static] |
Reimplemented from edm::EDFilter.
Definition at line 63 of file HLTRHemisphere.cc.
References edm::ParameterSetDescription::add(), and edm::ConfigurationDescriptions::add().
{ edm::ParameterSetDescription desc; desc.add<edm::InputTag>("inputTag",edm::InputTag("hltMCJetCorJetIcone5HF07")); desc.add<edm::InputTag>("muonTag",edm::InputTag("")); desc.add<bool>("doMuonCorrection",false); desc.add<double>("maxMuonEta",2.1); desc.add<double>("minJetPt",30.0); desc.add<double>("maxEta",3.0); desc.add<int>("maxNJ",7); desc.add<bool>("acceptNJ",true); descriptions.add("hltRHemisphere",desc); }
bool HLTRHemisphere::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDFilter.
Definition at line 82 of file HLTRHemisphere.cc.
References abs, accNJJets_, ComputeHemispheres(), doMuonCorrection_, edm::Event::getByLabel(), i, getHLTprescales::index, inputTag_, fwrapper::jets, max_Eta_, max_NJ_, min_Jet_Pt_, muonEta_, patZpeak::muons, muonTag_, n, reco::LeafCandidate::p4(), edm::Event::put(), and dt_dqm_sourceclient_common_cff::reco.
{ using namespace std; using namespace edm; using namespace reco; using namespace math; using namespace trigger; typedef XYZTLorentzVector LorentzVector; // get hold of collection of objects // Handle<CaloJetCollection> jets; Handle<View<Jet> > jets; iEvent.getByLabel (inputTag_,jets); // get hold of the muons, if necessary Handle<vector<reco::RecoChargedCandidate> > muons; if(doMuonCorrection_) iEvent.getByLabel( muonTag_,muons ); // The output Collection std::auto_ptr<vector<math::XYZTLorentzVector> > Hemispheres(new vector<math::XYZTLorentzVector> ); // look at all objects, check cuts and add to filter object int n(0); vector<math::XYZTLorentzVector> JETS; for (unsigned int i=0; i<jets->size(); i++) { if(std::abs(jets->at(i).eta()) < max_Eta_ && jets->at(i).pt() >= min_Jet_Pt_){ JETS.push_back(jets->at(i).p4()); n++; } } if(n>max_NJ_ && max_NJ_!=-1){ iEvent.put(Hemispheres); return accNJJets_; // too many jets, accept for timing } if(doMuonCorrection_){ const int nMu = 2; int muonIndex[nMu] = { -1, -1 }; std::vector<reco::RecoChargedCandidate>::const_iterator muonIt; int index = 0; int nPassMu = 0; for(muonIt = muons->begin(); muonIt!=muons->end(); muonIt++,index++){ if(std::abs(muonIt->eta()) > muonEta_ || muonIt->pt() < min_Jet_Pt_) continue; // skip muons out of eta range or too low pT if(nPassMu >= 2){ // if we have already accepted two muons, accept the event iEvent.put(Hemispheres); // too many muons, accept for timing return true; } muonIndex[nPassMu++] = index; } //muons as MET this->ComputeHemispheres(Hemispheres,JETS); //lead muon as jet if(nPassMu>0){ std::vector<math::XYZTLorentzVector> muonJets; reco::RecoChargedCandidate leadMu = muons->at(muonIndex[0]); muonJets.push_back(leadMu.p4()); Hemispheres->push_back(leadMu.p4()); this->ComputeHemispheres(Hemispheres,JETS,&muonJets); // lead muon as jet if(nPassMu>1){ // two passing muons muonJets.pop_back(); reco::RecoChargedCandidate secondMu = muons->at(muonIndex[1]); muonJets.push_back(secondMu.p4()); Hemispheres->push_back(secondMu.p4()); this->ComputeHemispheres(Hemispheres,JETS,&muonJets); // lead muon as v, second muon as jet muonJets.push_back(leadMu.p4()); this->ComputeHemispheres(Hemispheres,JETS,&muonJets); // both muon as jets } } }else{ // do MuonCorrection==false if(n<2) return false; // not enough jets and not adding in muons this->ComputeHemispheres(Hemispheres,JETS); // don't do the muon isolation, just run once and done } //Format: // 0 muon: 2 hemispheres (2) // 1 muon: 2 hemisheress + leadMuP4 + 2 hemispheres (5) // 2 muon: 2 hemispheres + leadMuP4 + 2 hemispheres + 2ndMuP4 + 4 Hemispheres (10) iEvent.put(Hemispheres); return true; }
bool HLTRHemisphere::accNJJets_ [private] |
Definition at line 35 of file HLTRHemisphere.h.
Referenced by filter(), and HLTRHemisphere().
bool HLTRHemisphere::doMuonCorrection_ [private] |
Definition at line 30 of file HLTRHemisphere.h.
Referenced by filter().
edm::InputTag HLTRHemisphere::inputTag_ [private] |
Definition at line 28 of file HLTRHemisphere.h.
Referenced by filter(), and HLTRHemisphere().
double HLTRHemisphere::max_Eta_ [private] |
Definition at line 33 of file HLTRHemisphere.h.
Referenced by filter(), and HLTRHemisphere().
int HLTRHemisphere::max_NJ_ [private] |
Definition at line 34 of file HLTRHemisphere.h.
Referenced by filter(), and HLTRHemisphere().
double HLTRHemisphere::min_Jet_Pt_ [private] |
Definition at line 32 of file HLTRHemisphere.h.
Referenced by filter(), and HLTRHemisphere().
double HLTRHemisphere::muonEta_ [private] |
Definition at line 31 of file HLTRHemisphere.h.
Referenced by filter().
edm::InputTag HLTRHemisphere::muonTag_ [private] |
Definition at line 29 of file HLTRHemisphere.h.
Referenced by filter().