ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/SelMods/src/HwwExampleAnalysisMod.cc
Revision: 1.6
Committed: Sat Oct 23 04:48:07 2010 UTC (14 years, 6 months ago) by ceballos
Content type: text/plain
Branch: MAIN
Changes since 1.5: +131 -21 lines
Log Message:
new cuts

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