ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/SelMods/src/HwwExampleAnalysisMod.cc
Revision: 1.14
Committed: Tue Apr 12 07:38:45 2011 UTC (14 years ago) by ceballos
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_020d, TMit_020d, Mit_020c, Mit_021pre1, Mit_020b
Changes since 1.13: +52 -76 lines
Log Message:
minor update

File Contents

# User Rev Content
1 ceballos 1.1 // $Id $
2    
3     #include "MitPhysics/SelMods/interface/HwwExampleAnalysisMod.h"
4     #include <TH1D.h>
5     #include <TH2D.h>
6     #include <TParameter.h>
7     #include "MitAna/DataUtil/interface/Debug.h"
8     #include "MitAna/DataTree/interface/Names.h"
9     #include "MitPhysics/Init/interface/ModNames.h"
10     #include "MitAna/DataCont/interface/ObjArray.h"
11     #include "MitCommon/MathTools/interface/MathUtils.h"
12 ceballos 1.14 #include "MitPhysics/Utils/interface/MetTools.h"
13 ceballos 1.1 #include "MitAna/DataTree/interface/ParticleCol.h"
14     #include "TFile.h"
15     #include "TTree.h"
16    
17     using namespace mithep;
18     ClassImp(mithep::HwwExampleAnalysisMod)
19    
20     //--------------------------------------------------------------------------------------------------
21     HwwExampleAnalysisMod::HwwExampleAnalysisMod(const char *name, const char *title) :
22     BaseMod(name,title),
23     fMuonBranchName(Names::gkMuonBrn),
24     fMetName("NoDefaultNameSet"),
25     fCleanJetsName("NoDefaultNameSet"),
26 ceballos 1.6 fCleanJetsNoPtCutName("NoDefaultNameSet"),
27 ceballos 1.5 fVertexName(ModNames::gkGoodVertexesName),
28 ceballos 1.14 fPFCandidatesName(Names::gkPFCandidatesBrn),
29 ceballos 1.1 fMuons(0),
30     fMet(0),
31 ceballos 1.6 fVertices(0),
32 ceballos 1.14 fPFCandidates(0),
33 ceballos 1.13 fNEventsSelected(0)
34 ceballos 1.1 {
35     // Constructor.
36     }
37    
38     //--------------------------------------------------------------------------------------------------
39     void HwwExampleAnalysisMod::Begin()
40     {
41     // Run startup code on the client machine. For this module, we dont do
42     // anything here.
43     }
44    
45     //--------------------------------------------------------------------------------------------------
46     void HwwExampleAnalysisMod::SlaveBegin()
47     {
48     // Run startup code on the computer (slave) doing the actual analysis. Here,
49     // we typically initialize histograms and other analysis objects and request
50     // branches. For this module, we request a branch of the MitTree.
51    
52     // Load Branches
53 ceballos 1.14 ReqBranch(fMuonBranchName, fMuons);
54     ReqBranch(fPFCandidatesName, fPFCandidates);
55 ceballos 1.1
56     //Create your histograms here
57    
58     //*************************************************************************************************
59     // Selection Histograms
60     //*************************************************************************************************
61 ceballos 1.7 AddTH1(fHWWSelection,"hHWWSelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
62     AddTH1(fHWWToEESelection,"hHWWToEESelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
63     AddTH1(fHWWToMuMuSelection,"hHWWToMuMuSelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
64     AddTH1(fHWWToEMuSelection,"hHWWToEMuSelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
65 ceballos 1.1
66     //***********************************************************************************************
67     // Histograms after preselection
68     //***********************************************************************************************
69     AddTH1(fLeptonEta ,"hLeptonEta",";LeptonEta;Number of Events",100,-5.,5.0);
70     AddTH1(fLeptonPtMax ,"hLeptonPtMax",";Lepton P_t Max;Number of Events",150,0.,150.);
71     AddTH1(fLeptonPtMin ,"hLeptonPtMin",";Lepton P_t Min;Number of Events",150,0.,150.);
72     AddTH1(fMetPtHist ,"hMetPtHist",";Met;Number of Events",150,0.,300.);
73     AddTH1(fMetPhiHist ,"hMetPhiHist",";#phi;Number of Events",28,-3.5,3.5);
74     AddTH1(fUncorrMetPtHist ,"hUncorrMetPtHist",";Met;Number of Events",150,0.,300.);
75     AddTH1(fUncorrMetPhiHist ,"hUncorrMetPhiHist",";#phi;Number of Events",28,-3.5,3.5);
76     AddTH1(fDeltaPhiLeptons ,"hDeltaPhiLeptons",";#Delta#phi_{ll};Number of Events",90,0,180);
77 ceballos 1.8 AddTH1(fDeltaEtaLeptons ,"hDeltaEtaLeptons",";#Delta#eta_{ll};Number of Events",100,-5.,5.0);
78 ceballos 1.1 AddTH1(fDileptonMass ,"hDileptonMass",";Mass_{ll};Number of Events",150,0.,300.);
79    
80     //***********************************************************************************************
81     // N-1 Histograms
82     //***********************************************************************************************
83     //All events
84     AddTH1(fLeptonPtMax_NMinusOne ,"hLeptonPtMax_NMinusOne",
85     ";Lepton P_t Max;Number of Events",150,0.,150.);
86     AddTH1(fLeptonPtMin_NMinusOne ,"hLeptonPtMin_NMinusOne",
87     ";Lepton P_t Min;Number of Events",150,0.,150.);
88     AddTH1(fMetPtHist_NMinusOne ,"hMetPtHist_NMinusOne",
89     ";Met;Number of Events",150,0.,300.);
90     AddTH1(fMetPhiHist_NMinusOne ,"hMetPhiHist_NMinusOne",
91     ";#phi;Number of Events",28,-3.5,3.5);
92     AddTH1(fMETdeltaPhilEtHist_NMinusOne ,"hMETdeltaPhilEtHist_NMinusOne",
93     ";METdeltaPhilEtHist;Number of Events",150,0.,300.);
94     AddTH1(fNCentralJets_NMinusOne ,"hNCentralJets_NMinusOne",
95     ";Number of Central Jets;Number of Events",6,-0.5,5.5);
96     AddTH1(fNSoftMuonsHist_NMinusOne ,"hNSoftMuonsHist_NMinusOne",
97     ";Number of Dirty Muons; Number of Events",6,-0.5,5.5);
98     AddTH1(fDeltaPhiLeptons_NMinusOne ,"hDeltaPhiLeptons_NMinusOne",
99     ";#Delta#phi_{ll};Number of Events",90,0,180);
100     AddTH1(fDeltaEtaLeptons_NMinusOne ,"hDeltaEtaLeptons_NMinusOne",
101     ";#Delta#eta_{ll};Number of Events",100,-5.0,5.0);
102     AddTH1(fDileptonMass_NMinusOne ,"hDileptonMass_NMinusOne",
103     ";Mass_{ll};Number of Events",150,0.,300.);
104     AddTH1(fMinDeltaPhiLeptonMet_NMinusOne ,"hMinDeltaPhiLeptonMet_NMinusOne",
105     ";Min #Delta#phi_{l,Met};Number of Events",90,0.,180);
106    
107     //***********************************************************************************************
108     // After all cuts Histograms
109     //***********************************************************************************************
110     AddTH1(fMinDeltaPhiLeptonMet_afterCuts ,"hMinDeltaPhiLeptonMet_afterCuts",
111     ";Min #Delta#phi_{l,Met};Number of Events",90,0.,180);
112     AddTH1(fMtLepton1_afterCuts ,"hMtLepton1_afterCuts",
113     ";M_t (Lepton1,Met);Number of Events",100,0.,200.);
114     AddTH1(fMtLepton2_afterCuts ,"hMtLepton2_afterCuts",
115     ";M_t (Lepton2,Met);Number of Events",100,0.,200.);
116     AddTH1(fMtHiggs_afterCuts ,"hMtHiggs_afterCuts",
117     ";M_t (l1+l2+Met);Number of Events",150,0.,300.);
118     AddTH1(fLeptonPtPlusMet_afterCuts ,"hLeptonPtPlusMet_afterCuts",
119     ";LeptonPtPlusMet;Number of Events",150,0., 300.);
120    
121     }
122    
123     //--------------------------------------------------------------------------------------------------
124     void HwwExampleAnalysisMod::Process()
125     {
126     // Process entries of the tree. For this module, we just load the branches and
127     LoadBranch(fMuonBranchName);
128 ceballos 1.14 LoadBranch(fPFCandidatesName);
129    
130 ceballos 1.1 //Obtain all the good objects from the event cleaning module
131 ceballos 1.5 fVertices = GetObjThisEvt<VertexOArr>(fVertexName);
132 ceballos 1.1 ObjArray<Muon> *CleanMuons = dynamic_cast<ObjArray<Muon>* >(FindObjThisEvt(ModNames::gkCleanMuonsName));
133 ceballos 1.6 ObjArray<Electron> *CleanElectrons = dynamic_cast<ObjArray<Electron>* >(FindObjThisEvt(ModNames::gkCleanElectronsName));
134 ceballos 1.1 ParticleOArr *CleanLeptons = dynamic_cast<mithep::ParticleOArr*>
135     (FindObjThisEvt(ModNames::gkMergedLeptonsName));
136     ObjArray<Jet> *CleanJets = dynamic_cast<ObjArray<Jet>* >
137     (FindObjThisEvt(fCleanJetsName.Data()));
138 ceballos 1.6 ObjArray<Jet> *CleanJetsNoPtCut = dynamic_cast<ObjArray<Jet>* >
139     (FindObjThisEvt(fCleanJetsNoPtCutName.Data()));
140 ceballos 1.2 TParameter<Double_t> *NNLOWeight = GetObjThisEvt<TParameter<Double_t> >("NNLOWeight");
141 ceballos 1.1
142     MetCol *met = dynamic_cast<ObjArray<Met>* >(FindObjThisEvt(fMetName));
143 ceballos 1.14 const Met *stdMet = 0;
144 ceballos 1.1 if (met) {
145 ceballos 1.14 stdMet = met->At(0);
146 ceballos 1.1 } else {
147     cout << "Error: Met Collection " << fMetName << " could not be loaded.\n";
148     return;
149     }
150    
151     //***********************************************************************************************
152     //Kinematic PreSelection
153     //***********************************************************************************************
154     // At least two leptons in the event
155     if (CleanLeptons->GetEntries() < 2) return;
156     // Pt1 > 20 && Pt2 > 10
157 ceballos 1.4 if(CleanLeptons->At(0)->Pt() <= 20 || CleanLeptons->At(1)->Pt() <= 10) return;
158 ceballos 1.1 // opposite charge leptons
159     if(CleanLeptons->At(0)->Charge() * CleanLeptons->At(1)->Charge() > 0) return;
160    
161     CompositeParticle *dilepton = new CompositeParticle();
162     dilepton->AddDaughter(CleanLeptons->At(0));
163     dilepton->AddDaughter(CleanLeptons->At(1));
164    
165     //***********************************************************************************************
166     //Get Dirty Muons: Non-isolated Muons (exclude the clean muons)
167     //***********************************************************************************************
168     ObjArray<Muon> *SoftMuons = new ObjArray<Muon>;
169     for (UInt_t i=0; i<fMuons->GetEntries(); ++i) {
170     const Muon *mu = fMuons->At(i);
171     if(!MuonTools::PassSoftMuonCut(mu, fVertices)) continue;
172    
173     bool isCleanMuon = kFALSE;
174     for (UInt_t j=0; j<CleanMuons->GetEntries(); j++) {
175     if(fMuons->At(i) == CleanMuons->At(j) &&
176     CleanMuons->At(j)->Pt() > 10) isCleanMuon = kTRUE;
177     }
178     if(isCleanMuon == kFALSE) SoftMuons->Add(mu);
179     }
180    
181     //***********************************************************************************************
182 ceballos 1.6 //|Z_vert-Z_l| maximum
183     //***********************************************************************************************
184 ceballos 1.14 std::vector<double> leptonsDz;
185 ceballos 1.6 double zDiffMax = 0.0;
186     if(fVertices->GetEntries() > 0) {
187     for (UInt_t j=0; j<CleanMuons->GetEntries(); j++) {
188 ceballos 1.14 double pDz = CleanMuons->At(j)->BestTrk()->DzCorrected(*fVertices->At(0));
189     leptonsDz.push_back(pDz);
190 ceballos 1.6 }
191     for (UInt_t j=0; j<CleanElectrons->GetEntries(); j++) {
192 ceballos 1.14 double pDz = CleanElectrons->At(j)->GsfTrk()->DzCorrected(*fVertices->At(0));
193     leptonsDz.push_back(pDz);
194     }
195     for(UInt_t t=0; t<leptonsDz.size(); t++) {
196     for(UInt_t i=t+1; i<leptonsDz.size(); i++) {
197     if(TMath::Abs(leptonsDz[t]-leptonsDz[i]) > zDiffMax) zDiffMax = TMath::Abs(leptonsDz[t]-leptonsDz[i]);
198 ceballos 1.7 }
199 ceballos 1.6 }
200 ceballos 1.14 leptonsDz.clear();
201 ceballos 1.6 }
202    
203     //***********************************************************************************************
204 ceballos 1.1 //Define Event Variables
205     //***********************************************************************************************
206     //delta phi between the 2 leptons in degrees
207     double deltaPhiLeptons = MathUtils::DeltaPhi(CleanLeptons->At(0)->Phi(),
208     CleanLeptons->At(1)->Phi())* 180.0 / TMath::Pi();
209    
210 ceballos 1.8 double deltaEtaLeptons = CleanLeptons->At(0)->Eta() - CleanLeptons->At(1)->Eta();
211 ceballos 1.1
212 ceballos 1.14 double deltaPhiDileptonMet = MathUtils::DeltaPhi(stdMet->Phi(),
213 ceballos 1.1 dilepton->Phi())*180.0 / TMath::Pi();
214    
215 ceballos 1.14 double mtHiggs = TMath::Sqrt(2.0*dilepton->Pt() * stdMet->Pt()*
216 ceballos 1.1 (1.0 - cos(deltaPhiDileptonMet * TMath::Pi() / 180.0)));
217    
218     //angle between MET and closest lepton
219 ceballos 1.14 double deltaPhiMetLepton[2] = {MathUtils::DeltaPhi(stdMet->Phi(), CleanLeptons->At(0)->Phi()),
220     MathUtils::DeltaPhi(stdMet->Phi(), CleanLeptons->At(1)->Phi())};
221 ceballos 1.1
222 ceballos 1.14 double mTW[2] = {TMath::Sqrt(2.0*CleanLeptons->At(0)->Pt()*stdMet->Pt()*
223 ceballos 1.1 (1.0 - cos(deltaPhiMetLepton[0]))),
224 ceballos 1.14 TMath::Sqrt(2.0*CleanLeptons->At(1)->Pt()*stdMet->Pt()*
225 ceballos 1.1 (1.0 - cos(deltaPhiMetLepton[1])))};
226    
227     double minDeltaPhiMetLepton = (deltaPhiMetLepton[0] < deltaPhiMetLepton[1])?
228     deltaPhiMetLepton[0]:deltaPhiMetLepton[1];
229    
230 ceballos 1.14 MetTools metTools(CleanMuons, CleanElectrons, fPFCandidates, fVertices->At(0), 0.2, 8.0, 5.0);
231     double pMET[2] = {metTools.GetProjectedMet(CleanLeptons,stdMet),
232     metTools.GetProjectedTrackMet(CleanLeptons)};
233    
234     double METdeltaPhilEt = TMath::Min(pMET[0],pMET[1]);
235 ceballos 1.1
236 ceballos 1.6 //count the number of central Jets for vetoing and b-tagging
237     vector<Jet*> sortedJetsAll;
238     vector<Jet*> sortedJets;
239     vector<Jet*> sortedJetsLowPt;
240     for(UInt_t i=0; i<CleanJetsNoPtCut->GetEntries(); i++){
241     Jet* jet_a = new Jet(CleanJetsNoPtCut->At(i)->Px(),
242     CleanJetsNoPtCut->At(i)->Py(),
243     CleanJetsNoPtCut->At(i)->Pz(),
244     CleanJetsNoPtCut->At(i)->E() );
245 ceballos 1.14 jet_a->SetMatchedMCFlavor(CleanJetsNoPtCut->At(i)->MatchedMCFlavor());
246     jet_a->SetCombinedSecondaryVertexBJetTagsDisc(CleanJetsNoPtCut->At(i)->CombinedSecondaryVertexBJetTagsDisc());
247     jet_a->SetCombinedSecondaryVertexMVABJetTagsDisc(CleanJetsNoPtCut->At(i)->CombinedSecondaryVertexMVABJetTagsDisc());
248     jet_a->SetJetProbabilityBJetTagsDisc(CleanJetsNoPtCut->At(i)->JetProbabilityBJetTagsDisc());
249     jet_a->SetJetBProbabilityBJetTagsDisc(CleanJetsNoPtCut->At(i)->JetBProbabilityBJetTagsDisc());
250     jet_a->SetTrackCountingHighEffBJetTagsDisc(CleanJetsNoPtCut->At(i)->TrackCountingHighEffBJetTagsDisc());
251     jet_a->SetTrackCountingHighPurBJetTagsDisc(CleanJetsNoPtCut->At(i)->TrackCountingHighPurBJetTagsDisc());
252     jet_a->SetSimpleSecondaryVertexBJetTagsDisc(CleanJetsNoPtCut->At(i)->SimpleSecondaryVertexBJetTagsDisc());
253     jet_a->SetSimpleSecondaryVertexHighEffBJetTagsDisc(CleanJetsNoPtCut->At(i)->SimpleSecondaryVertexHighEffBJetTagsDisc());
254     jet_a->SetSimpleSecondaryVertexHighPurBJetTagsDisc(CleanJetsNoPtCut->At(i)->SimpleSecondaryVertexHighPurBJetTagsDisc());
255 ceballos 1.6 sortedJetsAll.push_back(jet_a);
256     }
257    
258 ceballos 1.1 for(UInt_t i=0; i<CleanJets->GetEntries(); i++){
259 ceballos 1.9 if(TMath::Abs(CleanJets->At(i)->Eta()) < 5.0 &&
260 ceballos 1.14 CleanJets->At(i)->Pt() > 30.0){
261 ceballos 1.6 Jet* jet_b = new Jet(CleanJets->At(i)->Px(),
262     CleanJets->At(i)->Py(),
263     CleanJets->At(i)->Pz(),
264     CleanJets->At(i)->E() );
265     sortedJets.push_back(jet_b);
266     }
267     }
268    
269     for(UInt_t i=0; i<sortedJetsAll.size(); i++){
270     bool overlap = kFALSE;
271     for(UInt_t j=0; j<sortedJets.size(); j++){
272     if(sortedJetsAll[i]->Pt() == sortedJets[j]->Pt() ||
273     (sortedJetsAll[i]->CombinedSecondaryVertexBJetTagsDisc() == sortedJets[j]->CombinedSecondaryVertexBJetTagsDisc() &&
274     sortedJetsAll[i]->JetBProbabilityBJetTagsDisc() == sortedJets[j]->JetBProbabilityBJetTagsDisc() &&
275     sortedJetsAll[i]->TrackCountingHighPurBJetTagsDisc() == sortedJets[j]->TrackCountingHighPurBJetTagsDisc())
276     ) {
277     sortedJets[j]->SetMatchedMCFlavor(sortedJetsAll[i]->MatchedMCFlavor());
278     sortedJets[j]->SetCombinedSecondaryVertexBJetTagsDisc(sortedJetsAll[i]->CombinedSecondaryVertexBJetTagsDisc());
279     sortedJets[j]->SetCombinedSecondaryVertexMVABJetTagsDisc(sortedJetsAll[i]->CombinedSecondaryVertexMVABJetTagsDisc());
280     sortedJets[j]->SetJetProbabilityBJetTagsDisc(sortedJetsAll[i]->JetProbabilityBJetTagsDisc());
281     sortedJets[j]->SetJetBProbabilityBJetTagsDisc(sortedJetsAll[i]->JetBProbabilityBJetTagsDisc());
282     sortedJets[j]->SetTrackCountingHighEffBJetTagsDisc(sortedJetsAll[i]->TrackCountingHighEffBJetTagsDisc());
283     sortedJets[j]->SetTrackCountingHighPurBJetTagsDisc(sortedJetsAll[i]->TrackCountingHighPurBJetTagsDisc());
284     sortedJets[j]->SetSimpleSecondaryVertexBJetTagsDisc(sortedJetsAll[i]->SimpleSecondaryVertexBJetTagsDisc());
285     sortedJets[j]->SetSimpleSecondaryVertexHighEffBJetTagsDisc(sortedJetsAll[i]->SimpleSecondaryVertexHighEffBJetTagsDisc());
286     sortedJets[j]->SetSimpleSecondaryVertexHighPurBJetTagsDisc(sortedJetsAll[i]->SimpleSecondaryVertexHighPurBJetTagsDisc());
287     overlap = kTRUE;
288     break;
289     }
290     }
291     if(overlap == kFALSE){
292     sortedJetsLowPt.push_back(sortedJetsAll[i]);
293     }
294     }
295     double maxBtag = -99999.;
296     double imaxBtag = -1;
297     for(UInt_t i=0; i<sortedJetsLowPt.size(); i++){
298     if(sortedJetsLowPt[i]->TrackCountingHighEffBJetTagsDisc() > maxBtag){
299     maxBtag = sortedJetsLowPt[i]->TrackCountingHighEffBJetTagsDisc();
300     imaxBtag = i;
301 ceballos 1.1 }
302     }
303    
304     //Lepton Type
305     int finalstateType = -1;
306     if (CleanLeptons->At(0)->ObjType() == kMuon && CleanLeptons->At(1)->ObjType() == kMuon ){ // mumu
307     finalstateType = 10;
308     } else if(CleanLeptons->At(0)->ObjType() == kElectron && CleanLeptons->At(1)->ObjType() == kElectron ){ // ee
309     finalstateType = 11;
310     } else if((CleanLeptons->At(0)->ObjType() == kElectron && CleanLeptons->At(1)->ObjType() == kMuon) ||
311     (CleanLeptons->At(1)->ObjType() == kElectron && CleanLeptons->At(0)->ObjType() == kMuon)) {
312     finalstateType = 12;
313     } else {
314     cerr << "Error: finalstate lepton type not supported\n";
315     }
316    
317     //*********************************************************************************************
318     //Define Cuts
319     //*********************************************************************************************
320 ceballos 1.7 const int nCuts = 10;
321     bool passCut[nCuts] = {false, false, false, false, false,
322     false, false, false, false, false};
323 ceballos 1.1
324 ceballos 1.14 Bool_t PreselPtCut = kTRUE;
325     if(CleanLeptons->At(0)->Pt() <= 20) PreselPtCut = kFALSE;
326     if(CleanLeptons->At(1)->Pt() <= 10) PreselPtCut = kFALSE;
327     if(CleanLeptons->At(1)->ObjType() == kElectron && CleanLeptons->At(1)->Pt() <= 15) PreselPtCut = kFALSE;
328     if(PreselPtCut == kTRUE) passCut[0] = true;
329 ceballos 1.1
330 ceballos 1.14 if(zDiffMax < 100000.0) passCut[1] = true;
331 ceballos 1.1
332 ceballos 1.14 if(stdMet->Pt() > 20.0) passCut[2] = true;
333 ceballos 1.1
334 ceballos 1.7 if(dilepton->Mass() > 12.0) passCut[3] = true;
335    
336     if(sortedJets.size() < 1) passCut[6] = true;
337    
338     if(SoftMuons->GetEntries() == 0) passCut[7] = true;
339 ceballos 1.6
340 ceballos 1.7 if(CleanLeptons->GetEntries() == 2) passCut[8] = true;
341 ceballos 1.1
342 ceballos 1.7 if(maxBtag < 2.1) passCut[9] = true;
343 ceballos 1.1
344     if (finalstateType == 10 || finalstateType == 11){ // mumu/ee
345 ceballos 1.7 if(fabs(dilepton->Mass()-91.1876) > 15.0) passCut[4] = true;
346     if(METdeltaPhilEt > 35) passCut[5] = true;
347 ceballos 1.1 }
348     else if(finalstateType == 12) { // emu
349 ceballos 1.7 passCut[4] = true;
350     if(METdeltaPhilEt > 20) passCut[5] = true;
351 ceballos 1.1 }
352 ceballos 1.13
353 ceballos 1.1 //*********************************************************************************************
354     //Make Selection Histograms. Number of events passing each level of cut
355     //*********************************************************************************************
356     bool passAllCuts = true;
357     for(int c=0; c<nCuts; c++) passAllCuts = passAllCuts & passCut[c];
358 ceballos 1.13 if(passAllCuts) fNEventsSelected++;
359    
360 ceballos 1.1 //Cut Selection Histograms
361 ceballos 1.2 fHWWSelection->Fill(-1,NNLOWeight->GetVal());
362 ceballos 1.1 if (finalstateType == 10 )
363 ceballos 1.2 fHWWToMuMuSelection->Fill(-1,NNLOWeight->GetVal());
364 ceballos 1.1 else if(finalstateType == 11 )
365 ceballos 1.2 fHWWToEESelection->Fill(-1,NNLOWeight->GetVal());
366 ceballos 1.1 else if(finalstateType == 12 )
367 ceballos 1.2 fHWWToEMuSelection->Fill(-1,NNLOWeight->GetVal());
368 ceballos 1.1
369     for (int k=0;k<nCuts;k++) {
370     bool pass = true;
371     bool passPreviousCut = true;
372     for (int p=0;p<=k;p++) {
373     pass = (pass && passCut[p]);
374     if (p<k)
375     passPreviousCut = (passPreviousCut&& passCut[p]);
376     }
377    
378     if (pass) {
379 ceballos 1.2 fHWWSelection->Fill(k,NNLOWeight->GetVal());
380 ceballos 1.1 if (finalstateType == 10 )
381 ceballos 1.2 fHWWToMuMuSelection->Fill(k,NNLOWeight->GetVal());
382 ceballos 1.1 else if(finalstateType == 11)
383 ceballos 1.2 fHWWToEESelection->Fill(k,NNLOWeight->GetVal());
384 ceballos 1.1 else if(finalstateType == 12)
385 ceballos 1.2 fHWWToEMuSelection->Fill(k,NNLOWeight->GetVal());
386 ceballos 1.1 }
387     }
388    
389     //*****************************************************************************************
390     //Make Preselection Histograms
391     //*****************************************************************************************
392 ceballos 1.2 fLeptonEta->Fill(CleanLeptons->At(0)->Eta(),NNLOWeight->GetVal());
393     fLeptonEta->Fill(CleanLeptons->At(1)->Eta(),NNLOWeight->GetVal());
394     fLeptonPtMax->Fill(CleanLeptons->At(0)->Pt(),NNLOWeight->GetVal());
395     fLeptonPtMin->Fill(CleanLeptons->At(1)->Pt(),NNLOWeight->GetVal());
396 ceballos 1.14 fMetPtHist->Fill(stdMet->Pt(),NNLOWeight->GetVal());
397     fMetPhiHist->Fill(stdMet->Phi(),NNLOWeight->GetVal());
398 ceballos 1.2 fDeltaPhiLeptons->Fill(deltaPhiLeptons,NNLOWeight->GetVal());
399     fDeltaEtaLeptons->Fill(deltaEtaLeptons,NNLOWeight->GetVal());
400     fDileptonMass->Fill(dilepton->Mass(),NNLOWeight->GetVal());
401 ceballos 1.1
402     //*********************************************************************************************
403     // N-1 Histograms
404     //*********************************************************************************************
405     bool pass;;
406    
407     //N Jet Veto
408     pass = true;
409     for (int k=0;k<nCuts;k++) {
410 ceballos 1.11 if (k != 6) {
411 ceballos 1.1 pass = (pass && passCut[k]);
412     }
413     }
414     if (pass) {
415 ceballos 1.6 fNCentralJets_NMinusOne->Fill(sortedJets.size(),NNLOWeight->GetVal());
416 ceballos 1.1 }
417    
418     // Final Met Cut
419     pass = true;
420     for (int k=0;k<nCuts;k++) {
421 ceballos 1.11 if (k != 5) {
422 ceballos 1.1 pass = (pass && passCut[k]);
423     }
424     }
425     if (pass) {
426 ceballos 1.14 fMetPtHist_NMinusOne->Fill(stdMet->Pt(),NNLOWeight->GetVal());
427 ceballos 1.1 }
428    
429     // dilepton mass
430     pass = true;
431     for (int k=0;k<nCuts;k++) {
432 ceballos 1.11 if (k != 3 && k != 4)
433 ceballos 1.1 pass = (pass && passCut[k]);
434     }
435     if (pass) {
436 ceballos 1.2 fDileptonMass_NMinusOne->Fill(dilepton->Mass(),NNLOWeight->GetVal());
437 ceballos 1.1 }
438    
439     // Lepton Pt Max, Lepton Pt Min, DeltaPhiLeptons
440     pass = true;
441     for (int k=0;k<nCuts;k++) {
442 ceballos 1.12 if (k != 0)
443     pass = (pass && passCut[k]);
444 ceballos 1.1 }
445     if (pass) {
446 ceballos 1.2 fLeptonPtMax_NMinusOne->Fill(CleanLeptons->At(0)->Pt(),NNLOWeight->GetVal());
447     fLeptonPtMin_NMinusOne->Fill(CleanLeptons->At(1)->Pt(),NNLOWeight->GetVal());
448     fDeltaPhiLeptons_NMinusOne->Fill(deltaPhiLeptons,NNLOWeight->GetVal());
449 ceballos 1.1 }
450    
451     // NSoftMuons
452     pass = true;
453     for (int k=0;k<nCuts;k++) {
454 ceballos 1.11 if (k != 7)
455 ceballos 1.1 pass = (pass && passCut[k]);
456     }
457     if (pass) {
458 ceballos 1.2 fNSoftMuonsHist_NMinusOne->Fill(SoftMuons->GetEntries(),NNLOWeight->GetVal());
459 ceballos 1.1 }
460    
461     //*********************************************************************************************
462     //Plots after all Cuts
463     //*********************************************************************************************
464     if (passAllCuts) {
465 ceballos 1.2 fMinDeltaPhiLeptonMet_afterCuts->Fill(minDeltaPhiMetLepton,NNLOWeight->GetVal());
466     fMtLepton1_afterCuts->Fill(mTW[0],NNLOWeight->GetVal());
467     fMtLepton2_afterCuts->Fill(mTW[1],NNLOWeight->GetVal());
468     fMtHiggs_afterCuts->Fill(mtHiggs,NNLOWeight->GetVal());
469 ceballos 1.14 fLeptonPtPlusMet_afterCuts->Fill(CleanLeptons->At(0)->Pt()+CleanLeptons->At(1)->Pt()+stdMet->Pt(),NNLOWeight->GetVal());
470 ceballos 1.1 }
471    
472     delete dilepton;
473     delete SoftMuons;
474 ceballos 1.6 for(UInt_t i=0; i<sortedJets.size(); i++) delete sortedJets[i];
475     for(UInt_t i=0; i<sortedJetsAll.size(); i++) delete sortedJetsAll[i];
476 ceballos 1.1 return;
477     }
478    
479     //--------------------------------------------------------------------------------------------------
480     void HwwExampleAnalysisMod::SlaveTerminate()
481     {
482    
483     // Run finishing code on the computer (slave) that did the analysis. For this
484     // module, we dont do anything here.
485 ceballos 1.13 cout << "selected events on HwwExampleAnalysisMod: " << fNEventsSelected << endl;
486 ceballos 1.1
487     }
488     //--------------------------------------------------------------------------------------------------
489     void HwwExampleAnalysisMod::Terminate()
490     {
491     // Run finishing code on the client computer. For this module, we dont do
492     // anything here.
493    
494     }