ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/VHbbAnalysis/VHbbDataFormats/src/CutsAndHistos.cc
Revision: 1.3
Committed: Sun Jul 17 11:06:02 2011 UTC (13 years, 9 months ago) by arizzi
Content type: text/plain
Branch: MAIN
CVS Tags: Jul21st2011, Jul20th2011, Jul18th2011, Jul17th2011
Changes since 1.2: +6 -5 lines
Log Message:
tests using only the Candidate branch

File Contents

# User Rev Content
1 tboccali 1.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 bortigno 1.2 class mcHPtCut : public Cut {
13 tboccali 1.1 public:
14 bortigno 1.2 mcHPtCut(Double_t mcHPtCutMin_):
15     mcHPtCutMin(mcHPtCutMin_) {}
16     mcHPtCut(Double_t mcHPtCutMin_, Double_t mcHPtCutMax_):
17     mcHPtCutMin(mcHPtCutMin_),
18     mcHPtCutMax(mcHPtCutMax_) {}
19 tboccali 1.1 std::string name() {
20 bortigno 1.2 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 tboccali 1.1 }
28     bool pass(VHbbProxy &iProxy) {
29 bortigno 1.2 if(!mcHPtCutMax)
30     return ( iProxy.getVHbbEvent()->mcH.fourMomentum.Pt() > mcHPtCutMin );
31     else
32     return ( iProxy.getVHbbEvent()->mcH.fourMomentum.Pt() > mcHPtCutMin
33     && iProxy.getVHbbEvent()->mcH.fourMomentum.Pt() < mcHPtCutMax );
34 tboccali 1.1 }
35     private:
36 bortigno 1.2 Double_t mcHPtCutMin;
37     Double_t mcHPtCutMax;
38     std::ostringstream oss_mcHPtCutMin;
39     std::ostringstream oss_mcHPtCutMax;
40 tboccali 1.1 };
41    
42 bortigno 1.2
43     class mcWPtCut : public Cut {
44 tboccali 1.1 public:
45 bortigno 1.2 mcWPtCut(Double_t mcWPtCutMin_):
46     mcWPtCutMin(mcWPtCutMin_) {}
47     mcWPtCut(Double_t mcWPtCutMin_, Double_t mcWPtCutMax_):
48     mcWPtCutMin(mcWPtCutMin_),
49     mcWPtCutMax(mcWPtCutMax_) {}
50 tboccali 1.1 std::string name() {
51 bortigno 1.2 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 tboccali 1.1 }
59     bool pass(VHbbProxy &iProxy) {
60 bortigno 1.2 if(!mcWPtCutMax)
61     return ( iProxy.getVHbbEvent()->mcW.fourMomentum.Pt() > mcWPtCutMin );
62     else
63     return ( iProxy.getVHbbEvent()->mcW.fourMomentum.Pt() > mcWPtCutMin
64     && iProxy.getVHbbEvent()->mcW.fourMomentum.Pt() < mcWPtCutMax );
65 tboccali 1.1 }
66     private:
67 bortigno 1.2 Double_t mcWPtCutMin;
68     Double_t mcWPtCutMax;
69     std::ostringstream oss_mcWPtCutMin;
70     std::ostringstream oss_mcWPtCutMax;
71 tboccali 1.1 };
72    
73 bortigno 1.2
74    
75     class mcZPtCut : public Cut {
76 tboccali 1.1 public:
77 bortigno 1.2 mcZPtCut(Double_t mcZPtCutMin_):
78     mcZPtCutMin(mcZPtCutMin_) {}
79     mcZPtCut(Double_t mcZPtCutMin_, Double_t mcZPtCutMax_):
80     mcZPtCutMin(mcZPtCutMin_),
81     mcZPtCutMax(mcZPtCutMax_) {}
82 tboccali 1.1 std::string name() {
83 bortigno 1.2 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 tboccali 1.1 }
91     bool pass(VHbbProxy &iProxy) {
92 bortigno 1.2 if(!mcZPtCutMax)
93     return ( iProxy.getVHbbEvent()->mcZ.fourMomentum.Pt() > mcZPtCutMin );
94     else
95     return ( iProxy.getVHbbEvent()->mcZ.fourMomentum.Pt() > mcZPtCutMin
96     && iProxy.getVHbbEvent()->mcZ.fourMomentum.Pt() < mcZPtCutMax );
97     }
98 tboccali 1.1 private:
99 bortigno 1.2 Double_t mcZPtCutMin;
100     Double_t mcZPtCutMax;
101     std::ostringstream oss_mcZPtCutMin;
102     std::ostringstream oss_mcZPtCutMax;
103 tboccali 1.1 };
104    
105 bortigno 1.2
106     class SimpleJet1PtCut: public Cut {
107 tboccali 1.1 public:
108 bortigno 1.2 SimpleJet1PtCut(Double_t simpleJet1PtCutMin_):
109     simpleJet1PtCutMin(simpleJet1PtCutMin_) {}
110     SimpleJet1PtCut(Double_t simpleJet1PtCutMin_, Double_t simpleJet1PtCutMax_):
111     simpleJet1PtCutMin(simpleJet1PtCutMin_),
112     simpleJet1PtCutMax(simpleJet1PtCutMax_) {}
113 tboccali 1.1 std::string name() {
114 bortigno 1.2 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 tboccali 1.1 }
123     bool pass(VHbbProxy &iProxy) {
124     if(iProxy.getVHbbCandidate()->size() < 1)
125     return false;
126     else
127 bortigno 1.2 if(!simpleJet1PtCutMin)
128     return ((iProxy.getVHbbCandidate()->at(0)).H.jets.at(0).fourMomentum.Pt() > simpleJet1PtCutMin);
129     else
130     return ((iProxy.getVHbbCandidate()->at(0)).H.jets.at(0).fourMomentum.Pt() > simpleJet1PtCutMin
131     && (iProxy.getVHbbCandidate()->at(0)).H.jets.at(0).fourMomentum.Pt() < simpleJet1PtCutMax);
132     }
133 tboccali 1.1 private:
134 bortigno 1.2 Double_t simpleJet1PtCutMin;
135     Double_t simpleJet1PtCutMax;
136     std::ostringstream oss_simpleJet1PtCutMin;
137     std::ostringstream oss_simpleJet1PtCutMax;
138 tboccali 1.1 };
139    
140 bortigno 1.2 class VPtCut: public Cut {
141 tboccali 1.1 public:
142 bortigno 1.2 VPtCut(Double_t VptCutMin_):
143     VptCutMin(VptCutMin_) {}
144     VPtCut(Double_t VptCutMin_, Double_t VptCutMax_):
145     VptCutMin(VptCutMin_),
146     VptCutMax(VptCutMax_) {}
147 tboccali 1.1 std::string name() {
148 bortigno 1.2 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 tboccali 1.1 }
157     bool pass(VHbbProxy &iProxy) {
158     if(iProxy.getVHbbCandidate()->size() < 1)
159     return false;
160     else
161 bortigno 1.2 if(VptCutMax != 1e10)
162     return ((iProxy.getVHbbCandidate()->at(0)).V.fourMomentum.Pt() > VptCutMin
163     && (iProxy.getVHbbCandidate()->at(0)).V.fourMomentum.Pt() < VptCutMax);
164     else
165     return ((iProxy.getVHbbCandidate()->at(0)).V.fourMomentum.Pt() > VptCutMin );
166    
167     }
168 tboccali 1.1 private:
169 bortigno 1.2 Double_t VptCutMin;
170     Double_t VptCutMax;
171     std::ostringstream oss_VptMin;
172     std::ostringstream oss_VptMax;
173 tboccali 1.1 };
174    
175 bortigno 1.2
176     class HPtCut: public Cut {
177 tboccali 1.1 public:
178 bortigno 1.2 HPtCut(Double_t hPtCutMin_):
179     hPtCutMin(hPtCutMin_) {}
180     HPtCut(Double_t hPtCutMin_, Double_t hPtCutMax_):
181     hPtCutMin(hPtCutMin_),
182     hPtCutMax(hPtCutMax_) {}
183    
184 tboccali 1.1 std::string name() {
185 bortigno 1.2 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 tboccali 1.1 }
194     bool pass(VHbbProxy &iProxy) {
195     if(iProxy.getVHbbCandidate()->size() < 1)
196     return false;
197     else
198 bortigno 1.2 if(!hPtCutMax)
199     return ((iProxy.getVHbbCandidate()->at(0)).H.fourMomentum.Pt() > hPtCutMin);
200     else
201     return ((iProxy.getVHbbCandidate()->at(0)).H.fourMomentum.Pt() > hPtCutMin
202     && (iProxy.getVHbbCandidate()->at(0)).H.fourMomentum.Pt() < hPtCutMax);
203 tboccali 1.1 }
204     private:
205 bortigno 1.2 Double_t hPtCutMin;
206     Double_t hPtCutMax;
207     std::ostringstream oss_HPtCutMin;
208     std::ostringstream oss_HPtCutMax;
209 tboccali 1.1 };
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.fourMomentum.Pt() > Higgs_pt
260     && V.fourMomentum.Pt() > V_pt
261     && TMath::Abs( Geom::deltaPhi(H.fourMomentum.Phi(), V.fourMomentum.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 arizzi 1.3 if(iEvent)
332     {
333 tboccali 1.1 //from MC
334     McH_simHMass->Fill(iEvent->mcH.fourMomentum.M(), w);
335     McH_simHPt->Fill(iEvent->mcH.fourMomentum.Pt(), w);
336     McH_simZMass->Fill(iEvent->mcZ.fourMomentum.M(), w);
337     McH_simZPt->Fill(iEvent->mcZ.fourMomentum.Pt(), w);
338     McH_simWMass->Fill(iEvent->mcW.fourMomentum.M(), w);
339     McH_simWPt->Fill(iEvent->mcW.fourMomentum.Pt(), w);
340 arizzi 1.3 }
341 tboccali 1.1 }
342    
343    
344     TH1F * McH_simHMass;
345     TH1F * McH_simHPt;
346     TH1F * McH_simWMass;
347     TH1F * McH_simWPt;
348     TH1F * McH_simZMass;
349     TH1F * McH_simZPt;
350    
351     private:
352    
353     Int_t bin_mass;
354     Double_t min_mass;
355     Double_t max_mass;
356    
357     Int_t bin_pt;
358     Double_t min_pt;
359     Double_t max_pt;
360    
361     Int_t bin_hel;
362     Double_t min_hel;
363     Double_t max_hel;
364    
365    
366     };
367    
368    
369     class BTagHistos : public Histos {
370    
371     public:
372    
373     virtual void book(TFile &f, std::string suffix) {
374    
375     TDirectory *subDir;
376     if( ! f.GetDirectory(suffix.c_str()) )
377     subDir = f.mkdir(suffix.c_str());
378     else
379     subDir = f.GetDirectory(suffix.c_str());
380     subDir->cd();
381    
382     bin_btag_prob = 20;
383     min_btag_prob = 0;
384     max_btag_prob = 1;
385    
386     bin_btag_count = 10;
387     min_btag_count = 0;
388     max_btag_count = 10;
389    
390     //Candidates
391     BTagH_bJet1_csv = new TH1F(("BJet1_CSV"+suffix).c_str(),("BJet1 CSV ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
392     BTagH_bJet2_csv = new TH1F(("BJet2_CSV"+suffix).c_str(),("BJet2 CSV ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
393     BTagH_bJet1_csvmva = new TH1F(("BJet1_CSVMVA"+suffix).c_str(),("BJet1 CSVMVA ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
394     BTagH_bJet2_csvmva = new TH1F(("BJet2_CSVMVA"+suffix).c_str(),("BJet2 CSVMVA ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
395     BTagH_bJet1_ssvhe = new TH1F(("BJet1_SSVHE"+suffix).c_str(),("BJet1 SSVHE ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
396     BTagH_bJet2_ssvhe = new TH1F(("BJet2_SSVHE"+suffix).c_str(),("BJet2 SSVHE ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
397     BTagH_bJet1_jpb = new TH1F(("BJet1_JPB"+suffix).c_str(),("BJet1 JPB ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
398     BTagH_bJet2_jpb = new TH1F(("BJet2_JPB"+suffix).c_str(),("BJet2 JPB ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
399     BTagH_bJet1_tche = new TH1F(("BJet1_TCHE"+suffix).c_str(),("BJet1 TCHE ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
400     BTagH_bJet2_tche = new TH1F(("BJet2_TCHE"+suffix).c_str(),("BJet2 TCHE ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
401     BTagH_bJet1_jp = new TH1F(("BJet1_JP"+suffix).c_str(),("BJet1 JP ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
402     BTagH_bJet2_jp = new TH1F(("BJet2_JP"+suffix).c_str(),("BJet2 JP ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
403     BTagH_bJet1_tchp = new TH1F(("BJet1_TCHP"+suffix).c_str(),("BJet1 TCHP ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
404     BTagH_bJet2_tchp = new TH1F(("BJet2_TCHP"+suffix).c_str(),("BJet2 TCHP ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
405    
406     }
407    
408     virtual void fill(VHbbProxy &iProxy,float w) {
409    
410     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
411     const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
412    
413     //Candidates
414     if(iCand->size() > 0){
415     VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
416     VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
417    
418     BTagH_bJet1_csv->Fill(H.jets.at(0).csv, w);
419     BTagH_bJet2_csv->Fill(H.jets.at(1).csv, w);
420     BTagH_bJet1_csvmva->Fill(H.jets.at(0).csvmva, w);
421     BTagH_bJet2_csvmva->Fill(H.jets.at(1).csvmva, w);
422     BTagH_bJet1_ssvhe->Fill(H.jets.at(0).ssvhe, w);
423     BTagH_bJet2_ssvhe->Fill(H.jets.at(1).ssvhe, w);
424     BTagH_bJet1_tche->Fill(H.jets.at(0).tche, w);
425     BTagH_bJet2_tche->Fill(H.jets.at(1).tche, w);
426     BTagH_bJet1_tchp->Fill(H.jets.at(0).tchp, w);
427     BTagH_bJet2_tchp->Fill(H.jets.at(1).tchp, w);
428     BTagH_bJet1_jpb->Fill(H.jets.at(0).jpb, w);
429     BTagH_bJet2_jpb->Fill(H.jets.at(1).jpb, w);
430     BTagH_bJet1_jp->Fill(H.jets.at(0).jp, w);
431     BTagH_bJet2_jp->Fill(H.jets.at(1).jp, w);
432    
433     }
434     }
435    
436     TH1F * BTagH_bJet1_csv;
437     TH1F * BTagH_bJet2_csv;
438     TH1F * BTagH_bJet1_csvmva;
439     TH1F * BTagH_bJet2_csvmva;
440     TH1F * BTagH_bJet1_ssvhe;
441     TH1F * BTagH_bJet2_ssvhe;
442     TH1F * BTagH_bJet1_jpb;
443     TH1F * BTagH_bJet2_jpb;
444     TH1F * BTagH_bJet1_tche;
445     TH1F * BTagH_bJet2_tche;
446     TH1F * BTagH_bJet1_jp;
447     TH1F * BTagH_bJet2_jp;
448     TH1F * BTagH_bJet1_tchp;
449     TH1F * BTagH_bJet2_tchp;
450    
451     private:
452    
453     Int_t bin_btag_prob;
454     Double_t min_btag_prob;
455     Double_t max_btag_prob;
456    
457     Int_t bin_btag_count;
458     Double_t min_btag_count;
459     Double_t max_btag_count;
460    
461     };
462    
463     class StandardHistos : public Histos {
464    
465     public:
466    
467     virtual void book(TFile &f, std::string suffix) {
468    
469     TDirectory *subDir;
470     if( ! f.GetDirectory(suffix.c_str()) )
471     subDir = f.mkdir(suffix.c_str());
472     else
473     subDir = f.GetDirectory(suffix.c_str());
474     subDir->cd();
475    
476     bin_mass = 500;
477     min_mass = 0;
478     max_mass = 300;
479    
480     bin_pt = 500;
481     min_pt = 0;
482     max_pt = 500;
483    
484     bin_hel = 50;
485     min_hel = 0;
486     max_hel = 1;
487    
488     bin_btag = 20;
489     min_btag = 0;
490     max_btag = 1;
491    
492     bin_deltaR = bin_deltaPhi = bin_deltaEta = 20;
493     min_deltaR = min_deltaPhi = min_deltaEta = 0;
494     max_deltaR = max_deltaPhi = max_deltaEta = 5;
495    
496     //Candidates
497     StH_simpleJet1_pt = new TH1F(("SimpleJet1_pt"+suffix).c_str(),("Simple Jet1 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
498     StH_simpleJet2_pt = new TH1F(("SimpleJet2_pt"+suffix).c_str(),("Simple Jet2 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
499     StH_simpleJet1_bTag = new TH1F(("SimpleJet1_bTag"+suffix).c_str(),("Simple Jet1 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
500     StH_simpleJet2_bTag = new TH1F(("SimpleJet2_bTag"+suffix).c_str(),("Simple Jet2 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
501     StH_simpleJets_dR = new TH1F(("SimpleJets_dR"+suffix).c_str(),("Simple Jets deltaR ("+suffix+")").c_str(), bin_deltaR, min_deltaR, max_deltaR );
502     StH_simpleJets_dPhi = new TH1F(("SimpleJets_dPhi"+suffix).c_str(),("Simple Jets deltaPhi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
503     StH_simpleJets_dEta = new TH1F(("SimpleJets_dEta"+suffix).c_str(),("Simple Jets deltaEta ("+suffix+")").c_str(), bin_deltaEta, min_deltaEta, max_deltaEta );
504    
505     StH_HMass = new TH1F(("HiggsMass"+suffix).c_str(),(" Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
506     StH_HPt = new TH1F(("HiggsPt"+suffix).c_str(),(" Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
507     StH_HHel = new TH1F(("HiggsHel"+suffix).c_str(),("Higgs helicity angle ("+suffix+")").c_str(), bin_hel, min_hel, max_hel );
508     StH_HPullAngle = new TH1F(("HiggsPullAngle"+suffix).c_str(),("Higgs pull angle ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
509    
510     StH_ZMass = new TH1F(("ZMass"+suffix).c_str(),(" Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
511     StH_ZPt = new TH1F(("ZPt"+suffix).c_str(),(" Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
512     StH_ZH_dPhi = new TH1F(("ZH_dPhi"+suffix).c_str(),(" ZH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
513    
514     StH_WMass = new TH1F(("WMass"+suffix).c_str(),(" W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
515     StH_WPt = new TH1F(("WPt"+suffix).c_str(),(" W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
516     StH_WH_dPhi = new TH1F(("WH_dPhi"+suffix).c_str(),(" WH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
517    
518     }
519    
520     virtual void fill(VHbbProxy &iProxy,float w) {
521    
522     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
523     const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
524    
525     //Candidates
526     if(iCand->size() > 0){
527     VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
528     VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
529     VHbbCandidate::VectorCandidate V = iCand->at(0).V;
530    
531     StH_simpleJet1_pt->Fill(H.jets.at(0).fourMomentum.Pt(), w);
532     StH_simpleJet2_pt->Fill(H.jets.at(1).fourMomentum.Pt(), w);
533     StH_simpleJet1_bTag->Fill(H.jets.at(0).csv, w);
534     StH_simpleJet2_bTag->Fill(H.jets.at(1).csv, w);
535     StH_simpleJets_dR->Fill(H.jets.at(0).fourMomentum.DeltaR(H.jets.at(1).fourMomentum), w);
536     StH_simpleJets_dPhi->Fill(H.jets.at(0).fourMomentum.DeltaPhi(H.jets.at(1).fourMomentum), w);
537     StH_simpleJets_dEta->Fill(TMath::Abs(H.jets.at(0).fourMomentum.Eta()-H.jets.at(1).fourMomentum.Eta()), w);
538    
539     StH_HMass->Fill(H.fourMomentum.M(), w);
540     StH_HPt->Fill(H.fourMomentum.Pt(), w);
541     // StH_HHel->Fill(H.hel(), w);
542     StH_HPullAngle->Fill(H.deltaTheta, w);
543     if( iCandType == VHbbCandidate::Zmumu || iCandType == VHbbCandidate::Zee || iCandType == VHbbCandidate::Znn ){
544     StH_ZMass->Fill(V.fourMomentum.M(), w);
545     StH_ZPt->Fill(V.fourMomentum.Pt(), w);
546     StH_ZH_dPhi->Fill(V.fourMomentum.DeltaPhi(H.fourMomentum.Phi()), w);
547     }
548     else if(iCandType == VHbbCandidate::Wen || iCandType == VHbbCandidate::Wmun){
549     StH_WMass->Fill(V.fourMomentum.M(), w);
550     StH_WPt->Fill(V.fourMomentum.Pt(), w);
551     StH_WH_dPhi->Fill(V.fourMomentum.DeltaPhi(H.fourMomentum.Phi()), w);
552     }
553    
554     }
555     }
556    
557     TH1F * StH_simpleJet1_pt;
558     TH1F * StH_simpleJet2_pt;
559     TH1F * StH_simpleJet1_bTag;
560     TH1F * StH_simpleJet2_bTag;
561     TH1F * StH_simpleJets_dR;
562     TH1F * StH_simpleJets_dPhi;
563     TH1F * StH_simpleJets_dEta;
564    
565     TH1F * StH_HMass;
566     TH1F * StH_HPt;
567     TH1F * StH_HHel;
568     TH1F * StH_HPullAngle;
569     TH1F * StH_WMass;
570     TH1F * StH_WPt;
571     TH1F * StH_WH_dPhi;
572     TH1F * StH_ZMass;
573     TH1F * StH_ZPt;
574     TH1F * StH_ZH_dPhi;
575    
576     private:
577    
578     Int_t bin_btag;
579     Double_t min_btag;
580     Double_t max_btag;
581    
582     Int_t bin_deltaEta;
583     Double_t min_deltaEta;
584     Double_t max_deltaEta;
585    
586     Int_t bin_deltaPhi;
587     Double_t min_deltaPhi;
588     Double_t max_deltaPhi;
589    
590     Int_t bin_deltaR;
591     Double_t min_deltaR;
592     Double_t max_deltaR;
593    
594     Int_t bin_mass;
595     Double_t min_mass;
596     Double_t max_mass;
597    
598     Int_t bin_pt;
599     Double_t min_pt;
600     Double_t max_pt;
601    
602     Int_t bin_hel;
603     Double_t min_hel;
604     Double_t max_hel;
605    
606     };
607    
608     class HardJetHistos : public Histos {
609    
610     public:
611    
612     virtual void book(TFile &f, std::string suffix) {
613    
614     TDirectory *subDir;
615     if( ! f.GetDirectory(suffix.c_str()) )
616     subDir = f.mkdir(suffix.c_str());
617     else
618     subDir = f.GetDirectory(suffix.c_str());
619     subDir->cd();
620    
621     bin_mass = 500;
622     min_mass = 0;
623     max_mass = 300;
624    
625     bin_pt = 500;
626     min_pt = 0;
627     max_pt = 500;
628    
629     bin_hel = 50;
630     min_hel = 0;
631     max_hel = 1;
632    
633     bin_btag = 20;
634     min_btag = 0;
635     max_btag = 1;
636    
637     bin_deltaR = bin_deltaPhi = bin_deltaEta = 20;
638     min_deltaR = min_deltaPhi = min_deltaEta = 0;
639     max_deltaR = max_deltaPhi = max_deltaEta = 5;
640    
641     //Candidates
642    
643     HardJetH_subJet1_pt = new TH1F(("SubJet1_pt"+suffix).c_str(),("Sub Jet1 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
644     HardJetH_subJet2_pt = new TH1F(("SubJet2_pt"+suffix).c_str(),("Sub Jet2 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
645     HardJetH_subJet1_bTag = new TH1F(("SubJet1_bTag"+suffix).c_str(),("Sub Jet1 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
646     HardJetH_subJet2_bTag = new TH1F(("SubJet2_bTag"+suffix).c_str(),("Sub Jet2 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
647     HardJetH_subJets_dR = new TH1F(("SubJets_dR"+suffix).c_str(),("Sub Jets deltaR ("+suffix+")").c_str(), bin_deltaR, min_deltaR, max_deltaR );
648     HardJetH_subJets_dPhi = new TH1F(("SubJets_dPhi"+suffix).c_str(),("Sub Jets deltaPhi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
649     HardJetH_subJets_dEta = new TH1F(("SubJets_dEta"+suffix).c_str(),("Sub Jets deltaEta ("+suffix+")").c_str(), bin_deltaEta, min_deltaEta, max_deltaEta );
650    
651     HardJetH_HMass = new TH1F(("HiggsMass"+suffix).c_str(),(" Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
652     HardJetH_HPt = new TH1F(("HiggsPt"+suffix).c_str(),(" Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
653     HardJetH_HHel = new TH1F(("HiggsHel"+suffix).c_str(),("Higgs helicity angle ("+suffix+")").c_str(), bin_hel, min_hel, max_hel );
654    
655     HardJetH_ZMass = new TH1F(("ZMass"+suffix).c_str(),(" Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
656     HardJetH_ZPt = new TH1F(("ZPt"+suffix).c_str(),(" Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
657     HardJetH_ZH_dPhi = new TH1F(("ZH_dPhi"+suffix).c_str(),(" ZH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
658    
659     HardJetH_WMass = new TH1F(("WMass"+suffix).c_str(),(" W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
660     HardJetH_WPt = new TH1F(("WPt"+suffix).c_str(),(" W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
661     HardJetH_WH_dPhi = new TH1F(("WH_dPhi"+suffix).c_str(),(" WH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
662    
663     }
664    
665     virtual void fill(VHbbProxy &iProxy,float w) {
666    
667     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
668 arizzi 1.3 if(iEvent)
669     { const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
670 tboccali 1.1
671     //Candidates
672     if(iCand->size() > 0){
673     VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
674     VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
675     VHbbCandidate::VectorCandidate V = iCand->at(0).V;
676     std::vector<VHbbEvent::HardJet> iHardJets = iEvent->hardJets;
677     VHbbEvent::HardJet iHardJet = iHardJets.at(0);
678    
679     HardJetH_subJet1_pt->Fill(iHardJet.subFourMomentum.at(0).Pt(), w);
680     HardJetH_subJet2_pt->Fill(iHardJet.subFourMomentum.at(1).Pt(), w);
681     //SubJet information on btag missing
682     // HardJetH_subJet1_bTag->Fill(iHardJet.at(0).csv, w);
683     // HardJetH_subJet2_bTag->Fill(iHardJet.at(1).csv, w);
684     HardJetH_subJets_dR->Fill(iHardJet.subFourMomentum.at(0).DeltaR(iHardJet.subFourMomentum.at(1)), w);
685     HardJetH_subJets_dPhi->Fill(iHardJet.subFourMomentum.at(0).DeltaPhi(iHardJet.subFourMomentum.at(1)), w);
686     HardJetH_subJets_dEta->Fill(TMath::Abs(iHardJet.subFourMomentum.at(0).Eta()-iHardJet.subFourMomentum.at(1).Eta()), w);
687    
688     //Here there should be the higgs candidate from HardJet
689     // HardJetH_HMass->Fill(H.fourMomentum.M(), w);
690     // HardJetH_HPt->Fill(H.fourMomentum.Pt(), w);
691     // // HardJetH_HHel->Fill(H.hel(), w);
692     // if( iCandType == VHbbCandidate::Zmumu || iCandType == VHbbCandidate::Zee || iCandType == VHbbCandidate::Znn ){
693     // HardJetH_ZMass->Fill(V.fourMomentum.M(), w);
694     // HardJetH_ZPt->Fill(V.fourMomentum.Pt(), w);
695     // HardJetH_ZH_dPhi->Fill(V.fourMomentum.DeltaPhi(H.fourMomentum.Phi()), w);
696     // }
697     // else if(iCandType == VHbbCandidate::Wen || iCandType == VHbbCandidate::Wmun){
698     // HardJetH_WMass->Fill(V.fourMomentum.M(), w);
699     // HardJetH_WPt->Fill(V.fourMomentum.Pt(), w);
700     // HardJetH_WH_dPhi->Fill(V.fourMomentum.DeltaPhi(H.fourMomentum.Phi()), w);
701     // }
702 arizzi 1.3 }
703 tboccali 1.1 }
704     }
705    
706     TH1F * HardJetH_subJet1_pt;
707     TH1F * HardJetH_subJet2_pt;
708     TH1F * HardJetH_subJet1_bTag;
709     TH1F * HardJetH_subJet2_bTag;
710     TH1F * HardJetH_subJets_dR;
711     TH1F * HardJetH_subJets_dPhi;
712     TH1F * HardJetH_subJets_dEta;
713    
714     TH1F * HardJetH_HMass;
715     TH1F * HardJetH_HPt;
716     TH1F * HardJetH_HHel;
717     TH1F * HardJetH_WMass;
718     TH1F * HardJetH_WPt;
719     TH1F * HardJetH_WH_dPhi;
720     TH1F * HardJetH_ZMass;
721     TH1F * HardJetH_ZPt;
722     TH1F * HardJetH_ZH_dPhi;
723    
724     private:
725    
726     Int_t bin_btag;
727     Double_t min_btag;
728     Double_t max_btag;
729    
730     Int_t bin_deltaEta;
731     Double_t min_deltaEta;
732     Double_t max_deltaEta;
733    
734     Int_t bin_deltaPhi;
735     Double_t min_deltaPhi;
736     Double_t max_deltaPhi;
737    
738     Int_t bin_deltaR;
739     Double_t min_deltaR;
740     Double_t max_deltaR;
741    
742     Int_t bin_mass;
743     Double_t min_mass;
744     Double_t max_mass;
745    
746     Int_t bin_pt;
747     Double_t min_pt;
748     Double_t max_pt;
749    
750     Int_t bin_hel;
751     Double_t min_hel;
752     Double_t max_hel;
753    
754     };