ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/SelMods/src/HwwExampleAnalysisMod.cc
Revision: 1.19
Committed: Tue Oct 11 13:43:24 2011 UTC (13 years, 6 months ago) by ceballos
Content type: text/plain
Branch: MAIN
Changes since 1.18: +9 -3 lines
Log Message:
new

File Contents

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