CMS 3D CMS Logo

List of all members | Public Member Functions | Public Attributes | Private Attributes
HTMHTAnalyzer Class Reference

#include <HTMHTAnalyzer.h>

Inheritance diagram for HTMHTAnalyzer:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &, const edm::TriggerResults &)
 Get the analysis. More...
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 
 HTMHTAnalyzer (const edm::ParameterSet &)
 Constructor. More...
 
 ~HTMHTAnalyzer () override
 Destructor. More...
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

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

Author
K. Hatakeyama, Rockefeller University

Definition at line 31 of file HTMHTAnalyzer.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 29 of file HTMHTAnalyzer.cc.

29  {
30 
31  parameters = pSet;
32  _ptThreshold = 30.;
33 
34 }
double _ptThreshold
Definition: HTMHTAnalyzer.h:66
HTMHTAnalyzer::~HTMHTAnalyzer ( )
override

Destructor.

Definition at line 37 of file HTMHTAnalyzer.cc.

37 { }

Member Function Documentation

void HTMHTAnalyzer::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup,
const edm::TriggerResults triggerResults 
)

Get the analysis.

Definition at line 75 of file HTMHTAnalyzer.cc.

References edm::HLTGlobalStatus::accept(), gather_cfg::cout, edm::Event::getByLabel(), mps_fire::i, edm::HandleBase::isValid(), LogTrace, metname, nanoDQM_cfi::MHT, njet, edm::HLTGlobalStatus::size(), edm::TriggerNames::triggerIndex(), and edm::Event::triggerNames().

76  {
77 
78  LogTrace(metname)<<"[HTMHTAnalyzer] Analyze HT & MHT";
79 
80  jetME->Fill(4);
81 
82  // ==========================================================
83  // Trigger information
84  //
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_.empty()) _trig_JetMB=triggerResults.size()-1;
112 
113  if (_trig_JetMB==0) return;
114 
115  // ==========================================================
116 
117  // **** Get the Calo Jet container
119 
120  // **** Get the SISCone Jet container
121  iEvent.getByLabel(theJetCollectionForHTMHTLabel, jetcoll);
122 
123  if(!jetcoll.isValid()) return;
124 
125  // ==========================================================
126  // Reconstructed HT & MHT Information
127 
128  int njet=0;
129  double MHx=0.;
130  double MHy=0.;
131  double MHT=0.;
132  double MHTPhi=0.;
133  double HT=0.;
134 
135  for (reco::CaloJetCollection::const_iterator calojet = jetcoll->begin(); calojet!=jetcoll->end(); ++calojet){
136  if (calojet->pt()>_ptThreshold){
137  njet++;
138  MHx += -1.*calojet->px();
139  MHy += -1.*calojet->py();
140  HT += calojet->pt();
141  }
142  }
143 
144  TVector2 *tv2 = new TVector2(MHx,MHy);
145  MHT =tv2->Mod();
146  MHTPhi=tv2->Phi();
147 
148  //std::cout << "HTMHT " << MHT << std::endl;
149 
150  hNevents->Fill(1.);
151  hNJets->Fill(njet);
152  if (njet>0){
153  hMHx->Fill(MHx);
154  hMHy->Fill(MHy);
155  hMHT->Fill(MHT);
156  hMHTPhi->Fill(MHTPhi);
157  hHT->Fill(HT);
158  }
159 
160 }
MonitorElement * hHT
Definition: HTMHTAnalyzer.h:80
MonitorElement * hMHT
Definition: HTMHTAnalyzer.h:77
std::string metname
Definition: HTMHTAnalyzer.h:54
bool accept() const
Has at least one path accepted the event?
MonitorElement * hMHTPhi
Definition: HTMHTAnalyzer.h:78
int njet
Definition: HydjetWrapper.h:95
void Fill(long long x)
unsigned int triggerIndex(std::string const &name) const
Definition: TriggerNames.cc:32
double _ptThreshold
Definition: HTMHTAnalyzer.h:66
unsigned int size() const
Get number of paths stored.
MonitorElement * hNevents
Definition: HTMHTAnalyzer.h:71
std::vector< std::string > HLTPathsJetMBByName_
Definition: HTMHTAnalyzer.h:61
bool isValid() const
Definition: HandleBase.h:74
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:535
#define LogTrace(id)
MonitorElement * jetME
Definition: HTMHTAnalyzer.h:69
MonitorElement * hMHy
Definition: HTMHTAnalyzer.h:76
MonitorElement * hNJets
Definition: HTMHTAnalyzer.h:73
MonitorElement * hMHx
Definition: HTMHTAnalyzer.h:75
Definition: HT.h:21
edm::InputTag theJetCollectionForHTMHTLabel
Definition: HTMHTAnalyzer.h:58
edm::TriggerNames const & triggerNames(edm::TriggerResults const &triggerResults) const override
Definition: Event.cc:301
void HTMHTAnalyzer::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &   
)
override

Definition at line 41 of file HTMHTAnalyzer.cc.

References DQMStore::IBooker::book1D(), LogTrace, metname, MonitorElement::setBinLabel(), DQMStore::IBooker::setCurrentFolder(), and AlCaHLTBitMon_QueryRunRegistry::string.

43  {
44  evtCounter = 0;
45  metname = "HTMHTAnalyzer";
46 
47  // PFMET information
48  theJetCollectionForHTMHTLabel = parameters.getParameter<edm::InputTag>("JetCollectionForHTMHTLabel");
49  _source = parameters.getParameter<std::string>("Source");
50 
51  LogTrace(metname)<<"[HTMHTAnalyzer] Parameters initialization";
52  ibooker.setCurrentFolder("JetMET/MET/"+_source);
53 
54  HLTPathsJetMBByName_ = parameters.getParameter<std::vector<std::string > >("HLTPathsJetMB");
55 
56  // misc
57  _verbose = parameters.getParameter<int>("verbose");
58  _ptThreshold = parameters.getParameter<double>("ptThreshold");
59 
60  jetME = ibooker.book1D("metReco", "metReco", 4, 1, 5);
61  jetME->setBinLabel(4,"HTMHT",1);
62 
63  hNevents = ibooker.book1D("METTask_Nevents", "METTask_Nevents", 1,0,1);
64  hNJets = ibooker.book1D("METTask_NJets", "METTask_NJets", 100, 0, 100);
65  hMHx = ibooker.book1D("METTask_MHx", "METTask_MHx", 500,-500,500);
66  hMHy = ibooker.book1D("METTask_MHy", "METTask_MHy", 500,-500,500);
67  hMHT = ibooker.book1D("METTask_MHT", "METTask_MHT", 500,0,1000);
68  hMHTPhi = ibooker.book1D("METTask_MhTPhi", "METTask_MhTPhi", 80,-4,4);
69  hHT = ibooker.book1D("METTask_HT", "METTask_HT", 500,0,2000);
70 
71 }
MonitorElement * hHT
Definition: HTMHTAnalyzer.h:80
MonitorElement * hMHT
Definition: HTMHTAnalyzer.h:77
std::string metname
Definition: HTMHTAnalyzer.h:54
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:78
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::string _source
Definition: HTMHTAnalyzer.h:56
double _ptThreshold
Definition: HTMHTAnalyzer.h:66
MonitorElement * hNevents
Definition: HTMHTAnalyzer.h:71
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
std::vector< std::string > HLTPathsJetMBByName_
Definition: HTMHTAnalyzer.h:61
#define LogTrace(id)
MonitorElement * jetME
Definition: HTMHTAnalyzer.h:69
MonitorElement * hMHy
Definition: HTMHTAnalyzer.h:76
MonitorElement * hNJets
Definition: HTMHTAnalyzer.h:73
MonitorElement * hMHx
Definition: HTMHTAnalyzer.h:75
edm::InputTag theJetCollectionForHTMHTLabel
Definition: HTMHTAnalyzer.h:58

Member Data Documentation

double HTMHTAnalyzer::_ptThreshold
private

Definition at line 66 of file HTMHTAnalyzer.h.

std::string HTMHTAnalyzer::_source
private

Definition at line 56 of file HTMHTAnalyzer.h.

int HTMHTAnalyzer::_trig_JetMB
private

Definition at line 63 of file HTMHTAnalyzer.h.

int HTMHTAnalyzer::_verbose
private
int HTMHTAnalyzer::evtCounter

Definition at line 45 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::hHT
private

Definition at line 80 of file HTMHTAnalyzer.h.

std::vector<std::string > HTMHTAnalyzer::HLTPathsJetMBByName_
private

Definition at line 61 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::hMHT
private

Definition at line 77 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::hMHTPhi
private

Definition at line 78 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::hMHx
private

Definition at line 75 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::hMHy
private

Definition at line 76 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::hNevents
private

Definition at line 71 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::hNJets
private

Definition at line 73 of file HTMHTAnalyzer.h.

MonitorElement* HTMHTAnalyzer::jetME
private

Definition at line 69 of file HTMHTAnalyzer.h.

std::string HTMHTAnalyzer::metname
private

Definition at line 54 of file HTMHTAnalyzer.h.

edm::ParameterSet HTMHTAnalyzer::parameters
private
edm::InputTag HTMHTAnalyzer::theJetCollectionForHTMHTLabel
private

Definition at line 58 of file HTMHTAnalyzer.h.