#include <HTMHTAnalyzer.h>
Public Member Functions | |
void | analyze (const edm::Event &, const edm::EventSetup &, const edm::TriggerResults &) |
Get the analysis. | |
void | beginJob (DQMStore *dbe) |
Inizialize parameters for histo binning. | |
HTMHTAnalyzer (const edm::ParameterSet &) | |
Constructor. | |
virtual | ~HTMHTAnalyzer () |
Destructor. | |
Public Attributes | |
int | evtCounter |
Private Attributes | |
double | _ptThreshold |
std::string | _source |
int | _trig_JetMB |
int | _verbose |
MonitorElement * | hHT |
std::vector< std::string > | HLTPathsJetMBByName_ |
MonitorElement * | hMHT |
MonitorElement * | hMHTPhi |
MonitorElement * | hMHx |
MonitorElement * | hMHy |
MonitorElement * | hNevents |
MonitorElement * | hNJets |
MonitorElement * | jetME |
std::string | metname |
edm::ParameterSet | parameters |
edm::InputTag | theJetCollectionForHTMHTLabel |
DQM monitoring source for HTMHT
Definition at line 33 of file HTMHTAnalyzer.h.
HTMHTAnalyzer::HTMHTAnalyzer | ( | const edm::ParameterSet & | pSet | ) |
Constructor.
Definition at line 31 of file HTMHTAnalyzer.cc.
References ExpressReco_HICollisions_FallBack::parameters.
{ parameters = pSet; _ptThreshold = 30.; }
HTMHTAnalyzer::~HTMHTAnalyzer | ( | ) | [virtual] |
void HTMHTAnalyzer::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup, | ||
const edm::TriggerResults & | triggerResults | ||
) |
Get the analysis.
Definition at line 74 of file HTMHTAnalyzer.cc.
References edm::HLTGlobalStatus::accept(), gather_cfg::cout, edm::Event::getByLabel(), i, edm::HandleBase::isValid(), LogTrace, metname, njet, edm::HLTGlobalStatus::size(), edm::TriggerNames::triggerIndex(), and edm::Event::triggerNames().
{ LogTrace(metname)<<"[HTMHTAnalyzer] Analyze HT & MHT"; jetME->Fill(4); // ========================================================== // Trigger information // if(&triggerResults) { // // // Check how many HLT triggers are in triggerResults int ntrigs = triggerResults.size(); if (_verbose) std::cout << "ntrigs=" << ntrigs << std::endl; // // // If index=ntrigs, this HLT trigger doesn't exist in the HLT table for this data. const edm::TriggerNames & triggerNames = iEvent.triggerNames(triggerResults); // // // count number of requested Jet or MB HLT paths which have fired for (unsigned int i=0; i!=HLTPathsJetMBByName_.size(); i++) { unsigned int triggerIndex = triggerNames.triggerIndex(HLTPathsJetMBByName_[i]); if (triggerIndex<triggerResults.size()) { if (triggerResults.accept(triggerIndex)) { _trig_JetMB++; } } } // for empty input vectors (n==0), take all HLT triggers! if (HLTPathsJetMBByName_.size()==0) _trig_JetMB=triggerResults.size()-1; if (_trig_JetMB==0) return; } else { edm::LogInfo("CaloMetAnalyzer") << "TriggerResults::HLT not found, " "automatically select events"; //return; } // ========================================================== // **** Get the Calo Jet container edm::Handle<reco::CaloJetCollection> jetcoll; // **** Get the SISCone Jet container iEvent.getByLabel(theJetCollectionForHTMHTLabel, jetcoll); if(!jetcoll.isValid()) return; // ========================================================== // Reconstructed HT & MHT Information int njet=0; double MHx=0.; double MHy=0.; double MHT=0.; double MHTPhi=0.; double HT=0.; for (reco::CaloJetCollection::const_iterator calojet = jetcoll->begin(); calojet!=jetcoll->end(); ++calojet){ if (calojet->pt()>_ptThreshold){ njet++; MHx += -1.*calojet->px(); MHy += -1.*calojet->py(); HT += calojet->pt(); } } TVector2 *tv2 = new TVector2(MHx,MHy); MHT =tv2->Mod(); MHTPhi=tv2->Phi(); //std::cout << "HTMHT " << MHT << std::endl; hNevents->Fill(1.); hNJets->Fill(njet); if (njet>0){ hMHx->Fill(MHx); hMHy->Fill(MHy); hMHT->Fill(MHT); hMHTPhi->Fill(MHTPhi); hHT->Fill(HT); } }
void HTMHTAnalyzer::beginJob | ( | DQMStore * | dbe | ) | [virtual] |
Inizialize parameters for histo binning.
Implements JetAnalyzerBase.
Definition at line 41 of file HTMHTAnalyzer.cc.
References DQMStore::book1D(), LogTrace, metname, ExpressReco_HICollisions_FallBack::parameters, MonitorElement::setBinLabel(), and DQMStore::setCurrentFolder().
{ evtCounter = 0; metname = "HTMHTAnalyzer"; // PFMET information theJetCollectionForHTMHTLabel = parameters.getParameter<edm::InputTag>("JetCollectionForHTMHTLabel"); _source = parameters.getParameter<std::string>("Source"); LogTrace(metname)<<"[HTMHTAnalyzer] Parameters initialization"; dbe->setCurrentFolder("JetMET/MET/"+_source); HLTPathsJetMBByName_ = parameters.getParameter<std::vector<std::string > >("HLTPathsJetMB"); // misc _verbose = parameters.getParameter<int>("verbose"); _ptThreshold = parameters.getParameter<double>("ptThreshold"); jetME = dbe->book1D("metReco", "metReco", 4, 1, 5); jetME->setBinLabel(4,"HTMHT",1); hNevents = dbe->book1D("METTask_Nevents", "METTask_Nevents", 1,0,1); hNJets = dbe->book1D("METTask_NJets", "METTask_NJets", 100, 0, 100); hMHx = dbe->book1D("METTask_MHx", "METTask_MHx", 500,-500,500); hMHy = dbe->book1D("METTask_MHy", "METTask_MHy", 500,-500,500); hMHT = dbe->book1D("METTask_MHT", "METTask_MHT", 500,0,1000); hMHTPhi = dbe->book1D("METTask_MhTPhi", "METTask_MhTPhi", 80,-4,4); hHT = dbe->book1D("METTask_HT", "METTask_HT", 500,0,2000); }
double HTMHTAnalyzer::_ptThreshold [private] |
Definition at line 70 of file HTMHTAnalyzer.h.
std::string HTMHTAnalyzer::_source [private] |
Definition at line 60 of file HTMHTAnalyzer.h.
int HTMHTAnalyzer::_trig_JetMB [private] |
Definition at line 67 of file HTMHTAnalyzer.h.
int HTMHTAnalyzer::_verbose [private] |
Definition at line 56 of file HTMHTAnalyzer.h.
Definition at line 49 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::hHT [private] |
Definition at line 84 of file HTMHTAnalyzer.h.
std::vector<std::string > HTMHTAnalyzer::HLTPathsJetMBByName_ [private] |
Definition at line 65 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::hMHT [private] |
Definition at line 81 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::hMHTPhi [private] |
Definition at line 82 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::hMHx [private] |
Definition at line 79 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::hMHy [private] |
Definition at line 80 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::hNevents [private] |
Definition at line 75 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::hNJets [private] |
Definition at line 77 of file HTMHTAnalyzer.h.
MonitorElement* HTMHTAnalyzer::jetME [private] |
Definition at line 73 of file HTMHTAnalyzer.h.
std::string HTMHTAnalyzer::metname [private] |
Definition at line 58 of file HTMHTAnalyzer.h.
edm::ParameterSet HTMHTAnalyzer::parameters [private] |
Definition at line 54 of file HTMHTAnalyzer.h.
Definition at line 62 of file HTMHTAnalyzer.h.