CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Public Attributes | Private Attributes
HTMHTAnalyzer Class Reference

#include <HTMHTAnalyzer.h>

Inheritance diagram for HTMHTAnalyzer:
JetAnalyzerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &, const edm::TriggerResults &)
 Get the analysis. More...
 
void beginJob (DQMStore *dbe)
 Inizialize parameters for histo binning. More...
 
 HTMHTAnalyzer (const edm::ParameterSet &)
 Constructor. More...
 
virtual ~HTMHTAnalyzer ()
 Destructor. More...
 
- Public Member Functions inherited from JetAnalyzerBase
void analyze (const edm::Event &, const edm::EventSetup &, reco::CaloJet &caloJet)
 Get the analysis of the muon properties. More...
 
 JetAnalyzerBase ()
 Constructor. More...
 
virtual ~JetAnalyzerBase ()
 Destructor. More...
 

Public Attributes

int evtCounter
 

Private Attributes

double _ptThreshold
 
std::string _source
 
int _trig_JetMB
 
int _verbose
 
MonitorElementhHT
 
std::vector< std::string > HLTPathsJetMBByName_
 
MonitorElementhMHT
 
MonitorElementhMHTPhi
 
MonitorElementhMHx
 
MonitorElementhMHy
 
MonitorElementhNevents
 
MonitorElementhNJets
 
MonitorElementjetME
 
std::string metname
 
edm::ParameterSet parameters
 
edm::InputTag theJetCollectionForHTMHTLabel
 

Detailed Description

DQM monitoring source for HTMHT

Date:
2010/02/24 19:08:53
Revision:
1.6
Author
K. Hatakeyama, Rockefeller University

Definition at line 33 of file HTMHTAnalyzer.h.

Constructor & Destructor Documentation

HTMHTAnalyzer::HTMHTAnalyzer ( const edm::ParameterSet pSet)

Constructor.

Definition at line 31 of file HTMHTAnalyzer.cc.

References Parameters::parameters.

31  {
32 
33  parameters = pSet;
34  _ptThreshold = 30.;
35 
36 }
edm::ParameterSet parameters
Definition: HTMHTAnalyzer.h:54
double _ptThreshold
Definition: HTMHTAnalyzer.h:70
HTMHTAnalyzer::~HTMHTAnalyzer ( )
virtual

Destructor.

Definition at line 39 of file HTMHTAnalyzer.cc.

39 { }

Member Function Documentation

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().

75  {
76 
77  LogTrace(metname)<<"[HTMHTAnalyzer] Analyze HT & MHT";
78 
79  jetME->Fill(4);
80 
81  // ==========================================================
82  // Trigger information
83  //
84  if(&triggerResults) {
85 
87 
88  //
89  //
90  // Check how many HLT triggers are in triggerResults
91  int ntrigs = triggerResults.size();
92  if (_verbose) std::cout << "ntrigs=" << ntrigs << std::endl;
93 
94  //
95  //
96  // If index=ntrigs, this HLT trigger doesn't exist in the HLT table for this data.
97  const edm::TriggerNames & triggerNames = iEvent.triggerNames(triggerResults);
98 
99  //
100  //
101  // count number of requested Jet or MB HLT paths which have fired
102  for (unsigned int i=0; i!=HLTPathsJetMBByName_.size(); i++) {
103  unsigned int triggerIndex = triggerNames.triggerIndex(HLTPathsJetMBByName_[i]);
104  if (triggerIndex<triggerResults.size()) {
105  if (triggerResults.accept(triggerIndex)) {
106  _trig_JetMB++;
107  }
108  }
109  }
110  // for empty input vectors (n==0), take all HLT triggers!
111  if (HLTPathsJetMBByName_.size()==0) _trig_JetMB=triggerResults.size()-1;
112 
113  if (_trig_JetMB==0) return;
114 
115  } else {
116 
117  edm::LogInfo("CaloMetAnalyzer") << "TriggerResults::HLT not found, "
118  "automatically select events";
119  //return;
120 
121  }
122 
123  // ==========================================================
124 
125  // **** Get the Calo Jet container
127 
128  // **** Get the SISCone Jet container
129  iEvent.getByLabel(theJetCollectionForHTMHTLabel, jetcoll);
130 
131  if(!jetcoll.isValid()) return;
132 
133  // ==========================================================
134  // Reconstructed HT & MHT Information
135 
136  int njet=0;
137  double MHx=0.;
138  double MHy=0.;
139  double MHT=0.;
140  double MHTPhi=0.;
141  double HT=0.;
142 
143  for (reco::CaloJetCollection::const_iterator calojet = jetcoll->begin(); calojet!=jetcoll->end(); ++calojet){
144  if (calojet->pt()>_ptThreshold){
145  njet++;
146  MHx += -1.*calojet->px();
147  MHy += -1.*calojet->py();
148  HT += calojet->pt();
149  }
150  }
151 
152  TVector2 *tv2 = new TVector2(MHx,MHy);
153  MHT =tv2->Mod();
154  MHTPhi=tv2->Phi();
155 
156  //std::cout << "HTMHT " << MHT << std::endl;
157 
158  hNevents->Fill(1.);
159  hNJets->Fill(njet);
160  if (njet>0){
161  hMHx->Fill(MHx);
162  hMHy->Fill(MHy);
163  hMHT->Fill(MHT);
164  hMHTPhi->Fill(MHTPhi);
165  hHT->Fill(HT);
166  }
167 
168 }
MonitorElement * hHT
Definition: HTMHTAnalyzer.h:84
int i
Definition: DBlmapReader.cc:9
virtual edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const
Definition: Event.cc:207
MonitorElement * hMHT
Definition: HTMHTAnalyzer.h:81
std::string metname
Definition: HTMHTAnalyzer.h:58
bool accept() const
Has at least one path accepted the event?
MonitorElement * hMHTPhi
Definition: HTMHTAnalyzer.h:82
int njet
Definition: HydjetWrapper.h:91
void Fill(long long x)
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:32
double _ptThreshold
Definition: HTMHTAnalyzer.h:70
unsigned int size() const
Get number of paths stored.
MonitorElement * hNevents
Definition: HTMHTAnalyzer.h:75
std::vector< std::string > HLTPathsJetMBByName_
Definition: HTMHTAnalyzer.h:65
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
#define LogTrace(id)
MonitorElement * jetME
Definition: HTMHTAnalyzer.h:73
MonitorElement * hMHy
Definition: HTMHTAnalyzer.h:80
MonitorElement * hNJets
Definition: HTMHTAnalyzer.h:77
tuple cout
Definition: gather_cfg.py:121
MonitorElement * hMHx
Definition: HTMHTAnalyzer.h:79
Definition: HT.h:20
edm::InputTag theJetCollectionForHTMHTLabel
Definition: HTMHTAnalyzer.h:62
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, Parameters::parameters, MonitorElement::setBinLabel(), DQMStore::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

41  {
42 
43  evtCounter = 0;
44  metname = "HTMHTAnalyzer";
45 
46  // PFMET information
47  theJetCollectionForHTMHTLabel = parameters.getParameter<edm::InputTag>("JetCollectionForHTMHTLabel");
49 
50  LogTrace(metname)<<"[HTMHTAnalyzer] Parameters initialization";
51  dbe->setCurrentFolder("JetMET/MET/"+_source);
52 
53  HLTPathsJetMBByName_ = parameters.getParameter<std::vector<std::string > >("HLTPathsJetMB");
54 
55  // misc
56  _verbose = parameters.getParameter<int>("verbose");
57  _ptThreshold = parameters.getParameter<double>("ptThreshold");
58 
59  jetME = dbe->book1D("metReco", "metReco", 4, 1, 5);
60  jetME->setBinLabel(4,"HTMHT",1);
61 
62  hNevents = dbe->book1D("METTask_Nevents", "METTask_Nevents", 1,0,1);
63  hNJets = dbe->book1D("METTask_NJets", "METTask_NJets", 100, 0, 100);
64  hMHx = dbe->book1D("METTask_MHx", "METTask_MHx", 500,-500,500);
65  hMHy = dbe->book1D("METTask_MHy", "METTask_MHy", 500,-500,500);
66  hMHT = dbe->book1D("METTask_MHT", "METTask_MHT", 500,0,1000);
67  hMHTPhi = dbe->book1D("METTask_MhTPhi", "METTask_MhTPhi", 80,-4,4);
68  hHT = dbe->book1D("METTask_HT", "METTask_HT", 500,0,2000);
69 
70 }
MonitorElement * hHT
Definition: HTMHTAnalyzer.h:84
T getParameter(std::string const &) const
MonitorElement * hMHT
Definition: HTMHTAnalyzer.h:81
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:722
std::string metname
Definition: HTMHTAnalyzer.h:58
void setBinLabel(int bin, const std::string &label, int axis=1)
set bin label for x, y or z axis (axis=1, 2, 3 respectively)
MonitorElement * hMHTPhi
Definition: HTMHTAnalyzer.h:82
edm::ParameterSet parameters
Definition: HTMHTAnalyzer.h:54
std::string _source
Definition: HTMHTAnalyzer.h:60
double _ptThreshold
Definition: HTMHTAnalyzer.h:70
MonitorElement * hNevents
Definition: HTMHTAnalyzer.h:75
std::vector< std::string > HLTPathsJetMBByName_
Definition: HTMHTAnalyzer.h:65
#define LogTrace(id)
MonitorElement * jetME
Definition: HTMHTAnalyzer.h:73
MonitorElement * hMHy
Definition: HTMHTAnalyzer.h:80
MonitorElement * hNJets
Definition: HTMHTAnalyzer.h:77
MonitorElement * hMHx
Definition: HTMHTAnalyzer.h:79
edm::InputTag theJetCollectionForHTMHTLabel
Definition: HTMHTAnalyzer.h:62
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434

Member Data Documentation

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
int HTMHTAnalyzer::evtCounter

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
edm::InputTag HTMHTAnalyzer::theJetCollectionForHTMHTLabel
private

Definition at line 62 of file HTMHTAnalyzer.h.