ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/VHbbAnalysis/VHbbDataFormats/interface/Histos.h
Revision: 1.1
Committed: Mon Aug 22 11:53:09 2011 UTC (13 years, 8 months ago) by arizzi
Content type: text/plain
Branch: MAIN
Log Message:
move Histos and Cuts200X to .h remove CutAndHistos

File Contents

# User Rev Content
1 arizzi 1.1 #ifndef HISTOS_H
2     #define HISTOS_H
3    
4     #include "VHbbAnalysis/VHbbDataFormats/interface/VHbbEvent.h"
5     #include "VHbbAnalysis/VHbbDataFormats/interface/VHbbCandidate.h"
6     #include "VHbbAnalysis/VHbbDataFormats/interface/VHbbProxy.h"
7     #include "VHbbAnalysis/VHbbDataFormats/interface/CutsAndHistos.h"
8     #include <TH1F.h>
9     #include "DataFormats/GeometryVector/interface/VectorUtil.h"
10     #include "DataFormats/Math/interface/deltaR.h"
11     #include <sstream>
12     #include "TKey.h"
13     #include "TMVA/Reader.h"
14     #include "TMVA/Config.h"
15     #include "TMVA/Tools.h"
16     #include "TMVA/MethodCuts.h"
17    
18     class MCHistos : public Histos {
19    
20     public:
21    
22     virtual void book(TFile &f, std::string suffix) {
23    
24     TDirectory *subDir;
25    
26     if( ! f.GetDirectory(suffix.c_str()) )
27     subDir = f.mkdir(suffix.c_str());
28     else
29     subDir = f.GetDirectory(suffix.c_str());
30    
31     subDir->cd();
32    
33     bin_mass = 500;
34     min_mass = 0;
35     max_mass = 300;
36    
37     bin_pt = 500;
38     min_pt = 0;
39     max_pt = 500;
40    
41     bin_hel = 50;
42     min_hel = 0;
43     max_hel = 1;
44    
45     //from MC
46     McH_simHMass = new TH1F(("simHiggsMass"+suffix).c_str(),("Sim Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
47     McH_simHPt = new TH1F(("simHiggsPt"+suffix).c_str(),("Sim Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
48     McH_simZMass = new TH1F(("simZMass"+suffix).c_str(),("Sim Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
49     McH_simZPt = new TH1F(("simZPt"+suffix).c_str(),("Sim Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
50     McH_simWMass = new TH1F(("simWMass"+suffix).c_str(),("Sim W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
51     McH_simWPt = new TH1F(("simWPt"+suffix).c_str(),("Sim W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
52    
53     }
54    
55     virtual void fill(VHbbProxy &iProxy,float w) {
56    
57     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
58     const VHbbEventAuxInfo *iAuxInfo = iProxy.getVHbbEventAuxInfo();
59     if(iAuxInfo)
60     {
61     //from MC
62     if (iAuxInfo->mcH.size()!=0)
63     McH_simHMass->Fill(iAuxInfo->mcH[0].p4.M(), w);
64     if (iAuxInfo->mcH.size()!=0)
65     McH_simHPt->Fill(iAuxInfo->mcH[0].p4.Pt(), w);
66     if (iAuxInfo->mcZ.size()!=0)
67     McH_simZMass->Fill(iAuxInfo->mcZ[0].p4.M(), w);
68     if (iAuxInfo->mcZ.size()!=0)
69     McH_simZPt->Fill(iAuxInfo->mcZ[0].p4.Pt(), w);
70     if (iAuxInfo->mcW.size()!=0)
71     McH_simWMass->Fill(iAuxInfo->mcW[0].p4.M(), w);
72     if (iAuxInfo->mcW.size()!=0)
73     McH_simWPt->Fill(iAuxInfo->mcW[0].p4.Pt(), w);
74     }
75     }
76    
77    
78     TH1F * McH_simHMass;
79     TH1F * McH_simHPt;
80     TH1F * McH_simWMass;
81     TH1F * McH_simWPt;
82     TH1F * McH_simZMass;
83     TH1F * McH_simZPt;
84    
85     private:
86    
87     Int_t bin_mass;
88     Double_t min_mass;
89     Double_t max_mass;
90    
91     Int_t bin_pt;
92     Double_t min_pt;
93     Double_t max_pt;
94    
95     Int_t bin_hel;
96     Double_t min_hel;
97     Double_t max_hel;
98    
99    
100     };
101    
102    
103     class BTagHistos : public Histos {
104    
105     public:
106    
107     virtual void book(TFile &f, std::string suffix) {
108    
109     TDirectory *subDir;
110     if( ! f.GetDirectory(suffix.c_str()) )
111     subDir = f.mkdir(suffix.c_str());
112     else
113     subDir = f.GetDirectory(suffix.c_str());
114     subDir->cd();
115    
116     bin_btag_prob = 20;
117     min_btag_prob = 0;
118     max_btag_prob = 1;
119    
120     bin_btag_count = 10;
121     min_btag_count = 0;
122     max_btag_count = 10;
123    
124     //Candidates
125     BTagH_bJet1_csv = new TH1F(("BJet1_CSV"+suffix).c_str(),("BJet1 CSV ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
126     BTagH_bJet2_csv = new TH1F(("BJet2_CSV"+suffix).c_str(),("BJet2 CSV ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
127     BTagH_bJet1_csvmva = new TH1F(("BJet1_CSVMVA"+suffix).c_str(),("BJet1 CSVMVA ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
128     BTagH_bJet2_csvmva = new TH1F(("BJet2_CSVMVA"+suffix).c_str(),("BJet2 CSVMVA ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
129     BTagH_bJet1_ssvhe = new TH1F(("BJet1_SSVHE"+suffix).c_str(),("BJet1 SSVHE ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
130     BTagH_bJet2_ssvhe = new TH1F(("BJet2_SSVHE"+suffix).c_str(),("BJet2 SSVHE ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
131     BTagH_bJet1_jpb = new TH1F(("BJet1_JPB"+suffix).c_str(),("BJet1 JPB ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
132     BTagH_bJet2_jpb = new TH1F(("BJet2_JPB"+suffix).c_str(),("BJet2 JPB ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
133     BTagH_bJet1_tche = new TH1F(("BJet1_TCHE"+suffix).c_str(),("BJet1 TCHE ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
134     BTagH_bJet2_tche = new TH1F(("BJet2_TCHE"+suffix).c_str(),("BJet2 TCHE ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
135     BTagH_bJet1_jp = new TH1F(("BJet1_JP"+suffix).c_str(),("BJet1 JP ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
136     BTagH_bJet2_jp = new TH1F(("BJet2_JP"+suffix).c_str(),("BJet2 JP ("+suffix+")").c_str(), bin_btag_prob, min_btag_prob, max_btag_prob );
137     BTagH_bJet1_tchp = new TH1F(("BJet1_TCHP"+suffix).c_str(),("BJet1 TCHP ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
138     BTagH_bJet2_tchp = new TH1F(("BJet2_TCHP"+suffix).c_str(),("BJet2 TCHP ("+suffix+")").c_str(), bin_btag_count, min_btag_count, max_btag_count );
139    
140     }
141    
142     virtual void fill(VHbbProxy &iProxy,float w) {
143    
144     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
145     const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
146    
147     //Candidates
148     if(iCand->size() > 0){
149     VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
150     VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
151    
152     BTagH_bJet1_csv->Fill(H.jets.at(0).csv, w);
153     BTagH_bJet2_csv->Fill(H.jets.at(1).csv, w);
154     BTagH_bJet1_csvmva->Fill(H.jets.at(0).csvmva, w);
155     BTagH_bJet2_csvmva->Fill(H.jets.at(1).csvmva, w);
156     BTagH_bJet1_ssvhe->Fill(H.jets.at(0).ssvhe, w);
157     BTagH_bJet2_ssvhe->Fill(H.jets.at(1).ssvhe, w);
158     BTagH_bJet1_tche->Fill(H.jets.at(0).tche, w);
159     BTagH_bJet2_tche->Fill(H.jets.at(1).tche, w);
160     BTagH_bJet1_tchp->Fill(H.jets.at(0).tchp, w);
161     BTagH_bJet2_tchp->Fill(H.jets.at(1).tchp, w);
162     BTagH_bJet1_jpb->Fill(H.jets.at(0).jpb, w);
163     BTagH_bJet2_jpb->Fill(H.jets.at(1).jpb, w);
164     BTagH_bJet1_jp->Fill(H.jets.at(0).jp, w);
165     BTagH_bJet2_jp->Fill(H.jets.at(1).jp, w);
166    
167     }
168     }
169    
170     TH1F * BTagH_bJet1_csv;
171     TH1F * BTagH_bJet2_csv;
172     TH1F * BTagH_bJet1_csvmva;
173     TH1F * BTagH_bJet2_csvmva;
174     TH1F * BTagH_bJet1_ssvhe;
175     TH1F * BTagH_bJet2_ssvhe;
176     TH1F * BTagH_bJet1_jpb;
177     TH1F * BTagH_bJet2_jpb;
178     TH1F * BTagH_bJet1_tche;
179     TH1F * BTagH_bJet2_tche;
180     TH1F * BTagH_bJet1_jp;
181     TH1F * BTagH_bJet2_jp;
182     TH1F * BTagH_bJet1_tchp;
183     TH1F * BTagH_bJet2_tchp;
184    
185     private:
186    
187     Int_t bin_btag_prob;
188     Double_t min_btag_prob;
189     Double_t max_btag_prob;
190    
191     Int_t bin_btag_count;
192     Double_t min_btag_count;
193     Double_t max_btag_count;
194    
195     };
196    
197    
198     class CountHisto : public Histos {
199     virtual void book(TFile &f, std::string suffix) {
200    
201     TDirectory *subDir;
202     if( ! f.GetDirectory(suffix.c_str()) )
203     subDir = f.mkdir(suffix.c_str());
204     else
205     subDir = f.GetDirectory(suffix.c_str());
206     subDir->cd();
207     count = new TH1F(("Count"+suffix).c_str(),("Count ("+suffix+")").c_str(), 1,0,2 );
208    
209     }
210     virtual void fill(VHbbProxy &iProxy,float w) {
211     count->Fill(1,w);
212     }
213    
214     TH1F * count;
215     };
216    
217    
218     class StandardHistos : public Histos {
219    
220     public:
221    
222     virtual void book(TFile &f, std::string suffix) {
223    
224     TDirectory *subDir;
225     if( ! f.GetDirectory(suffix.c_str()) )
226     subDir = f.mkdir(suffix.c_str());
227     else
228     subDir = f.GetDirectory(suffix.c_str());
229     subDir->cd();
230    
231     bin_mass = 500;
232     min_mass = 0;
233     max_mass = 300;
234    
235     bin_pt = 500;
236     min_pt = 0;
237     max_pt = 500;
238    
239     bin_hel = 50;
240     min_hel = 0;
241     max_hel = 1;
242    
243     bin_btag = 20;
244     min_btag = 0;
245     max_btag = 1;
246    
247     bin_deltaR = bin_deltaPhi = bin_deltaEta = 20;
248     min_deltaR = min_deltaPhi = min_deltaEta = 0;
249     max_deltaR = max_deltaPhi = max_deltaEta = 5;
250    
251     //Candidates
252     StH_simpleJet1_pt = new TH1F(("SimpleJet1_pt"+suffix).c_str(),("Simple Jet1 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
253     StH_simpleJet2_pt = new TH1F(("SimpleJet2_pt"+suffix).c_str(),("Simple Jet2 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
254     StH_simpleJet1_bTag = new TH1F(("SimpleJet1_bTag"+suffix).c_str(),("Simple Jet1 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
255     StH_simpleJet2_bTag = new TH1F(("SimpleJet2_bTag"+suffix).c_str(),("Simple Jet2 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
256     StH_simpleJets_dR = new TH1F(("SimpleJets_dR"+suffix).c_str(),("Simple Jets deltaR ("+suffix+")").c_str(), bin_deltaR, min_deltaR, max_deltaR );
257     StH_simpleJets_dPhi = new TH1F(("SimpleJets_dPhi"+suffix).c_str(),("Simple Jets deltaPhi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
258     StH_simpleJets_dEta = new TH1F(("SimpleJets_dEta"+suffix).c_str(),("Simple Jets deltaEta ("+suffix+")").c_str(), bin_deltaEta, min_deltaEta, max_deltaEta );
259    
260     StH_HMass = new TH1F(("HiggsMass"+suffix).c_str(),(" Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
261     StH_HPt = new TH1F(("HiggsPt"+suffix).c_str(),(" Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
262     StH_HHel = new TH1F(("HiggsHel"+suffix).c_str(),("Higgs helicity angle ("+suffix+")").c_str(), bin_hel, min_hel, max_hel );
263     StH_HPullAngle = new TH1F(("HiggsPullAngle"+suffix).c_str(),("Higgs pull angle ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
264    
265     StH_ZMass = new TH1F(("ZMass"+suffix).c_str(),(" Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
266     StH_ZPt = new TH1F(("ZPt"+suffix).c_str(),(" Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
267     StH_ZH_dPhi = new TH1F(("ZH_dPhi"+suffix).c_str(),(" ZH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
268    
269     StH_WMass = new TH1F(("WMass"+suffix).c_str(),(" W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
270     StH_WPt = new TH1F(("WPt"+suffix).c_str(),(" W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
271     StH_WH_dPhi = new TH1F(("WH_dPhi"+suffix).c_str(),(" WH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
272    
273     }
274    
275     virtual void fill(VHbbProxy &iProxy,float w) {
276    
277     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
278     const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
279    
280     //Candidates
281     if(iCand->size() > 0){
282     VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
283     VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
284     VHbbCandidate::VectorCandidate V = iCand->at(0).V;
285    
286     StH_simpleJet1_pt->Fill(H.jets.at(0).p4.Pt(), w);
287     StH_simpleJet2_pt->Fill(H.jets.at(1).p4.Pt(), w);
288     StH_simpleJet1_bTag->Fill(H.jets.at(0).csv, w);
289     StH_simpleJet2_bTag->Fill(H.jets.at(1).csv, w);
290     StH_simpleJets_dR->Fill(H.jets.at(0).p4.DeltaR(H.jets.at(1).p4), w);
291     StH_simpleJets_dPhi->Fill(H.jets.at(0).p4.DeltaPhi(H.jets.at(1).p4), w);
292     StH_simpleJets_dEta->Fill(TMath::Abs(H.jets.at(0).p4.Eta()-H.jets.at(1).p4.Eta()), w);
293    
294     StH_HMass->Fill(H.p4.M(), w);
295     StH_HPt->Fill(H.p4.Pt(), w);
296     // StH_HHel->Fill(H.hel(), w);
297     StH_HPullAngle->Fill(H.deltaTheta, w);
298    
299     if( iCandType == VHbbCandidate::Zmumu || iCandType == VHbbCandidate::Zee || iCandType == VHbbCandidate::Znn ){
300     StH_ZMass->Fill(V.p4.M(), w);
301     StH_ZPt->Fill(V.p4.Pt(), w);
302     StH_ZH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
303     }
304     else if(iCandType == VHbbCandidate::Wen || iCandType == VHbbCandidate::Wmun){
305     StH_WMass->Fill(V.p4.M(), w);
306     StH_WPt->Fill(V.p4.Pt(), w);
307     StH_WH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
308     }
309    
310     }
311     }
312    
313     TH1F * StH_simpleJet1_pt;
314     TH1F * StH_simpleJet2_pt;
315     TH1F * StH_simpleJet1_bTag;
316     TH1F * StH_simpleJet2_bTag;
317     TH1F * StH_simpleJets_dR;
318     TH1F * StH_simpleJets_dPhi;
319     TH1F * StH_simpleJets_dEta;
320    
321     TH1F * StH_HMass;
322     TH1F * StH_HPt;
323     TH1F * StH_HHel;
324     TH1F * StH_HPullAngle;
325     TH1F * StH_WMass;
326     TH1F * StH_WPt;
327     TH1F * StH_WH_dPhi;
328     TH1F * StH_ZMass;
329     TH1F * StH_ZPt;
330     TH1F * StH_ZH_dPhi;
331    
332     private:
333    
334     Int_t bin_btag;
335     Double_t min_btag;
336     Double_t max_btag;
337    
338     Int_t bin_deltaEta;
339     Double_t min_deltaEta;
340     Double_t max_deltaEta;
341    
342     Int_t bin_deltaPhi;
343     Double_t min_deltaPhi;
344     Double_t max_deltaPhi;
345    
346     Int_t bin_deltaR;
347     Double_t min_deltaR;
348     Double_t max_deltaR;
349    
350     Int_t bin_mass;
351     Double_t min_mass;
352     Double_t max_mass;
353    
354     Int_t bin_pt;
355     Double_t min_pt;
356     Double_t max_pt;
357    
358     Int_t bin_hel;
359     Double_t min_hel;
360     Double_t max_hel;
361    
362     };
363    
364    
365     class MVAHistos : public Histos {
366    
367     public:
368     MVAHistos( const std::string &channel_):
369     channel(channel_) {}
370    
371     virtual void book(TFile &f, std::string suffix) {
372    
373     // std::clog << "Start booking" << std::endl;
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_BDT = 100;
383     min_BDT = -0.8;
384     max_BDT = 0.8;
385    
386     MVAH_BDT = new TH1F(("BDT"+suffix).c_str(),(" BDT output ("+suffix+")").c_str(), bin_BDT, min_BDT, max_BDT );
387    
388     }
389    
390     virtual void fill(VHbbProxy &iProxy,float w) {
391    
392     // std::clog << "Start filling" << std::endl;
393    
394     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
395     const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
396    
397     TMVA::Reader * reader = new TMVA::Reader();
398    
399     float bbMass,bbPt,btag1,btag2,NaddJet,DeltaRbb,helicity,DeltaPhiVH,bPt1,bPt2,VMass,VPt,pullAngle,DeltaEtabb,deltaPhipfMETjet1,deltaPhipfMETjet2,pfMET,pfMETsig;
400    
401     if(channel == "Zmm"){
402     reader->AddVariable( "bbMass", &bbMass );
403     reader->AddVariable( "VMass", &VMass );
404     reader->AddVariable( "bbPt", &bbPt );
405     reader->AddVariable( "VPt", &VPt );
406     reader->AddVariable( "btag1", &btag1 );
407     reader->AddVariable( "btag2", &btag2 );
408     reader->AddVariable( "DeltaPhiVH", &DeltaPhiVH);
409     reader->AddVariable( "DeltaEtabb", &DeltaEtabb);
410     }
411     else if(channel == "Zee"){
412     reader->AddVariable( "bbMass", &bbMass );
413     reader->AddVariable( "VMass", &VMass );
414     reader->AddVariable( "bbPt", &bbPt );
415     reader->AddVariable( "VPt", &VPt );
416     reader->AddVariable( "btag1", &btag1 );
417     reader->AddVariable( "btag2", &btag2 );
418     reader->AddVariable( "DeltaPhiVH", &DeltaPhiVH);
419     reader->AddVariable( "DeltaEtabb", &DeltaEtabb);
420     }
421     else if(channel == "Znn"){
422     reader->AddVariable( "bbMass", &bbMass );
423     reader->AddVariable( "bbPt", &bbPt );
424     reader->AddVariable( "pfMET", &pfMET );
425     reader->AddVariable( "btag1", &btag1 );
426     reader->AddVariable( "btag2", &btag2 );
427     reader->AddVariable( "DeltaPhiVH", &DeltaPhiVH);
428     }
429     else if(channel == "We"){
430     reader->AddVariable( "bbMass", &bbMass );
431     reader->AddVariable( "bbPt", &bbPt );
432     reader->AddVariable( "VPt", &VPt );
433     reader->AddVariable( "btag1", &btag1 );
434     reader->AddVariable( "btag2", &btag2 );
435     reader->AddVariable( "DeltaPhiVH", &DeltaPhiVH);
436     reader->AddVariable( "DeltaEtabb", &DeltaEtabb);
437     }
438     else if(channel == "Wm"){
439     reader->AddVariable( "bbMass", &bbMass );
440     reader->AddVariable( "bbPt", &bbPt );
441     reader->AddVariable( "VPt", &VPt );
442     reader->AddVariable( "btag1", &btag1 );
443     reader->AddVariable( "btag2", &btag2 );
444     reader->AddVariable( "DeltaPhiVH", &DeltaPhiVH);
445     reader->AddVariable( "DeltaEtabb", &DeltaEtabb);
446     }
447    
448     // --- Book the MVA methods
449    
450     std::vector<std::string> methods;
451     methods.push_back("BDT");
452     TString dir = "weights/";
453     TString prefix = "TMVAClassification";
454     // std::clog << "Booking MVA method..." << std::endl;
455     for (size_t it = 0; it < methods.size(); ++it) {
456     TString methodName = methods.at(it)+ " method";
457     TString weightfile = dir + prefix + "_" + methods.at(it) + "." + channel + TString("_weight.xml");
458     reader->BookMVA( methodName, weightfile );
459     }
460    
461     //Candidates
462     if(iCand->size() > 0 and iCand->at(0).H.jets.size() > 1){
463     VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
464     VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
465     VHbbCandidate::VectorCandidate V = iCand->at(0).V;
466    
467     bbMass = iCand->at(0).H.p4.M();
468     bbPt = iCand->at(0).H.p4.Pt();
469     btag1 = iCand->at(0).H.jets[0].csv;
470     btag2 = iCand->at(0).H.jets[1].csv;
471     NaddJet = iCand->at(0).additionalJets.size();
472     DeltaRbb = deltaR(iCand->at(0).H.jets[0].p4.Eta(),iCand->at(0).H.jets[0].p4.Phi(),iCand->at(0).H.jets[1].p4.Eta(),iCand->at(0).H.jets[1].p4.Phi());
473     helicity = iCand->at(0).H.helicities[0];
474     DeltaPhiVH = Geom::deltaPhi(iCand->at(0).H.p4.Phi(),iCand->at(0).V.p4.Phi()) ;
475     bPt1 = iCand->at(0).H.jets[0].p4.Pt();
476     bPt2 = iCand->at(0).H.jets[1].p4.Pt();
477     VMass = iCand->at(0).V.p4.M();
478     VPt = iCand->at(0).V.p4.Pt();
479     pullAngle = iCand->at(0).H.deltaTheta;
480     DeltaEtabb = TMath::Abs( iCand->at(0).H.jets[0].p4.Eta() - iCand->at(0).H.jets[1].p4.Eta() );
481     deltaPhipfMETjet1 = Geom::deltaPhi( iCand->at(0).V.mets.at(0).p4.Phi(), iCand->at(0).H.jets[0].p4.Phi() );
482     deltaPhipfMETjet2 = Geom::deltaPhi( iCand->at(0).V.mets.at(0).p4.Phi(), iCand->at(0).H.jets[1].p4.Phi() );
483     pfMET = iCand->at(0).V.mets.at(0).p4.Pt();
484     pfMETsig = iCand->at(0).V.mets.at(0).metSig;
485    
486    
487     MVAH_BDT->Fill( reader->EvaluateMVA( "BDT method" ) );
488    
489     }
490     }
491    
492     TH1F * MVAH_BDT;
493    
494     private:
495    
496     const std::string &channel;
497     Int_t bin_BDT;
498     Double_t min_BDT;
499     Double_t max_BDT;
500    
501     };
502    
503    
504     class HardJetHistos : public Histos {
505    
506     public:
507    
508     virtual void book(TFile &f, std::string suffix) {
509    
510     TDirectory *subDir;
511     if( ! f.GetDirectory(suffix.c_str()) )
512     subDir = f.mkdir(suffix.c_str());
513     else
514     subDir = f.GetDirectory(suffix.c_str());
515     subDir->cd();
516    
517     bin_mass = 500;
518     min_mass = 0;
519     max_mass = 300;
520    
521     bin_pt = 500;
522     min_pt = 0;
523     max_pt = 500;
524    
525     bin_hel = 50;
526     min_hel = 0;
527     max_hel = 1;
528    
529     bin_btag = 20;
530     min_btag = 0;
531     max_btag = 1;
532    
533     bin_deltaR = bin_deltaPhi = bin_deltaEta = 20;
534     min_deltaR = min_deltaPhi = min_deltaEta = 0;
535     max_deltaR = max_deltaPhi = max_deltaEta = 5;
536    
537     //Candidates
538    
539     HardJetH_subJet1_pt = new TH1F(("SubJet1_pt"+suffix).c_str(),("Sub Jet1 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
540     HardJetH_subJet2_pt = new TH1F(("SubJet2_pt"+suffix).c_str(),("Sub Jet2 pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
541     HardJetH_subJet1_bTag = new TH1F(("SubJet1_bTag"+suffix).c_str(),("Sub Jet1 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
542     HardJetH_subJet2_bTag = new TH1F(("SubJet2_bTag"+suffix).c_str(),("Sub Jet2 bTag ("+suffix+")").c_str(), bin_btag, min_btag, max_btag );
543     HardJetH_subJets_dR = new TH1F(("SubJets_dR"+suffix).c_str(),("Sub Jets deltaR ("+suffix+")").c_str(), bin_deltaR, min_deltaR, max_deltaR );
544     HardJetH_subJets_dPhi = new TH1F(("SubJets_dPhi"+suffix).c_str(),("Sub Jets deltaPhi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
545     HardJetH_subJets_dEta = new TH1F(("SubJets_dEta"+suffix).c_str(),("Sub Jets deltaEta ("+suffix+")").c_str(), bin_deltaEta, min_deltaEta, max_deltaEta );
546    
547     HardJetH_HMass = new TH1F(("HiggsMass"+suffix).c_str(),(" Higgs Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
548     HardJetH_HPt = new TH1F(("HiggsPt"+suffix).c_str(),(" Higgs Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
549     HardJetH_HHel = new TH1F(("HiggsHel"+suffix).c_str(),("Higgs helicity angle ("+suffix+")").c_str(), bin_hel, min_hel, max_hel );
550    
551     HardJetH_ZMass = new TH1F(("ZMass"+suffix).c_str(),(" Z Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
552     HardJetH_ZPt = new TH1F(("ZPt"+suffix).c_str(),(" Z Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
553     HardJetH_ZH_dPhi = new TH1F(("ZH_dPhi"+suffix).c_str(),(" ZH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
554    
555     HardJetH_WMass = new TH1F(("WMass"+suffix).c_str(),(" W Mass ("+suffix+")").c_str(), bin_mass, min_mass, max_mass );
556     HardJetH_WPt = new TH1F(("WPt"+suffix).c_str(),(" W Pt ("+suffix+")").c_str(), bin_pt, min_pt, max_pt );
557     HardJetH_WH_dPhi = new TH1F(("WH_dPhi"+suffix).c_str(),(" WH delta Phi ("+suffix+")").c_str(), bin_deltaPhi, min_deltaPhi, max_deltaPhi );
558    
559     }
560    
561     virtual void fill(VHbbProxy &iProxy,float w) {
562    
563     const VHbbEvent *iEvent = iProxy.getVHbbEvent();
564     if(iEvent)
565     { const std::vector<VHbbCandidate> *iCand = iProxy.getVHbbCandidate();
566    
567     //Candidates
568     if(iCand->size() > 0){
569     VHbbCandidate::CandidateType iCandType = iCand->at(0).candidateType;
570     VHbbCandidate::HiggsCandidate H = iCand->at(0).H;
571     VHbbCandidate::VectorCandidate V = iCand->at(0).V;
572     std::vector<VHbbEvent::HardJet> iHardJets = iEvent->hardJets;
573     VHbbEvent::HardJet iHardJet = iHardJets.at(0);
574    
575     HardJetH_subJet1_pt->Fill(iHardJet.subFourMomentum.at(0).Pt(), w);
576     HardJetH_subJet2_pt->Fill(iHardJet.subFourMomentum.at(1).Pt(), w);
577     //SubJet information on btag missing
578     // HardJetH_subJet1_bTag->Fill(iHardJet.at(0).csv, w);
579     // HardJetH_subJet2_bTag->Fill(iHardJet.at(1).csv, w);
580     HardJetH_subJets_dR->Fill(iHardJet.subFourMomentum.at(0).DeltaR(iHardJet.subFourMomentum.at(1)), w);
581     HardJetH_subJets_dPhi->Fill(iHardJet.subFourMomentum.at(0).DeltaPhi(iHardJet.subFourMomentum.at(1)), w);
582     HardJetH_subJets_dEta->Fill(TMath::Abs(iHardJet.subFourMomentum.at(0).Eta()-iHardJet.subFourMomentum.at(1).Eta()), w);
583    
584     //Here there should be the higgs candidate from HardJet
585     // HardJetH_HMass->Fill(H.p4.M(), w);
586     // HardJetH_HPt->Fill(H.p4.Pt(), w);
587     // // HardJetH_HHel->Fill(H.hel(), w);
588     // if( iCandType == VHbbCandidate::Zmumu || iCandType == VHbbCandidate::Zee || iCandType == VHbbCandidate::Znn ){
589     // HardJetH_ZMass->Fill(V.p4.M(), w);
590     // HardJetH_ZPt->Fill(V.p4.Pt(), w);
591     // HardJetH_ZH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
592     // }
593     // else if(iCandType == VHbbCandidate::Wen || iCandType == VHbbCandidate::Wmun){
594     // HardJetH_WMass->Fill(V.p4.M(), w);
595     // HardJetH_WPt->Fill(V.p4.Pt(), w);
596     // HardJetH_WH_dPhi->Fill(V.p4.DeltaPhi(H.p4.Phi()), w);
597     // }
598     }
599     }
600     }
601    
602     TH1F * HardJetH_subJet1_pt;
603     TH1F * HardJetH_subJet2_pt;
604     TH1F * HardJetH_subJet1_bTag;
605     TH1F * HardJetH_subJet2_bTag;
606     TH1F * HardJetH_subJets_dR;
607     TH1F * HardJetH_subJets_dPhi;
608     TH1F * HardJetH_subJets_dEta;
609    
610     TH1F * HardJetH_HMass;
611     TH1F * HardJetH_HPt;
612     TH1F * HardJetH_HHel;
613     TH1F * HardJetH_WMass;
614     TH1F * HardJetH_WPt;
615     TH1F * HardJetH_WH_dPhi;
616     TH1F * HardJetH_ZMass;
617     TH1F * HardJetH_ZPt;
618     TH1F * HardJetH_ZH_dPhi;
619    
620     private:
621    
622     Int_t bin_btag;
623     Double_t min_btag;
624     Double_t max_btag;
625    
626     Int_t bin_deltaEta;
627     Double_t min_deltaEta;
628     Double_t max_deltaEta;
629    
630     Int_t bin_deltaPhi;
631     Double_t min_deltaPhi;
632     Double_t max_deltaPhi;
633    
634     Int_t bin_deltaR;
635     Double_t min_deltaR;
636     Double_t max_deltaR;
637    
638     Int_t bin_mass;
639     Double_t min_mass;
640     Double_t max_mass;
641    
642     Int_t bin_pt;
643     Double_t min_pt;
644     Double_t max_pt;
645    
646     Int_t bin_hel;
647     Double_t min_hel;
648     Double_t max_hel;
649    
650     };
651    
652    
653     #endif