ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/VHbbAnalysis/VHbbDataFormats/src/CutsAndHistos.cc
Revision: 1.5
Committed: Mon Jul 25 08:55:41 2011 UTC (13 years, 9 months ago) by tboccali
Content type: text/plain
Branch: MAIN
CVS Tags: Jul25th2011
Changes since 1.4: +48 -48 lines
Log Message:
as asked by andrea fourMomentum->p4

File Contents

# Content
1 #include "VHbbAnalysis/VHbbDataFormats/interface/VHbbEvent.h"
2 #include "VHbbAnalysis/VHbbDataFormats/interface/VHbbCandidate.h"
3 #include "VHbbAnalysis/VHbbDataFormats/interface/VHbbProxy.h"
4 #include "VHbbAnalysis/VHbbDataFormats/interface/CutsAndHistos.h"
5 #include <TH1F.h>
6 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
7 #include <sstream>
8 #include "TKey.h"
9
10 //class LeptonSelection : public Cut
11
12 class mcHPtCut : public Cut {
13 public:
14 mcHPtCut(Double_t mcHPtCutMin_):
15 mcHPtCutMin(mcHPtCutMin_) {}
16 mcHPtCut(Double_t mcHPtCutMin_, Double_t mcHPtCutMax_):
17 mcHPtCutMin(mcHPtCutMin_),
18 mcHPtCutMax(mcHPtCutMax_) {}
19 std::string name() {
20 oss_mcHPtCutMin << mcHPtCutMin;
21 if(!mcHPtCutMax)
22 return "mcHpT"+oss_mcHPtCutMin.str();
23 else{
24 oss_mcHPtCutMax << mcHPtCutMax;
25 return "mcHpT"+oss_mcHPtCutMin.str()+"To"+oss_mcHPtCutMax.str();
26 }
27 }
28 bool pass(VHbbProxy &iProxy) {
29 if(!mcHPtCutMax)
30 return ( iProxy.getVHbbEventAuxInfo()->mcH.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcH)[0]).p4.Pt() > mcHPtCutMin );
31 else
32 return ( iProxy.getVHbbEventAuxInfo()->mcH.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcH)[0]).p4.Pt() > mcHPtCutMin
33 && iProxy.getVHbbEventAuxInfo()->mcH.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcH)[0]).p4.Pt() < mcHPtCutMax );
34 }
35 private:
36 Double_t mcHPtCutMin;
37 Double_t mcHPtCutMax;
38 std::ostringstream oss_mcHPtCutMin;
39 std::ostringstream oss_mcHPtCutMax;
40 };
41
42
43 class mcWPtCut : public Cut {
44 public:
45 mcWPtCut(Double_t mcWPtCutMin_):
46 mcWPtCutMin(mcWPtCutMin_) {}
47 mcWPtCut(Double_t mcWPtCutMin_, Double_t mcWPtCutMax_):
48 mcWPtCutMin(mcWPtCutMin_),
49 mcWPtCutMax(mcWPtCutMax_) {}
50 std::string name() {
51 oss_mcWPtCutMin << mcWPtCutMin;
52 if(!mcWPtCutMax)
53 return "mcWpT"+oss_mcWPtCutMin.str();
54 else{
55 oss_mcWPtCutMax << mcWPtCutMax;
56 return "mcWpT"+oss_mcWPtCutMin.str()+"To"+oss_mcWPtCutMax.str();
57 }
58 }
59 bool pass(VHbbProxy &iProxy) {
60 if(!mcWPtCutMax)
61 return ( iProxy.getVHbbEventAuxInfo()->mcW.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcW)[0]).p4.Pt() > mcWPtCutMin );
62 else
63 return ( iProxy.getVHbbEventAuxInfo()->mcW.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcW)[0]).p4.Pt() > mcWPtCutMin
64 && iProxy.getVHbbEventAuxInfo()->mcW.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcW)[0]).p4.Pt() < mcWPtCutMax );
65 }
66 private:
67 Double_t mcWPtCutMin;
68 Double_t mcWPtCutMax;
69 std::ostringstream oss_mcWPtCutMin;
70 std::ostringstream oss_mcWPtCutMax;
71 };
72
73
74
75 class mcZPtCut : public Cut {
76 public:
77 mcZPtCut(Double_t mcZPtCutMin_):
78 mcZPtCutMin(mcZPtCutMin_) {}
79 mcZPtCut(Double_t mcZPtCutMin_, Double_t mcZPtCutMax_):
80 mcZPtCutMin(mcZPtCutMin_),
81 mcZPtCutMax(mcZPtCutMax_) {}
82 std::string name() {
83 oss_mcZPtCutMin << mcZPtCutMin;
84 if(!mcZPtCutMax)
85 return "mcZpT"+oss_mcZPtCutMin.str();
86 else{
87 oss_mcZPtCutMax << mcZPtCutMax;
88 return "mcZpT"+oss_mcZPtCutMin.str()+"To"+oss_mcZPtCutMax.str();
89 }
90 }
91 bool pass(VHbbProxy &iProxy) {
92 if(!mcZPtCutMax)
93 return ( iProxy.getVHbbEventAuxInfo()->mcZ.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcZ)[0]).p4.Pt() > mcZPtCutMin );
94 else
95 return ( iProxy.getVHbbEventAuxInfo()->mcZ.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcZ)[0]).p4.Pt() > mcZPtCutMin
96 && iProxy.getVHbbEventAuxInfo()->mcZ.size() !=0 && ((iProxy.getVHbbEventAuxInfo()->mcZ)[0]).p4.Pt() < mcZPtCutMax );
97 }
98 private:
99 Double_t mcZPtCutMin;
100 Double_t mcZPtCutMax;
101 std::ostringstream oss_mcZPtCutMin;
102 std::ostringstream oss_mcZPtCutMax;
103 };
104
105
106 class SimpleJet1PtCut: public Cut {
107 public:
108 SimpleJet1PtCut(Double_t simpleJet1PtCutMin_):
109 simpleJet1PtCutMin(simpleJet1PtCutMin_) {}
110 SimpleJet1PtCut(Double_t simpleJet1PtCutMin_, Double_t simpleJet1PtCutMax_):
111 simpleJet1PtCutMin(simpleJet1PtCutMin_),
112 simpleJet1PtCutMax(simpleJet1PtCutMax_) {}
113 std::string name() {
114 oss_simpleJet1PtCutMin << simpleJet1PtCutMin;
115 if(!simpleJet1PtCutMax){
116 return "SimpleJet1PtCut"+oss_simpleJet1PtCutMin.str();
117 }
118 else{
119 oss_simpleJet1PtCutMax << simpleJet1PtCutMax;
120 return "SimpleJet1PtCut"+oss_simpleJet1PtCutMin.str()+"To"+oss_simpleJet1PtCutMax.str();
121 }
122 }
123 bool pass(VHbbProxy &iProxy) {
124 if(iProxy.getVHbbCandidate()->size() < 1)
125 return false;
126 else
127 if(!simpleJet1PtCutMin)
128 return ((iProxy.getVHbbCandidate()->at(0)).H.jets.at(0).p4.Pt() > simpleJet1PtCutMin);
129 else
130 return ((iProxy.getVHbbCandidate()->at(0)).H.jets.at(0).p4.Pt() > simpleJet1PtCutMin
131 && (iProxy.getVHbbCandidate()->at(0)).H.jets.at(0).p4.Pt() < simpleJet1PtCutMax);
132 }
133 private:
134 Double_t simpleJet1PtCutMin;
135 Double_t simpleJet1PtCutMax;
136 std::ostringstream oss_simpleJet1PtCutMin;
137 std::ostringstream oss_simpleJet1PtCutMax;
138 };
139
140 class VPtCut: public Cut {
141 public:
142 VPtCut(Double_t VptCutMin_):
143 VptCutMin(VptCutMin_) {}
144 VPtCut(Double_t VptCutMin_, Double_t VptCutMax_):
145 VptCutMin(VptCutMin_),
146 VptCutMax(VptCutMax_) {}
147 std::string name() {
148 oss_VptMin << VptCutMin;
149 if(! VptCutMax ){
150 return "VpT"+oss_VptMin.str();
151 }
152 else{
153 oss_VptMax << VptCutMax;
154 return "VpT"+oss_VptMin.str()+"To"+oss_VptMax.str();
155 }
156 }
157 bool pass(VHbbProxy &iProxy) {
158 if(iProxy.getVHbbCandidate()->size() < 1)
159 return false;
160 else
161 if(VptCutMax != 1e10)
162 return ((iProxy.getVHbbCandidate()->at(0)).V.p4.Pt() > VptCutMin
163 && (iProxy.getVHbbCandidate()->at(0)).V.p4.Pt() < VptCutMax);
164 else
165 return ((iProxy.getVHbbCandidate()->at(0)).V.p4.Pt() > VptCutMin );
166
167 }
168 private:
169 Double_t VptCutMin;
170 Double_t VptCutMax;
171 std::ostringstream oss_VptMin;
172 std::ostringstream oss_VptMax;
173 };
174
175
176 class HPtCut: public Cut {
177 public:
178 HPtCut(Double_t hPtCutMin_):
179 hPtCutMin(hPtCutMin_) {}
180 HPtCut(Double_t hPtCutMin_, Double_t hPtCutMax_):
181 hPtCutMin(hPtCutMin_),
182 hPtCutMax(hPtCutMax_) {}
183
184 std::string name() {
185 oss_HPtCutMin << hPtCutMin;
186 if(!hPtCutMax){
187 return "HpT"+oss_HPtCutMin.str();
188 }
189 else{
190 oss_HPtCutMax << hPtCutMax;
191 return "HpT"+oss_HPtCutMin.str()+"To"+oss_HPtCutMax.str();
192 }
193 }
194 bool pass(VHbbProxy &iProxy) {
195 if(iProxy.getVHbbCandidate()->size() < 1)
196 return false;
197 else
198 if(!hPtCutMax)
199 return ((iProxy.getVHbbCandidate()->at(0)).H.p4.Pt() > hPtCutMin);
200 else
201 return ((iProxy.getVHbbCandidate()->at(0)).H.p4.Pt() > hPtCutMin
202 && (iProxy.getVHbbCandidate()->at(0)).H.p4.Pt() < hPtCutMax);
203 }
204 private:
205 Double_t hPtCutMin;
206 Double_t hPtCutMax;
207 std::ostringstream oss_HPtCutMin;
208 std::ostringstream oss_HPtCutMax;
209 };
210
211
212 class SignalRegion: public Cut {
213 std::string name() {return "SignalRegion";};
214
215 Bool_t pass(VHbbProxy &iProxy){
216
217 const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
218
219 if(iCand->size() < 1)
220 return false;
221
222 VHbbCandidate::VectorCandidate V = iProxy.getVHbbCandidate()->at(0).V;
223 VHbbCandidate::HiggsCandidate H = iProxy.getVHbbCandidate()->at(0).H;
224
225 if(iCand->at(0).candidateType == VHbbCandidate::Zmumu || iCand->at(0).candidateType == VHbbCandidate::Zee){
226 btag_csv_min = 0.5;
227 btag_csv_max = 0.9;
228 Higgs_pt = 150;
229 V_pt = 150;
230 VH_deltaPhi = 2.70;
231 nOfAdditionalJets = 2;
232 pullAngle = 1.57;
233 helicityAngle = 0.8;
234 }
235 else if(iCand->at(0).candidateType == VHbbCandidate::Wmun || iCand->at(0).candidateType == VHbbCandidate::Wen){
236 btag_csv_min = 0.5;
237 btag_csv_max = 0.9;
238 Higgs_pt = 150;
239 V_pt = 150;
240 VH_deltaPhi = 2.95;
241 nOfAdditionalJets = 2;
242 pullAngle = 1.57;
243 helicityAngle = 0.8;
244 }
245 else if(iCand->at(0).candidateType == VHbbCandidate::Znn){
246 btag_csv_min = 0.5;
247 btag_csv_max = 0.9;
248 Higgs_pt = 150;
249 V_pt = 150;
250 VH_deltaPhi = 2.95;
251 nOfAdditionalJets = 2;
252 pullAngle = 1.57;
253 helicityAngle = 0.8;
254 }
255 else
256 std::cerr << "No vector boson reconstructed. No histos will be filled." << std::endl;
257
258 Bool_t go = false;
259 if( H.p4.Pt() > Higgs_pt
260 && V.p4.Pt() > V_pt
261 && TMath::Abs( Geom::deltaPhi(H.p4.Phi(), V.p4.Phi()) ) > VH_deltaPhi
262 && ( H.jets.at(0).csv > btag_csv_min && H.jets.at(1).csv > btag_csv_min )
263 && ( H.jets.at(0).csv > btag_csv_max || H.jets.at(1).csv > btag_csv_max )
264 && iCand->at(0).additionalJets.size() < nOfAdditionalJets
265 // && V.additionalLeptons.size() < nOfAddionalLeptons
266 && TMath::Abs(H.deltaTheta) < pullAngle
267 // && TMath::Abs(H.helicityAngle) < helicityAngle
268 )
269 go = true;
270 return go;
271
272 }
273
274 private:
275
276 Double_t btag_csv_min;
277 Double_t btag_csv_max;
278 Double_t Higgs_pt;
279 Double_t V_pt;
280 Double_t VH_deltaPhi;
281 unsigned int nOfAdditionalJets;
282 unsigned int nOfAdditionalLeptons;
283 Double_t pullAngle;
284 Double_t helicityAngle;
285
286 };
287
288
289
290
291 class MCHistos : public Histos {
292
293 public:
294
295 virtual void book(TFile &f, std::string suffix) {
296
297 TDirectory *subDir;
298
299 if( ! f.GetDirectory(suffix.c_str()) )
300 subDir = f.mkdir(suffix.c_str());
301 else
302 subDir = f.GetDirectory(suffix.c_str());
303
304 subDir->cd();
305
306 bin_mass = 500;
307 min_mass = 0;
308 max_mass = 300;
309
310 bin_pt = 500;
311 min_pt = 0;
312 max_pt = 500;
313
314 bin_hel = 50;
315 min_hel = 0;
316 max_hel = 1;
317
318 //from MC
319 McH_simHMass = new TH1F(("simHiggsMass"+suffix).c_str(),("Sim Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
320 McH_simHPt = new TH1F(("simHiggsPt"+suffix).c_str(),("Sim Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
321 McH_simZMass = new TH1F(("simZMass"+suffix).c_str(),("Sim Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
322 McH_simZPt = new TH1F(("simZPt"+suffix).c_str(),("Sim Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
323 McH_simWMass = new TH1F(("simWMass"+suffix).c_str(),("Sim W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
324 McH_simWPt = new TH1F(("simWPt"+suffix).c_str(),("Sim W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
325
326 }
327
328 virtual void fill(VHbbProxy &iProxy,float w) {
329
330 // const VHbbEvent *iEvent = iProxy.getVHbbEvent();
331 const VHbbEventAuxInfo *iAuxInfo = iProxy.getVHbbEventAuxInfo();
332 if(iAuxInfo) {
333 //from MC
334 if (iAuxInfo->mcH.size() )
335 McH_simHMass->Fill(iAuxInfo->mcH[0].p4.M(), w);
336 if (iAuxInfo->mcH.size() )
337 McH_simHPt->Fill(iAuxInfo->mcH[0].p4.Pt(), w);
338 if (iAuxInfo->mcZ.size() )
339 McH_simZMass->Fill(iAuxInfo->mcZ[0].p4.M(), w);
340 if (iAuxInfo->mcZ.size() )
341 McH_simZPt->Fill(iAuxInfo->mcZ[0].p4.Pt(), w);
342 if (iAuxInfo->mcW.size() )
343 McH_simWMass->Fill(iAuxInfo->mcW[0].p4.M(), w);
344 if (iAuxInfo->mcW.size() )
345 McH_simWPt->Fill(iAuxInfo->mcW[0].p4.Pt(), w);
346 }
347 }
348
349
350 TH1F * McH_simHMass;
351 TH1F * McH_simHPt;
352 TH1F * McH_simWMass;
353 TH1F * McH_simWPt;
354 TH1F * McH_simZMass;
355 TH1F * McH_simZPt;
356
357 private:
358
359 Int_t bin_mass;
360 Double_t min_mass;
361 Double_t max_mass;
362
363 Int_t bin_pt;
364 Double_t min_pt;
365 Double_t max_pt;
366
367 Int_t bin_hel;
368 Double_t min_hel;
369 Double_t max_hel;
370
371
372 };
373
374
375 class BTagHistos : public Histos {
376
377 public:
378
379 virtual void book(TFile &f, std::string suffix) {
380
381 TDirectory *subDir;
382 if( ! f.GetDirectory(suffix.c_str()) )
383 subDir = f.mkdir(suffix.c_str());
384 else
385 subDir = f.GetDirectory(suffix.c_str());
386 subDir->cd();
387
388 bin_btag_prob = 20;
389 min_btag_prob = 0;
390 max_btag_prob = 1;
391
392 bin_btag_count = 10;
393 min_btag_count = 0;
394 max_btag_count = 10;
395
396 //Candidates
397 BTagH_bJet1_csv = new TH1F(("BJet1_CSV"+suffix).c_str(),("BJet1 CSV ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
398 BTagH_bJet2_csv = new TH1F(("BJet2_CSV"+suffix).c_str(),("BJet2 CSV ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
399 BTagH_bJet1_csvmva = new TH1F(("BJet1_CSVMVA"+suffix).c_str(),("BJet1 CSVMVA ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
400 BTagH_bJet2_csvmva = new TH1F(("BJet2_CSVMVA"+suffix).c_str(),("BJet2 CSVMVA ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
401 BTagH_bJet1_ssvhe = new TH1F(("BJet1_SSVHE"+suffix).c_str(),("BJet1 SSVHE ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
402 BTagH_bJet2_ssvhe = new TH1F(("BJet2_SSVHE"+suffix).c_str(),("BJet2 SSVHE ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
403 BTagH_bJet1_jpb = new TH1F(("BJet1_JPB"+suffix).c_str(),("BJet1 JPB ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
404 BTagH_bJet2_jpb = new TH1F(("BJet2_JPB"+suffix).c_str(),("BJet2 JPB ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
405 BTagH_bJet1_tche = new TH1F(("BJet1_TCHE"+suffix).c_str(),("BJet1 TCHE ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
406 BTagH_bJet2_tche = new TH1F(("BJet2_TCHE"+suffix).c_str(),("BJet2 TCHE ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
407 BTagH_bJet1_jp = new TH1F(("BJet1_JP"+suffix).c_str(),("BJet1 JP ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
408 BTagH_bJet2_jp = new TH1F(("BJet2_JP"+suffix).c_str(),("BJet2 JP ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
409 BTagH_bJet1_tchp = new TH1F(("BJet1_TCHP"+suffix).c_str(),("BJet1 TCHP ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
410 BTagH_bJet2_tchp = new TH1F(("BJet2_TCHP"+suffix).c_str(),("BJet2 TCHP ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
411
412 }
413
414 virtual void fill(VHbbProxy &iProxy,float w) {
415
416 const VHbbEvent *iEvent = iProxy.getVHbbEvent();
417 const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
418
419 //Candidates
420 if(iCand->size() > 0){
421 VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
422 VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
423
424 BTagH_bJet1_csv->Fill(H.jets.at(0).csv, w);
425 BTagH_bJet2_csv->Fill(H.jets.at(1).csv, w);
426 BTagH_bJet1_csvmva->Fill(H.jets.at(0).csvmva, w);
427 BTagH_bJet2_csvmva->Fill(H.jets.at(1).csvmva, w);
428 BTagH_bJet1_ssvhe->Fill(H.jets.at(0).ssvhe, w);
429 BTagH_bJet2_ssvhe->Fill(H.jets.at(1).ssvhe, w);
430 BTagH_bJet1_tche->Fill(H.jets.at(0).tche, w);
431 BTagH_bJet2_tche->Fill(H.jets.at(1).tche, w);
432 BTagH_bJet1_tchp->Fill(H.jets.at(0).tchp, w);
433 BTagH_bJet2_tchp->Fill(H.jets.at(1).tchp, w);
434 BTagH_bJet1_jpb->Fill(H.jets.at(0).jpb, w);
435 BTagH_bJet2_jpb->Fill(H.jets.at(1).jpb, w);
436 BTagH_bJet1_jp->Fill(H.jets.at(0).jp, w);
437 BTagH_bJet2_jp->Fill(H.jets.at(1).jp, w);
438
439 }
440 }
441
442 TH1F * BTagH_bJet1_csv;
443 TH1F * BTagH_bJet2_csv;
444 TH1F * BTagH_bJet1_csvmva;
445 TH1F * BTagH_bJet2_csvmva;
446 TH1F * BTagH_bJet1_ssvhe;
447 TH1F * BTagH_bJet2_ssvhe;
448 TH1F * BTagH_bJet1_jpb;
449 TH1F * BTagH_bJet2_jpb;
450 TH1F * BTagH_bJet1_tche;
451 TH1F * BTagH_bJet2_tche;
452 TH1F * BTagH_bJet1_jp;
453 TH1F * BTagH_bJet2_jp;
454 TH1F * BTagH_bJet1_tchp;
455 TH1F * BTagH_bJet2_tchp;
456
457 private:
458
459 Int_t bin_btag_prob;
460 Double_t min_btag_prob;
461 Double_t max_btag_prob;
462
463 Int_t bin_btag_count;
464 Double_t min_btag_count;
465 Double_t max_btag_count;
466
467 };
468
469 class StandardHistos : public Histos {
470
471 public:
472
473 virtual void book(TFile &f, std::string suffix) {
474
475 TDirectory *subDir;
476 if( ! f.GetDirectory(suffix.c_str()) )
477 subDir = f.mkdir(suffix.c_str());
478 else
479 subDir = f.GetDirectory(suffix.c_str());
480 subDir->cd();
481
482 bin_mass = 500;
483 min_mass = 0;
484 max_mass = 300;
485
486 bin_pt = 500;
487 min_pt = 0;
488 max_pt = 500;
489
490 bin_hel = 50;
491 min_hel = 0;
492 max_hel = 1;
493
494 bin_btag = 20;
495 min_btag = 0;
496 max_btag = 1;
497
498 bin_deltaR = bin_deltaPhi = bin_deltaEta = 20;
499 min_deltaR = min_deltaPhi = min_deltaEta = 0;
500 max_deltaR = max_deltaPhi = max_deltaEta = 5;
501
502 //Candidates
503 StH_simpleJet1_pt = new TH1F(("SimpleJet1_pt"+suffix).c_str(),("Simple Jet1 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
504 StH_simpleJet2_pt = new TH1F(("SimpleJet2_pt"+suffix).c_str(),("Simple Jet2 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
505 StH_simpleJet1_bTag = new TH1F(("SimpleJet1_bTag"+suffix).c_str(),("Simple Jet1 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
506 StH_simpleJet2_bTag = new TH1F(("SimpleJet2_bTag"+suffix).c_str(),("Simple Jet2 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
507 StH_simpleJets_dR = new TH1F(("SimpleJets_dR"+suffix).c_str(),("Simple Jets deltaR ("+suffix+")").c_str(), bin_deltaR, min_deltaR, max_deltaR );
508 StH_simpleJets_dPhi = new TH1F(("SimpleJets_dPhi"+suffix).c_str(),("Simple Jets deltaPhi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
509 StH_simpleJets_dEta = new TH1F(("SimpleJets_dEta"+suffix).c_str(),("Simple Jets deltaEta ("+suffix+")").c_str(), bin_deltaEta, min_deltaEta, max_deltaEta );
510
511 StH_HMass = new TH1F(("HiggsMass"+suffix).c_str(),(" Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
512 StH_HPt = new TH1F(("HiggsPt"+suffix).c_str(),(" Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
513 StH_HHel = new TH1F(("HiggsHel"+suffix).c_str(),("Higgs helicity angle ("+suffix+")").c_str(), bin_hel, min_hel, max_hel );
514 StH_HPullAngle = new TH1F(("HiggsPullAngle"+suffix).c_str(),("Higgs pull angle ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
515
516 StH_ZMass = new TH1F(("ZMass"+suffix).c_str(),(" Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
517 StH_ZPt = new TH1F(("ZPt"+suffix).c_str(),(" Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
518 StH_ZH_dPhi = new TH1F(("ZH_dPhi"+suffix).c_str(),(" ZH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
519
520 StH_WMass = new TH1F(("WMass"+suffix).c_str(),(" W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
521 StH_WPt = new TH1F(("WPt"+suffix).c_str(),(" W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
522 StH_WH_dPhi = new TH1F(("WH_dPhi"+suffix).c_str(),(" WH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
523
524 }
525
526 virtual void fill(VHbbProxy &iProxy,float w) {
527
528 const VHbbEvent *iEvent = iProxy.getVHbbEvent();
529 const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
530
531 //Candidates
532 if(iCand->size() > 0){
533 VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
534 VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
535 VHbbCandidate::VectorCandidate V = iCand->at(0).V;
536
537 StH_simpleJet1_pt->Fill(H.jets.at(0).p4.Pt(), w);
538 StH_simpleJet2_pt->Fill(H.jets.at(1).p4.Pt(), w);
539 StH_simpleJet1_bTag->Fill(H.jets.at(0).csv, w);
540 StH_simpleJet2_bTag->Fill(H.jets.at(1).csv, w);
541 StH_simpleJets_dR->Fill(H.jets.at(0).p4.DeltaR(H.jets.at(1).p4), w);
542 StH_simpleJets_dPhi->Fill(H.jets.at(0).p4.DeltaPhi(H.jets.at(1).p4), w);
543 StH_simpleJets_dEta->Fill(TMath::Abs(H.jets.at(0).p4.Eta()-H.jets.at(1).p4.Eta()), w);
544
545 StH_HMass->Fill(H.p4.M(), w);
546 StH_HPt->Fill(H.p4.Pt(), w);
547 // StH_HHel->Fill(H.hel(), w);
548 StH_HPullAngle->Fill(H.deltaTheta, w);
549 if( iCandType == VHbbCandidate::Zmumu || iCandType == VHbbCandidate::Zee || iCandType == VHbbCandidate::Znn ){
550 StH_ZMass->Fill(V.p4.M(), w);
551 StH_ZPt->Fill(V.p4.Pt(), w);
552 StH_ZH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
553 }
554 else if(iCandType == VHbbCandidate::Wen || iCandType == VHbbCandidate::Wmun){
555 StH_WMass->Fill(V.p4.M(), w);
556 StH_WPt->Fill(V.p4.Pt(), w);
557 StH_WH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
558 }
559
560 }
561 }
562
563 TH1F * StH_simpleJet1_pt;
564 TH1F * StH_simpleJet2_pt;
565 TH1F * StH_simpleJet1_bTag;
566 TH1F * StH_simpleJet2_bTag;
567 TH1F * StH_simpleJets_dR;
568 TH1F * StH_simpleJets_dPhi;
569 TH1F * StH_simpleJets_dEta;
570
571 TH1F * StH_HMass;
572 TH1F * StH_HPt;
573 TH1F * StH_HHel;
574 TH1F * StH_HPullAngle;
575 TH1F * StH_WMass;
576 TH1F * StH_WPt;
577 TH1F * StH_WH_dPhi;
578 TH1F * StH_ZMass;
579 TH1F * StH_ZPt;
580 TH1F * StH_ZH_dPhi;
581
582 private:
583
584 Int_t bin_btag;
585 Double_t min_btag;
586 Double_t max_btag;
587
588 Int_t bin_deltaEta;
589 Double_t min_deltaEta;
590 Double_t max_deltaEta;
591
592 Int_t bin_deltaPhi;
593 Double_t min_deltaPhi;
594 Double_t max_deltaPhi;
595
596 Int_t bin_deltaR;
597 Double_t min_deltaR;
598 Double_t max_deltaR;
599
600 Int_t bin_mass;
601 Double_t min_mass;
602 Double_t max_mass;
603
604 Int_t bin_pt;
605 Double_t min_pt;
606 Double_t max_pt;
607
608 Int_t bin_hel;
609 Double_t min_hel;
610 Double_t max_hel;
611
612 };
613
614 class HardJetHistos : public Histos {
615
616 public:
617
618 virtual void book(TFile &f, std::string suffix) {
619
620 TDirectory *subDir;
621 if( ! f.GetDirectory(suffix.c_str()) )
622 subDir = f.mkdir(suffix.c_str());
623 else
624 subDir = f.GetDirectory(suffix.c_str());
625 subDir->cd();
626
627 bin_mass = 500;
628 min_mass = 0;
629 max_mass = 300;
630
631 bin_pt = 500;
632 min_pt = 0;
633 max_pt = 500;
634
635 bin_hel = 50;
636 min_hel = 0;
637 max_hel = 1;
638
639 bin_btag = 20;
640 min_btag = 0;
641 max_btag = 1;
642
643 bin_deltaR = bin_deltaPhi = bin_deltaEta = 20;
644 min_deltaR = min_deltaPhi = min_deltaEta = 0;
645 max_deltaR = max_deltaPhi = max_deltaEta = 5;
646
647 //Candidates
648
649 HardJetH_subJet1_pt = new TH1F(("SubJet1_pt"+suffix).c_str(),("Sub Jet1 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
650 HardJetH_subJet2_pt = new TH1F(("SubJet2_pt"+suffix).c_str(),("Sub Jet2 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
651 HardJetH_subJet1_bTag = new TH1F(("SubJet1_bTag"+suffix).c_str(),("Sub Jet1 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
652 HardJetH_subJet2_bTag = new TH1F(("SubJet2_bTag"+suffix).c_str(),("Sub Jet2 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
653 HardJetH_subJets_dR = new TH1F(("SubJets_dR"+suffix).c_str(),("Sub Jets deltaR ("+suffix+")").c_str(), bin_deltaR, min_deltaR, max_deltaR );
654 HardJetH_subJets_dPhi = new TH1F(("SubJets_dPhi"+suffix).c_str(),("Sub Jets deltaPhi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
655 HardJetH_subJets_dEta = new TH1F(("SubJets_dEta"+suffix).c_str(),("Sub Jets deltaEta ("+suffix+")").c_str(), bin_deltaEta, min_deltaEta, max_deltaEta );
656
657 HardJetH_HMass = new TH1F(("HiggsMass"+suffix).c_str(),(" Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
658 HardJetH_HPt = new TH1F(("HiggsPt"+suffix).c_str(),(" Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
659 HardJetH_HHel = new TH1F(("HiggsHel"+suffix).c_str(),("Higgs helicity angle ("+suffix+")").c_str(), bin_hel, min_hel, max_hel );
660
661 HardJetH_ZMass = new TH1F(("ZMass"+suffix).c_str(),(" Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
662 HardJetH_ZPt = new TH1F(("ZPt"+suffix).c_str(),(" Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
663 HardJetH_ZH_dPhi = new TH1F(("ZH_dPhi"+suffix).c_str(),(" ZH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
664
665 HardJetH_WMass = new TH1F(("WMass"+suffix).c_str(),(" W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
666 HardJetH_WPt = new TH1F(("WPt"+suffix).c_str(),(" W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
667 HardJetH_WH_dPhi = new TH1F(("WH_dPhi"+suffix).c_str(),(" WH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
668
669 }
670
671 virtual void fill(VHbbProxy &iProxy,float w) {
672
673 const VHbbEvent *iEvent = iProxy.getVHbbEvent();
674 if(iEvent)
675 { const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
676
677 //Candidates
678 if(iCand->size() > 0){
679 VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
680 VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
681 VHbbCandidate::VectorCandidate V = iCand->at(0).V;
682 std::vector<VHbbEvent::HardJet> iHardJets = iEvent->hardJets;
683 VHbbEvent::HardJet iHardJet = iHardJets.at(0);
684
685 HardJetH_subJet1_pt->Fill(iHardJet.subFourMomentum.at(0).Pt(), w);
686 HardJetH_subJet2_pt->Fill(iHardJet.subFourMomentum.at(1).Pt(), w);
687 //SubJet information on btag missing
688 // HardJetH_subJet1_bTag->Fill(iHardJet.at(0).csv, w);
689 // HardJetH_subJet2_bTag->Fill(iHardJet.at(1).csv, w);
690 HardJetH_subJets_dR->Fill(iHardJet.subFourMomentum.at(0).DeltaR(iHardJet.subFourMomentum.at(1)), w);
691 HardJetH_subJets_dPhi->Fill(iHardJet.subFourMomentum.at(0).DeltaPhi(iHardJet.subFourMomentum.at(1)), w);
692 HardJetH_subJets_dEta->Fill(TMath::Abs(iHardJet.subFourMomentum.at(0).Eta()-iHardJet.subFourMomentum.at(1).Eta()), w);
693
694 //Here there should be the higgs candidate from HardJet
695 // HardJetH_HMass->Fill(H.p4.M(), w);
696 // HardJetH_HPt->Fill(H.p4.Pt(), w);
697 // // HardJetH_HHel->Fill(H.hel(), w);
698 // if( iCandType == VHbbCandidate::Zmumu || iCandType == VHbbCandidate::Zee || iCandType == VHbbCandidate::Znn ){
699 // HardJetH_ZMass->Fill(V.p4.M(), w);
700 // HardJetH_ZPt->Fill(V.p4.Pt(), w);
701 // HardJetH_ZH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
702 // }
703 // else if(iCandType == VHbbCandidate::Wen || iCandType == VHbbCandidate::Wmun){
704 // HardJetH_WMass->Fill(V.p4.M(), w);
705 // HardJetH_WPt->Fill(V.p4.Pt(), w);
706 // HardJetH_WH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
707 // }
708 }
709 }
710 }
711
712 TH1F * HardJetH_subJet1_pt;
713 TH1F * HardJetH_subJet2_pt;
714 TH1F * HardJetH_subJet1_bTag;
715 TH1F * HardJetH_subJet2_bTag;
716 TH1F * HardJetH_subJets_dR;
717 TH1F * HardJetH_subJets_dPhi;
718 TH1F * HardJetH_subJets_dEta;
719
720 TH1F * HardJetH_HMass;
721 TH1F * HardJetH_HPt;
722 TH1F * HardJetH_HHel;
723 TH1F * HardJetH_WMass;
724 TH1F * HardJetH_WPt;
725 TH1F * HardJetH_WH_dPhi;
726 TH1F * HardJetH_ZMass;
727 TH1F * HardJetH_ZPt;
728 TH1F * HardJetH_ZH_dPhi;
729
730 private:
731
732 Int_t bin_btag;
733 Double_t min_btag;
734 Double_t max_btag;
735
736 Int_t bin_deltaEta;
737 Double_t min_deltaEta;
738 Double_t max_deltaEta;
739
740 Int_t bin_deltaPhi;
741 Double_t min_deltaPhi;
742 Double_t max_deltaPhi;
743
744 Int_t bin_deltaR;
745 Double_t min_deltaR;
746 Double_t max_deltaR;
747
748 Int_t bin_mass;
749 Double_t min_mass;
750 Double_t max_mass;
751
752 Int_t bin_pt;
753 Double_t min_pt;
754 Double_t max_pt;
755
756 Int_t bin_hel;
757 Double_t min_hel;
758 Double_t max_hel;
759
760 };