CMS 3D CMS Logo

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

#include <QcdHighPtDQM.h>

Inheritance diagram for QcdHighPtDQM:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &, const edm::EventSetup &)
 Get the analysis. More...
 
void beginJob ()
 Inizialize parameters for histo binning. More...
 
void endJob (void)
 
 QcdHighPtDQM (const edm::ParameterSet &)
 Constructor. More...
 
virtual ~QcdHighPtDQM ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndex indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndex > &) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Member Functions

float moverl (const reco::CaloMETCollection &metcollection, float &ljpt)
 
float movers (const reco::CaloMETCollection &metcollection)
 

Private Attributes

edm::InputTag jetLabel_
 
std::map< std::string,
MonitorElement * > 
MEcontainer_
 
edm::InputTag metLabel1_
 
edm::InputTag metLabel2_
 
edm::InputTag metLabel3_
 
edm::InputTag metLabel4_
 
DQMStoretheDbe
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
CurrentProcessingContext const * currentContext () const
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

DQM Physics Module for High Pt QCD group

Based on DQM/SiPixel and DQM/Physics code Version 1.0, 7/7/09 By Keith Rose

Definition at line 26 of file QcdHighPtDQM.h.

Constructor & Destructor Documentation

QcdHighPtDQM::QcdHighPtDQM ( const edm::ParameterSet iConfig)

Constructor.

Definition at line 35 of file QcdHighPtDQM.cc.

References cppFunctionSkipper::operator, and theDbe.

35  :
36  jetLabel_(iConfig.getUntrackedParameter<edm::InputTag>("jetTag")),
37  metLabel1_(iConfig.getUntrackedParameter<edm::InputTag>("metTag1")),
38  metLabel2_(iConfig.getUntrackedParameter<edm::InputTag>("metTag2")),
39  metLabel3_(iConfig.getUntrackedParameter<edm::InputTag>("metTag3")),
41 {
42 
44 
45 }
T getUntrackedParameter(std::string const &, T const &) const
edm::InputTag metLabel2_
Definition: QcdHighPtDQM.h:55
DQMStore * theDbe
Definition: QcdHighPtDQM.h:50
edm::InputTag jetLabel_
Definition: QcdHighPtDQM.h:53
edm::InputTag metLabel1_
Definition: QcdHighPtDQM.h:54
edm::InputTag metLabel3_
Definition: QcdHighPtDQM.h:56
edm::InputTag metLabel4_
Definition: QcdHighPtDQM.h:57
QcdHighPtDQM::~QcdHighPtDQM ( )
virtual

Destructor.

Definition at line 47 of file QcdHighPtDQM.cc.

47  {
48 
49 
50 }

Member Function Documentation

void QcdHighPtDQM::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
virtual

Get the analysis.

Implements edm::EDAnalyzer.

Definition at line 138 of file QcdHighPtDQM.cc.

References reco::LeafCandidate::energy(), edm::Event::getByLabel(), jetLabel_, fwrapper::jets, M_PI, MEcontainer_, CaloMET_cfi::met, CaloMET_cfi::metHO, metLabel1_, metLabel2_, metLabel3_, metLabel4_, CaloMET_cfi::metNoHF, CaloMET_cfi::metNoHFHO, movers(), reco::LeafCandidate::p4(), and reco::LeafCandidate::pt().

138  {
139 
140  //Get Jets
142  iEvent.getByLabel(jetLabel_,jetHandle);
143  const CaloJetCollection & jets = *jetHandle;
144  CaloJetCollection::const_iterator jet_iter;
145 
146  //Get MET collections
148  iEvent.getByLabel(metLabel1_, metHandle);
149  const CaloMETCollection &met = *metHandle;
150 
151  edm::Handle<CaloMETCollection> metHOHandle;
152  iEvent.getByLabel(metLabel2_, metHOHandle);
153  const CaloMETCollection &metHO = *metHOHandle;
154 
155  edm::Handle<CaloMETCollection> metNoHFHandle;
156  iEvent.getByLabel(metLabel3_, metNoHFHandle);
157  const CaloMETCollection &metNoHF = *metNoHFHandle;
158 
159  edm::Handle<CaloMETCollection> metNoHFHOHandle;
160  iEvent.getByLabel(metLabel4_, metNoHFHOHandle);
161  const CaloMETCollection &metNoHFHO = *metNoHFHOHandle;
162 
163  //initialize leading jet value and jet multiplicity counter
164  int njets = 0;
165  int njets30 = 0;
166  float leading_jetpt = 0;
167  float leading_jeteta = 0;
168 
169  //initialize variables for picking out leading 2 barrel jets
170  reco::CaloJet leadingbarreljet;
171  reco::CaloJet secondbarreljet;
172  int nbarreljets = 0;
173 
174  //get bins in eta.
175  //Bins correspond to calotower regions.
176 
177  const float etabins[83] = {-5.191, -4.889, -4.716, -4.538, -4.363, -4.191, -4.013, -3.839, -3.664, -3.489, -3.314, -3.139, -2.964, -2.853, -2.650, -2.500, -2.322, -2.172, -2.043, -1.930, -1.830, -1.740, -1.653, -1.566, -1.479, -1.392, -1.305, -1.218, -1.131, -1.044, -.957, -.879, -.783, -.696, -.609, -.522, -.435, -.348, -.261, -.174, -.087, 0, .087, .174, .261, .348, .435, .522, .609, .696, .783, .879, .957, 1.044, 1.131, 1.218, 1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830, 1.930, 2.043, 2.172, 2.322, 2.500, 2.650, 2.853, 2.964, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.889, 5.191};
178 
179  for(jet_iter = jets.begin(); jet_iter!= jets.end(); ++jet_iter){
180  njets++;
181 
182  //get Jet stats
183  float jet_pt = jet_iter->pt();
184  float jet_eta = jet_iter->eta();
185  float jet_phi = jet_iter->phi();
186 
187  //fill jet Pt and jet EMF
188  MEcontainer_["inclusive_jet_pt"]->Fill(jet_pt);
189  MEcontainer_["inclusive_jet_EMF"]->Fill(jet_iter->emEnergyFraction());
190 
191  // pick out up to the first 2 leading barrel jets
192  // for use in calculating dijet mass in barrel region
193  // also fill jet Pt histogram for barrel
194 
195  if(jet_eta <= 1.3){
196  MEcontainer_["inclusive_jet_pt_barrel"]->Fill(jet_pt);
197  if(nbarreljets == 0){
198  leadingbarreljet = jets[(njets-1)];
199  nbarreljets++;
200  }
201  else if(nbarreljets == 1){
202  secondbarreljet = jets[(njets-1)];
203  nbarreljets++;
204  }
205 
206  }
207 
208  // fill jet Pt for endcap and forward regions
209  else if(jet_eta <= 3.0 && jet_eta > 1.3){MEcontainer_["inclusive_jet_pt_endcap"]->Fill(jet_pt);}
210 
211  else if(jet_eta <= 5.0 && jet_eta > 3.0){MEcontainer_["inclusive_jet_pt_forward"]->Fill(jet_pt);}
212 
213 
214  // count jet multiplicity for jets with Pt > 30
215  if((jet_pt) > 30) njets30++;
216 
217  // check leading jet quantities
218  if(jet_pt > leading_jetpt){
219  leading_jetpt = jet_pt;
220  leading_jeteta = jet_eta;
221  }
222 
223  //fill eta-phi plot
224  for (int eit = 0; eit < 81; eit++){
225  for(int pit = 0; pit < 72; pit++){
226  float low_eta = etabins[eit];
227  float high_eta = etabins[eit+1];
228  float low_phi = (-M_PI) + pit*(M_PI/36);
229  float high_phi = low_phi + (M_PI/36);
230  if(jet_eta > low_eta && jet_eta < high_eta && jet_phi > low_phi && jet_phi < high_phi){
231  MEcontainer_["etaphi"]->Fill((eit - 41), jet_phi);}
232 
233 
234  }
235  }
236  }
237 
238  // after iterating over all jets, fill leading jet quantity histograms
239  // and jet multiplicity histograms
240 
241  MEcontainer_["leading_jet_pt"]->Fill(leading_jetpt);
242 
243  if(leading_jeteta <= 1.3){MEcontainer_["leading_jet_pt_barrel"]->Fill(leading_jetpt);}
244 
245  else if(leading_jeteta <= 3.0 && leading_jeteta > 1.3){MEcontainer_["leading_jet_pt_endcap"]->Fill(leading_jetpt);}
246 
247  else if(leading_jeteta <= 5.0 && leading_jeteta > 3.0){MEcontainer_["leading_jet_pt_forward"]->Fill(leading_jetpt);}
248 
249  MEcontainer_["njets"]->Fill(njets);
250 
251  MEcontainer_["njets30"]->Fill(njets30);
252 
253  // fill MET over Sum ET and Leading jet PT for all MET flavors
254  MEcontainer_["movers_met"]->Fill(movers(met));
255  MEcontainer_["moverl_met"]->Fill(movers(met), leading_jetpt);
256  MEcontainer_["movers_metho"]->Fill(movers(metHO));
257  MEcontainer_["moverl_metho"]->Fill(movers(metHO), leading_jetpt);
258  MEcontainer_["movers_metnohf"]->Fill(movers(metNoHF));
259  MEcontainer_["moverl_metnohf"]->Fill(movers(metNoHF), leading_jetpt);
260  MEcontainer_["movers_metnohfho"]->Fill(movers(metNoHFHO));
261  MEcontainer_["moverl_metnohfho"]->Fill(movers(metNoHFHO), leading_jetpt);
262 
263 
264  // fetch first 3 jet EMF
265 
266  if(jets.size() >= 1){
267  MEcontainer_["leading_jet_EMF"]->Fill(jets[0].emEnergyFraction());
268  if(jets.size() >=2){
269  MEcontainer_["second_jet_EMF"]->Fill(jets[1].emEnergyFraction());
270  if(jets.size() >= 3){
271  MEcontainer_["third_jet_EMF"]->Fill(jets[2].emEnergyFraction());
272  }
273  }
274  }
275 
276  // if 2 nontrivial barrel jets, reconstruct dijet mass
277 
278  if(nbarreljets == 2){
279  if(leadingbarreljet.energy() > 0 && secondbarreljet.energy() > 0){
280  math::XYZTLorentzVector DiJet = leadingbarreljet.p4() + secondbarreljet.p4();
281  float dijet_mass = DiJet.mass();
282  MEcontainer_["dijet_mass"]->Fill(dijet_mass);
283  }
284  }
285 
286 
287 }
virtual double energy() const GCC11_FINAL
energy
edm::InputTag metLabel2_
Definition: QcdHighPtDQM.h:55
Jets made from CaloTowers.
Definition: CaloJet.h:30
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
edm::InputTag jetLabel_
Definition: QcdHighPtDQM.h:53
edm::InputTag metLabel1_
Definition: QcdHighPtDQM.h:54
tuple metNoHF
Definition: CaloMET_cfi.py:50
std::map< std::string, MonitorElement * > MEcontainer_
Definition: QcdHighPtDQM.h:60
Collection of Calo MET.
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
tuple metNoHFHO
Definition: CaloMET_cfi.py:62
vector< PseudoJet > jets
edm::InputTag metLabel3_
Definition: QcdHighPtDQM.h:56
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
float movers(const reco::CaloMETCollection &metcollection)
#define M_PI
Definition: BFit3D.cc:3
edm::InputTag metLabel4_
Definition: QcdHighPtDQM.h:57
tuple metHO
Definition: CaloMET_cfi.py:30
virtual float pt() const GCC11_FINAL
transverse momentum
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
void QcdHighPtDQM::beginJob ( void  )
virtual

Inizialize parameters for histo binning.

Reimplemented from edm::EDAnalyzer.

Definition at line 53 of file QcdHighPtDQM.cc.

References DQMStore::book1D(), DQMStore::book2D(), M_PI, MEcontainer_, DQMStore::setCurrentFolder(), and theDbe.

53  {
54 
55 
56  //Book MEs
57 
58  theDbe->setCurrentFolder("Physics/QcdHighPt");
59 
60  MEcontainer_["dijet_mass"] = theDbe->book1D("dijet_mass", "dijet resonance invariant mass, barrel region", 100, 0, 1000);
61  MEcontainer_["njets"] = theDbe->book1D("njets", "jet multiplicity", 10, 0, 10);
62  MEcontainer_["etaphi"] = theDbe->book2D("etaphi", "eta/phi distribution", 83, -42, 42, 72, -M_PI, M_PI);
63  MEcontainer_["njets30"] = theDbe->book1D("njets30", "jet multiplicity, pt > 30 GeV", 10, 0, 10);
64 
65  //book histograms for inclusive jet quantities
66  MEcontainer_["inclusive_jet_pt"] = theDbe->book1D("inclusive_jet_pt", "inclusive jet Pt spectrum", 100, 0, 1000);
67  MEcontainer_["inclusive_jet_pt_barrel"] = theDbe->book1D("inclusive_jet_pt_barrel", "inclusive jet Pt, eta < 1.3", 100, 0, 1000);
68  MEcontainer_["inclusive_jet_pt_forward"] = theDbe->book1D("inclusive_jet_pt_forward", "inclusive jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
69  MEcontainer_["inclusive_jet_pt_endcap"] = theDbe->book1D("inclusive_jet_pt_endcap", "inclusive jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
70 
71  //book histograms for leading jet quantities
72  MEcontainer_["leading_jet_pt"] = theDbe->book1D("leading_jet_pt", "leading jet Pt", 100, 0, 1000);
73  MEcontainer_["leading_jet_pt_barrel"] = theDbe->book1D("leading_jet_pt_barrel", "leading jet Pt, eta < 1.3", 100, 0, 1000);
74  MEcontainer_["leading_jet_pt_forward"] = theDbe->book1D("leading_jet_pt_forward", "leading jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
75  MEcontainer_["leading_jet_pt_endcap"] = theDbe->book1D("leading_jet_pt_endcap", "leading jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
76 
77  //book histograms for met over sum et and met over leading jet pt for various
78  //flavors of MET
79  MEcontainer_["movers_met"] = theDbe->book1D("movers_met", "MET over Sum ET for basic MET collection", 50, 0, 1);
80  MEcontainer_["moverl_met"] = theDbe->book1D("moverl_met", "MET over leading jet Pt for basic MET collection", 50, 0, 2);
81 
82  MEcontainer_["movers_metho"] = theDbe->book1D("movers_metho", "MET over Sum ET for MET HO collection", 50, 0, 1);
83  MEcontainer_["moverl_metho"] = theDbe->book1D("moverl_metho", "MET over leading jet Pt for MET HO collection", 50, 0, 2);
84 
85  MEcontainer_["movers_metnohf"] = theDbe->book1D("movers_metnohf", "MET over Sum ET for MET no HF collection", 50, 0, 1);
86  MEcontainer_["moverl_metnohf"] = theDbe->book1D("moverl_metnohf", "MET over leading jet Pt for MET no HF collection", 50, 0, 2);
87 
88  MEcontainer_["movers_metnohfho"] = theDbe->book1D("movers_metnohfho", "MET over Sum ET for MET no HF HO collection", 50, 0, 1);
89  MEcontainer_["moverl_metnohfho"] = theDbe->book1D("moverl_metnohfho", "MET over leading jet Pt for MET no HF HO collection", 50, 0, 2);
90 
91  //book histograms for EMF fraction for all jets and first 3 jets
92  MEcontainer_["inclusive_jet_EMF"] = theDbe->book1D("inclusive_jet_EMF", "inclusive jet EMF", 50, -1, 1);
93  MEcontainer_["leading_jet_EMF"] = theDbe->book1D("leading_jet_EMF", "leading jet EMF", 50, -1, 1);
94  MEcontainer_["second_jet_EMF"] = theDbe->book1D("second_jet_EMF", "second jet EMF", 50, -1, 1);
95  MEcontainer_["third_jet_EMF"] = theDbe->book1D("third_jet_EMF", "third jet EMF", 50, -1, 1);
96 
97 
98 
99 
100 }
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:722
DQMStore * theDbe
Definition: QcdHighPtDQM.h:50
std::map< std::string, MonitorElement * > MEcontainer_
Definition: QcdHighPtDQM.h:60
#define M_PI
Definition: BFit3D.cc:3
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:850
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434
void QcdHighPtDQM::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 102 of file QcdHighPtDQM.cc.

102  {
103 
104 }
float QcdHighPtDQM::moverl ( const reco::CaloMETCollection metcollection,
float &  ljpt 
)
private

Definition at line 124 of file QcdHighPtDQM.cc.

References CaloMET_cfi::met, and mathSSE::sqrt().

124  {
125  float metoverl = 0;
126  CaloMETCollection::const_iterator met_iter;
127  for(met_iter = metcollection.begin(); met_iter!= metcollection.end(); ++met_iter)
128  {
129  float mex = met_iter->momentum().x();
130  float mey = met_iter->momentum().y();
131  float met = sqrt(mex*mex+mey*mey);
132  metoverl = (met / ljpt);
133  }
134  return metoverl;
135 }
T sqrt(T t)
Definition: SSEVec.h:48
float QcdHighPtDQM::movers ( const reco::CaloMETCollection metcollection)
private

Definition at line 107 of file QcdHighPtDQM.cc.

References CaloMET_cfi::met, and mathSSE::sqrt().

Referenced by analyze().

107  {
108 
109  float metovers = 0;
110  CaloMETCollection::const_iterator met_iter;
111  for(met_iter = metcollection.begin(); met_iter!= metcollection.end(); ++met_iter)
112  {
113  float mex = met_iter->momentum().x();
114  float mey = met_iter->momentum().y();
115  float met = sqrt(mex*mex+mey*mey);
116  float sumet = met_iter->sumEt();
117  metovers = (met / sumet);
118  }
119  return metovers;
120 }
T sqrt(T t)
Definition: SSEVec.h:48

Member Data Documentation

edm::InputTag QcdHighPtDQM::jetLabel_
private

Definition at line 53 of file QcdHighPtDQM.h.

Referenced by analyze().

std::map<std::string, MonitorElement*> QcdHighPtDQM::MEcontainer_
private

Definition at line 60 of file QcdHighPtDQM.h.

Referenced by analyze(), and beginJob().

edm::InputTag QcdHighPtDQM::metLabel1_
private

Definition at line 54 of file QcdHighPtDQM.h.

Referenced by analyze().

edm::InputTag QcdHighPtDQM::metLabel2_
private

Definition at line 55 of file QcdHighPtDQM.h.

Referenced by analyze().

edm::InputTag QcdHighPtDQM::metLabel3_
private

Definition at line 56 of file QcdHighPtDQM.h.

Referenced by analyze().

edm::InputTag QcdHighPtDQM::metLabel4_
private

Definition at line 57 of file QcdHighPtDQM.h.

Referenced by analyze().

DQMStore* QcdHighPtDQM::theDbe
private

Definition at line 50 of file QcdHighPtDQM.h.

Referenced by beginJob(), and QcdHighPtDQM().