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