CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
QcdHighPtDQM.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author S. Bolognesi, Eric - CERN
5  */
6 
8 
13 
17 
20 #include <vector>
21 
22 #include <string>
23 #include <cmath>
24 
25 using namespace std;
26 using namespace edm;
27 using namespace reco;
28 using namespace math;
29 
30 
31 //Get Jets and MET (no MET plots yet pending converging w/JetMET group)
32 
34  jetToken_(consumes<CaloJetCollection>
35  (iConfig.getUntrackedParameter<edm::InputTag>("jetTag"))),
36  metToken1_(consumes<CaloMETCollection>
37  (iConfig.getUntrackedParameter<edm::InputTag>("metTag1"))),
38  metToken2_(consumes<CaloMETCollection>
39  (iConfig.getUntrackedParameter<edm::InputTag>("metTag2"))),
40  metToken3_(consumes<CaloMETCollection>
41  (iConfig.getUntrackedParameter<edm::InputTag>("metTag3"))),
42  metToken4_(consumes<CaloMETCollection>
43  (iConfig.getUntrackedParameter<edm::InputTag>("metTag4")))
44 {
45 
47 
48 }
49 
51 
52 
53 }
54 
55 
57 
58 
59  //Book MEs
60 
61  theDbe->setCurrentFolder("Physics/QcdHighPt");
62 
63  MEcontainer_["dijet_mass"] = theDbe->book1D("dijet_mass", "dijet resonance invariant mass, barrel region", 100, 0, 1000);
64  MEcontainer_["njets"] = theDbe->book1D("njets", "jet multiplicity", 10, 0, 10);
65  MEcontainer_["etaphi"] = theDbe->book2D("etaphi", "eta/phi distribution", 83, -42, 42, 72, -M_PI, M_PI);
66  MEcontainer_["njets30"] = theDbe->book1D("njets30", "jet multiplicity, pt > 30 GeV", 10, 0, 10);
67 
68  //book histograms for inclusive jet quantities
69  MEcontainer_["inclusive_jet_pt"] = theDbe->book1D("inclusive_jet_pt", "inclusive jet Pt spectrum", 100, 0, 1000);
70  MEcontainer_["inclusive_jet_pt_barrel"] = theDbe->book1D("inclusive_jet_pt_barrel", "inclusive jet Pt, eta < 1.3", 100, 0, 1000);
71  MEcontainer_["inclusive_jet_pt_forward"] = theDbe->book1D("inclusive_jet_pt_forward", "inclusive jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
72  MEcontainer_["inclusive_jet_pt_endcap"] = theDbe->book1D("inclusive_jet_pt_endcap", "inclusive jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
73 
74  //book histograms for leading jet quantities
75  MEcontainer_["leading_jet_pt"] = theDbe->book1D("leading_jet_pt", "leading jet Pt", 100, 0, 1000);
76  MEcontainer_["leading_jet_pt_barrel"] = theDbe->book1D("leading_jet_pt_barrel", "leading jet Pt, eta < 1.3", 100, 0, 1000);
77  MEcontainer_["leading_jet_pt_forward"] = theDbe->book1D("leading_jet_pt_forward", "leading jet Pt, 3.0 < eta < 5.0", 100, 0, 1000);
78  MEcontainer_["leading_jet_pt_endcap"] = theDbe->book1D("leading_jet_pt_endcap", "leading jet Pt, 1.3 < eta < 3.0", 100, 0, 1000);
79 
80  //book histograms for met over sum et and met over leading jet pt for various
81  //flavors of MET
82  MEcontainer_["movers_met"] = theDbe->book1D("movers_met", "MET over Sum ET for basic MET collection", 50, 0, 1);
83  MEcontainer_["moverl_met"] = theDbe->book1D("moverl_met", "MET over leading jet Pt for basic MET collection", 50, 0, 2);
84 
85  MEcontainer_["movers_metho"] = theDbe->book1D("movers_metho", "MET over Sum ET for MET HO collection", 50, 0, 1);
86  MEcontainer_["moverl_metho"] = theDbe->book1D("moverl_metho", "MET over leading jet Pt for MET HO collection", 50, 0, 2);
87 
88  MEcontainer_["movers_metnohf"] = theDbe->book1D("movers_metnohf", "MET over Sum ET for MET no HF collection", 50, 0, 1);
89  MEcontainer_["moverl_metnohf"] = theDbe->book1D("moverl_metnohf", "MET over leading jet Pt for MET no HF collection", 50, 0, 2);
90 
91  MEcontainer_["movers_metnohfho"] = theDbe->book1D("movers_metnohfho", "MET over Sum ET for MET no HF HO collection", 50, 0, 1);
92  MEcontainer_["moverl_metnohfho"] = theDbe->book1D("moverl_metnohfho", "MET over leading jet Pt for MET no HF HO collection", 50, 0, 2);
93 
94  //book histograms for EMF fraction for all jets and first 3 jets
95  MEcontainer_["inclusive_jet_EMF"] = theDbe->book1D("inclusive_jet_EMF", "inclusive jet EMF", 50, -1, 1);
96  MEcontainer_["leading_jet_EMF"] = theDbe->book1D("leading_jet_EMF", "leading jet EMF", 50, -1, 1);
97  MEcontainer_["second_jet_EMF"] = theDbe->book1D("second_jet_EMF", "second jet EMF", 50, -1, 1);
98  MEcontainer_["third_jet_EMF"] = theDbe->book1D("third_jet_EMF", "third jet EMF", 50, -1, 1);
99 
100 
101 
102 
103 }
104 
106 
107 }
108 
109 // method to calculate MET over Sum ET from a particular MET collection
110 float QcdHighPtDQM::movers(const CaloMETCollection &metcollection) {
111 
112  float metovers = 0;
113  CaloMETCollection::const_iterator met_iter;
114  for(met_iter = metcollection.begin(); met_iter!= metcollection.end(); ++met_iter)
115  {
116  float mex = met_iter->momentum().x();
117  float mey = met_iter->momentum().y();
118  float met = sqrt(mex*mex+mey*mey);
119  float sumet = met_iter->sumEt();
120  metovers = (met / sumet);
121  }
122  return metovers;
123 }
124 
125 
126 // method to calculate MET over Leading jet PT for a particular MET collection
127 float QcdHighPtDQM::moverl(const CaloMETCollection &metcollection, float &ljpt) {
128  float metoverl = 0;
129  CaloMETCollection::const_iterator met_iter;
130  for(met_iter = metcollection.begin(); met_iter!= metcollection.end(); ++met_iter)
131  {
132  float mex = met_iter->momentum().x();
133  float mey = met_iter->momentum().y();
134  float met = sqrt(mex*mex+mey*mey);
135  metoverl = (met / ljpt);
136  }
137  return metoverl;
138 }
139 
140 
141 void QcdHighPtDQM::analyze(const Event& iEvent, const EventSetup& iSetup) {
142 
143  //Get Jets
145  iEvent.getByToken(jetToken_,jetHandle);
146  const CaloJetCollection & jets = *jetHandle;
147  CaloJetCollection::const_iterator jet_iter;
148 
149  //Get MET collections
151  iEvent.getByToken(metToken1_, metHandle);
152  const CaloMETCollection &met = *metHandle;
153 
154  edm::Handle<CaloMETCollection> metHOHandle;
155  iEvent.getByToken(metToken2_, metHOHandle);
156  const CaloMETCollection &metHO = *metHOHandle;
157 
158  edm::Handle<CaloMETCollection> metNoHFHandle;
159  iEvent.getByToken(metToken3_, metNoHFHandle);
160  const CaloMETCollection &metNoHF = *metNoHFHandle;
161 
162  edm::Handle<CaloMETCollection> metNoHFHOHandle;
163  iEvent.getByToken(metToken4_, metNoHFHOHandle);
164  const CaloMETCollection &metNoHFHO = *metNoHFHOHandle;
165 
166  //initialize leading jet value and jet multiplicity counter
167  int njets = 0;
168  int njets30 = 0;
169  float leading_jetpt = 0;
170  float leading_jeteta = 0;
171 
172  //initialize variables for picking out leading 2 barrel jets
173  reco::CaloJet leadingbarreljet;
174  reco::CaloJet secondbarreljet;
175  int nbarreljets = 0;
176 
177  //get bins in eta.
178  //Bins correspond to calotower regions.
179 
180  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};
181 
182  for(jet_iter = jets.begin(); jet_iter!= jets.end(); ++jet_iter){
183  njets++;
184 
185  //get Jet stats
186  float jet_pt = jet_iter->pt();
187  float jet_eta = jet_iter->eta();
188  float jet_phi = jet_iter->phi();
189 
190  //fill jet Pt and jet EMF
191  MEcontainer_["inclusive_jet_pt"]->Fill(jet_pt);
192  MEcontainer_["inclusive_jet_EMF"]->Fill(jet_iter->emEnergyFraction());
193 
194  // pick out up to the first 2 leading barrel jets
195  // for use in calculating dijet mass in barrel region
196  // also fill jet Pt histogram for barrel
197 
198  if(jet_eta <= 1.3){
199  MEcontainer_["inclusive_jet_pt_barrel"]->Fill(jet_pt);
200  if(nbarreljets == 0){
201  leadingbarreljet = jets[(njets-1)];
202  nbarreljets++;
203  }
204  else if(nbarreljets == 1){
205  secondbarreljet = jets[(njets-1)];
206  nbarreljets++;
207  }
208 
209  }
210 
211  // fill jet Pt for endcap and forward regions
212  else if(jet_eta <= 3.0 && jet_eta > 1.3){MEcontainer_["inclusive_jet_pt_endcap"]->Fill(jet_pt);}
213 
214  else if(jet_eta <= 5.0 && jet_eta > 3.0){MEcontainer_["inclusive_jet_pt_forward"]->Fill(jet_pt);}
215 
216 
217  // count jet multiplicity for jets with Pt > 30
218  if((jet_pt) > 30) njets30++;
219 
220  // check leading jet quantities
221  if(jet_pt > leading_jetpt){
222  leading_jetpt = jet_pt;
223  leading_jeteta = jet_eta;
224  }
225 
226  //fill eta-phi plot
227  for (int eit = 0; eit < 81; eit++){
228  for(int pit = 0; pit < 72; pit++){
229  float low_eta = etabins[eit];
230  float high_eta = etabins[eit+1];
231  float low_phi = (-M_PI) + pit*(M_PI/36);
232  float high_phi = low_phi + (M_PI/36);
233  if(jet_eta > low_eta && jet_eta < high_eta && jet_phi > low_phi && jet_phi < high_phi){
234  MEcontainer_["etaphi"]->Fill((eit - 41), jet_phi);}
235 
236 
237  }
238  }
239  }
240 
241  // after iterating over all jets, fill leading jet quantity histograms
242  // and jet multiplicity histograms
243 
244  MEcontainer_["leading_jet_pt"]->Fill(leading_jetpt);
245 
246  if(leading_jeteta <= 1.3){MEcontainer_["leading_jet_pt_barrel"]->Fill(leading_jetpt);}
247 
248  else if(leading_jeteta <= 3.0 && leading_jeteta > 1.3){MEcontainer_["leading_jet_pt_endcap"]->Fill(leading_jetpt);}
249 
250  else if(leading_jeteta <= 5.0 && leading_jeteta > 3.0){MEcontainer_["leading_jet_pt_forward"]->Fill(leading_jetpt);}
251 
252  MEcontainer_["njets"]->Fill(njets);
253 
254  MEcontainer_["njets30"]->Fill(njets30);
255 
256  // fill MET over Sum ET and Leading jet PT for all MET flavors
257  MEcontainer_["movers_met"]->Fill(movers(met));
258  MEcontainer_["moverl_met"]->Fill(movers(met), leading_jetpt);
259  MEcontainer_["movers_metho"]->Fill(movers(metHO));
260  MEcontainer_["moverl_metho"]->Fill(movers(metHO), leading_jetpt);
261  MEcontainer_["movers_metnohf"]->Fill(movers(metNoHF));
262  MEcontainer_["moverl_metnohf"]->Fill(movers(metNoHF), leading_jetpt);
263  MEcontainer_["movers_metnohfho"]->Fill(movers(metNoHFHO));
264  MEcontainer_["moverl_metnohfho"]->Fill(movers(metNoHFHO), leading_jetpt);
265 
266 
267  // fetch first 3 jet EMF
268 
269  if(jets.size() >= 1){
270  MEcontainer_["leading_jet_EMF"]->Fill(jets[0].emEnergyFraction());
271  if(jets.size() >=2){
272  MEcontainer_["second_jet_EMF"]->Fill(jets[1].emEnergyFraction());
273  if(jets.size() >= 3){
274  MEcontainer_["third_jet_EMF"]->Fill(jets[2].emEnergyFraction());
275  }
276  }
277  }
278 
279  // if 2 nontrivial barrel jets, reconstruct dijet mass
280 
281  if(nbarreljets == 2){
282  if(leadingbarreljet.energy() > 0 && secondbarreljet.energy() > 0){
283  math::XYZTLorentzVector DiJet = leadingbarreljet.p4() + secondbarreljet.p4();
284  float dijet_mass = DiJet.mass();
285  MEcontainer_["dijet_mass"]->Fill(dijet_mass);
286  }
287  }
288 
289 
290 }
291 
292 
293 
294 
295 
296 
297 
298 
299 
300 
301 
302 
303 
304 
305 // Local Variables:
306 // show-trailing-whitespace: t
307 // truncate-lines: t
308 // End:
edm::EDGetTokenT< reco::CaloMETCollection > metToken1_
Definition: QcdHighPtDQM.h:54
virtual double energy() const GCC11_FINAL
energy
float moverl(const reco::CaloMETCollection &metcollection, float &ljpt)
void analyze(const edm::Event &, const edm::EventSetup &)
Get the analysis.
tuple met
____________________________________________________________________________||
Definition: CaloMET_cfi.py:4
Jets made from CaloTowers.
Definition: CaloJet.h:29
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:928
virtual const LorentzVector & p4() const GCC11_FINAL
four-momentum Lorentz vector
DQMStore * theDbe
Definition: QcdHighPtDQM.h:50
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
edm::EDGetTokenT< reco::CaloMETCollection > metToken2_
Definition: QcdHighPtDQM.h:55
tuple metNoHF
____________________________________________________________________________||
Definition: CaloMET_cfi.py:34
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:29
int iEvent
Definition: GenABIO.cc:230
T sqrt(T t)
Definition: SSEVec.h:48
tuple metNoHFHO
____________________________________________________________________________||
Definition: CaloMET_cfi.py:44
vector< PseudoJet > jets
void beginJob()
Inizialize parameters for histo binning.
Definition: QcdHighPtDQM.cc:56
edm::EDGetTokenT< reco::CaloMETCollection > metToken4_
Definition: QcdHighPtDQM.h:57
edm::EDGetTokenT< reco::CaloMETCollection > metToken3_
Definition: QcdHighPtDQM.h:56
#define M_PI
float movers(const reco::CaloMETCollection &metcollection)
void endJob(void)
QcdHighPtDQM(const edm::ParameterSet &)
Constructor.
Definition: QcdHighPtDQM.cc:33
tuple metHO
____________________________________________________________________________||
Definition: CaloMET_cfi.py:14
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:1056
virtual float pt() const GCC11_FINAL
transverse momentum
edm::EDGetTokenT< reco::CaloJetCollection > jetToken_
Definition: QcdHighPtDQM.h:53
virtual ~QcdHighPtDQM()
Destructor.
Definition: QcdHighPtDQM.cc:50
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:639
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects