ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/SelMods/src/HwwExampleAnalysisMod.cc
Revision: 1.7
Committed: Tue Oct 26 09:58:30 2010 UTC (14 years, 6 months ago) by ceballos
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_016
Changes since 1.6: +37 -20 lines
Log Message:
new v1 ww selection

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 {
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 ReqBranch(fMuonBranchName, fMuons);
52 ReqBranch(fCaloJetName0, fCaloJet0);
53
54 //Create your histograms here
55
56 //*************************************************************************************************
57 // Selection Histograms
58 //*************************************************************************************************
59 AddTH1(fHWWSelection,"hHWWSelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
60 AddTH1(fHWWToEESelection,"hHWWToEESelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
61 AddTH1(fHWWToMuMuSelection,"hHWWToMuMuSelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
62 AddTH1(fHWWToEMuSelection,"hHWWToEMuSelection", ";Cut Number;Number of Events", 11, -1.5, 9.5);
63
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 LoadBranch(fCaloJetName0);
127
128 //Obtain all the good objects from the event cleaning module
129 fVertices = GetObjThisEvt<VertexOArr>(fVertexName);
130 ObjArray<Muon> *CleanMuons = dynamic_cast<ObjArray<Muon>* >(FindObjThisEvt(ModNames::gkCleanMuonsName));
131 ObjArray<Electron> *CleanElectrons = dynamic_cast<ObjArray<Electron>* >(FindObjThisEvt(ModNames::gkCleanElectronsName));
132 ParticleOArr *CleanLeptons = dynamic_cast<mithep::ParticleOArr*>
133 (FindObjThisEvt(ModNames::gkMergedLeptonsName));
134 ObjArray<Jet> *CleanJets = dynamic_cast<ObjArray<Jet>* >
135 (FindObjThisEvt(fCleanJetsName.Data()));
136 ObjArray<Jet> *CleanJetsNoPtCut = dynamic_cast<ObjArray<Jet>* >
137 (FindObjThisEvt(fCleanJetsNoPtCutName.Data()));
138 TParameter<Double_t> *NNLOWeight = GetObjThisEvt<TParameter<Double_t> >("NNLOWeight");
139
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 if(CleanLeptons->At(0)->Pt() <= 20 || CleanLeptons->At(1)->Pt() <= 10) return;
156 // 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 //|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 double pDz = 0.0;
186 for(uint i0 = 0; i0 < fVertices->GetEntries(); i0++) {
187 if(fVertices->At(i0)->NTracks() > 0){
188 pDz = TMath::Abs(CleanMuons->At(j)->BestTrk()->DzCorrected(*fVertices->At(i0)));
189 break;
190 }
191 }
192 if(pDz > zDiffMax) pDz = zDiffMax;
193 }
194 for (UInt_t j=0; j<CleanElectrons->GetEntries(); j++) {
195 double pDz = 0.0;
196 for(uint i0 = 0; i0 < fVertices->GetEntries(); i0++) {
197 if(fVertices->At(i0)->NTracks() > 0){
198 pDz = TMath::Abs(CleanElectrons->At(j)->GsfTrk()->DzCorrected(*fVertices->At(i0)));
199 break;
200 }
201 }
202 if(pDz > zDiffMax) pDz = zDiffMax;
203 }
204 }
205
206 //***********************************************************************************************
207 //Define Event Variables
208 //***********************************************************************************************
209 //delta phi between the 2 leptons in degrees
210 double deltaPhiLeptons = MathUtils::DeltaPhi(CleanLeptons->At(0)->Phi(),
211 CleanLeptons->At(1)->Phi())* 180.0 / TMath::Pi();
212
213 double deltaEtaLeptons = abs(CleanLeptons->At(0)->Eta() - CleanLeptons->At(1)->Eta()) * 180.0 / TMath::Pi();
214
215 double deltaPhiDileptonMet = MathUtils::DeltaPhi(caloMet->Phi(),
216 dilepton->Phi())*180.0 / TMath::Pi();
217
218 double mtHiggs = TMath::Sqrt(2.0*dilepton->Pt() * caloMet->Pt()*
219 (1.0 - cos(deltaPhiDileptonMet * TMath::Pi() / 180.0)));
220
221 //angle between MET and closest lepton
222 double deltaPhiMetLepton[2] = {MathUtils::DeltaPhi(caloMet->Phi(), CleanLeptons->At(0)->Phi()),
223 MathUtils::DeltaPhi(caloMet->Phi(), CleanLeptons->At(1)->Phi())};
224
225 double mTW[2] = {TMath::Sqrt(2.0*CleanLeptons->At(0)->Pt()*caloMet->Pt()*
226 (1.0 - cos(deltaPhiMetLepton[0]))),
227 TMath::Sqrt(2.0*CleanLeptons->At(1)->Pt()*caloMet->Pt()*
228 (1.0 - cos(deltaPhiMetLepton[1])))};
229
230 double minDeltaPhiMetLepton = (deltaPhiMetLepton[0] < deltaPhiMetLepton[1])?
231 deltaPhiMetLepton[0]:deltaPhiMetLepton[1];
232
233 double METdeltaPhilEt = caloMet->Pt();
234 if(minDeltaPhiMetLepton < TMath::Pi()/2.)
235 METdeltaPhilEt = METdeltaPhilEt * sin(minDeltaPhiMetLepton);
236
237 //count the number of central Jets for vetoing and b-tagging
238 vector<Jet*> sortedJetsAll;
239 vector<Jet*> sortedJets;
240 vector<Jet*> sortedJetsLowPt;
241 for(UInt_t i=0; i<CleanJetsNoPtCut->GetEntries(); i++){
242 Jet* jet_a = new Jet(CleanJetsNoPtCut->At(i)->Px(),
243 CleanJetsNoPtCut->At(i)->Py(),
244 CleanJetsNoPtCut->At(i)->Pz(),
245 CleanJetsNoPtCut->At(i)->E() );
246
247 int nCloseStdJet = -1;
248 double deltaRMin = 999.;
249 for(UInt_t nj=0; nj<fCaloJet0->GetEntries(); nj++){
250 const CaloJet *jet = fCaloJet0->At(nj);
251 Double_t deltaR = MathUtils::DeltaR(jet_a->Mom(),jet->Mom());
252 if(deltaR < deltaRMin) {
253 nCloseStdJet = nj;
254 deltaRMin = deltaR;
255 }
256 }
257 if(nCloseStdJet >= 0 && deltaRMin < 0.5){
258 jet_a->SetMatchedMCFlavor(fCaloJet0->At(nCloseStdJet)->MatchedMCFlavor());
259 jet_a->SetCombinedSecondaryVertexBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->CombinedSecondaryVertexBJetTagsDisc());
260 jet_a->SetCombinedSecondaryVertexMVABJetTagsDisc(fCaloJet0->At(nCloseStdJet)->CombinedSecondaryVertexMVABJetTagsDisc());
261 jet_a->SetJetProbabilityBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->JetProbabilityBJetTagsDisc());
262 jet_a->SetJetBProbabilityBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->JetBProbabilityBJetTagsDisc());
263 jet_a->SetTrackCountingHighEffBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->TrackCountingHighEffBJetTagsDisc());
264 jet_a->SetTrackCountingHighPurBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->TrackCountingHighPurBJetTagsDisc());
265 jet_a->SetSimpleSecondaryVertexBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->SimpleSecondaryVertexBJetTagsDisc());
266 jet_a->SetSimpleSecondaryVertexHighEffBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->SimpleSecondaryVertexHighEffBJetTagsDisc());
267 jet_a->SetSimpleSecondaryVertexHighPurBJetTagsDisc(fCaloJet0->At(nCloseStdJet)->SimpleSecondaryVertexHighPurBJetTagsDisc());
268 }
269 else {
270 jet_a->SetMatchedMCFlavor(CleanJetsNoPtCut->At(i)->MatchedMCFlavor());
271 jet_a->SetCombinedSecondaryVertexBJetTagsDisc(CleanJetsNoPtCut->At(i)->CombinedSecondaryVertexBJetTagsDisc());
272 jet_a->SetCombinedSecondaryVertexMVABJetTagsDisc(CleanJetsNoPtCut->At(i)->CombinedSecondaryVertexMVABJetTagsDisc());
273 jet_a->SetJetProbabilityBJetTagsDisc(CleanJetsNoPtCut->At(i)->JetProbabilityBJetTagsDisc());
274 jet_a->SetJetBProbabilityBJetTagsDisc(CleanJetsNoPtCut->At(i)->JetBProbabilityBJetTagsDisc());
275 jet_a->SetTrackCountingHighEffBJetTagsDisc(CleanJetsNoPtCut->At(i)->TrackCountingHighEffBJetTagsDisc());
276 jet_a->SetTrackCountingHighPurBJetTagsDisc(CleanJetsNoPtCut->At(i)->TrackCountingHighPurBJetTagsDisc());
277 jet_a->SetSimpleSecondaryVertexBJetTagsDisc(CleanJetsNoPtCut->At(i)->SimpleSecondaryVertexBJetTagsDisc());
278 jet_a->SetSimpleSecondaryVertexHighEffBJetTagsDisc(CleanJetsNoPtCut->At(i)->SimpleSecondaryVertexHighEffBJetTagsDisc());
279 jet_a->SetSimpleSecondaryVertexHighPurBJetTagsDisc(CleanJetsNoPtCut->At(i)->SimpleSecondaryVertexHighPurBJetTagsDisc());
280 }
281 sortedJetsAll.push_back(jet_a);
282 }
283
284 for(UInt_t i=0; i<CleanJets->GetEntries(); i++){
285 if(TMath::Abs(CleanJets->At(i)->RawMom().Eta()) < 5.0 &&
286 CleanJets->At(i)->RawMom().Pt() > 25.0){
287 Jet* jet_b = new Jet(CleanJets->At(i)->Px(),
288 CleanJets->At(i)->Py(),
289 CleanJets->At(i)->Pz(),
290 CleanJets->At(i)->E() );
291 sortedJets.push_back(jet_b);
292 }
293 }
294
295 for(UInt_t i=0; i<sortedJetsAll.size(); i++){
296 bool overlap = kFALSE;
297 for(UInt_t j=0; j<sortedJets.size(); j++){
298 if(sortedJetsAll[i]->Pt() == sortedJets[j]->Pt() ||
299 (sortedJetsAll[i]->CombinedSecondaryVertexBJetTagsDisc() == sortedJets[j]->CombinedSecondaryVertexBJetTagsDisc() &&
300 sortedJetsAll[i]->JetBProbabilityBJetTagsDisc() == sortedJets[j]->JetBProbabilityBJetTagsDisc() &&
301 sortedJetsAll[i]->TrackCountingHighPurBJetTagsDisc() == sortedJets[j]->TrackCountingHighPurBJetTagsDisc())
302 ) {
303 sortedJets[j]->SetMatchedMCFlavor(sortedJetsAll[i]->MatchedMCFlavor());
304 sortedJets[j]->SetCombinedSecondaryVertexBJetTagsDisc(sortedJetsAll[i]->CombinedSecondaryVertexBJetTagsDisc());
305 sortedJets[j]->SetCombinedSecondaryVertexMVABJetTagsDisc(sortedJetsAll[i]->CombinedSecondaryVertexMVABJetTagsDisc());
306 sortedJets[j]->SetJetProbabilityBJetTagsDisc(sortedJetsAll[i]->JetProbabilityBJetTagsDisc());
307 sortedJets[j]->SetJetBProbabilityBJetTagsDisc(sortedJetsAll[i]->JetBProbabilityBJetTagsDisc());
308 sortedJets[j]->SetTrackCountingHighEffBJetTagsDisc(sortedJetsAll[i]->TrackCountingHighEffBJetTagsDisc());
309 sortedJets[j]->SetTrackCountingHighPurBJetTagsDisc(sortedJetsAll[i]->TrackCountingHighPurBJetTagsDisc());
310 sortedJets[j]->SetSimpleSecondaryVertexBJetTagsDisc(sortedJetsAll[i]->SimpleSecondaryVertexBJetTagsDisc());
311 sortedJets[j]->SetSimpleSecondaryVertexHighEffBJetTagsDisc(sortedJetsAll[i]->SimpleSecondaryVertexHighEffBJetTagsDisc());
312 sortedJets[j]->SetSimpleSecondaryVertexHighPurBJetTagsDisc(sortedJetsAll[i]->SimpleSecondaryVertexHighPurBJetTagsDisc());
313 overlap = kTRUE;
314 break;
315 }
316 }
317 if(overlap == kFALSE){
318 sortedJetsLowPt.push_back(sortedJetsAll[i]);
319 }
320 }
321 double maxBtag = -99999.;
322 double imaxBtag = -1;
323 for(UInt_t i=0; i<sortedJetsLowPt.size(); i++){
324 if(sortedJetsLowPt[i]->TrackCountingHighEffBJetTagsDisc() > maxBtag){
325 maxBtag = sortedJetsLowPt[i]->TrackCountingHighEffBJetTagsDisc();
326 imaxBtag = i;
327 }
328 }
329
330 //Lepton Type
331 int finalstateType = -1;
332 if (CleanLeptons->At(0)->ObjType() == kMuon && CleanLeptons->At(1)->ObjType() == kMuon ){ // mumu
333 finalstateType = 10;
334 } else if(CleanLeptons->At(0)->ObjType() == kElectron && CleanLeptons->At(1)->ObjType() == kElectron ){ // ee
335 finalstateType = 11;
336 } else if((CleanLeptons->At(0)->ObjType() == kElectron && CleanLeptons->At(1)->ObjType() == kMuon) ||
337 (CleanLeptons->At(1)->ObjType() == kElectron && CleanLeptons->At(0)->ObjType() == kMuon)) {
338 finalstateType = 12;
339 } else {
340 cerr << "Error: finalstate lepton type not supported\n";
341 }
342
343 //*********************************************************************************************
344 //Define Cuts
345 //*********************************************************************************************
346 const int nCuts = 10;
347 bool passCut[nCuts] = {false, false, false, false, false,
348 false, false, false, false, false};
349
350 if(CleanLeptons->At(0)->Pt() > 20.0 &&
351 CleanLeptons->At(1)->Pt() >= 20.0) passCut[0] = true;
352
353 if(zDiffMax < 1.0) passCut[1] = true;
354
355 if(caloMet->Pt() > 20.0) passCut[2] = true;
356
357 if(dilepton->Mass() > 12.0) passCut[3] = true;
358
359 if(sortedJets.size() < 1) passCut[6] = true;
360
361 if(SoftMuons->GetEntries() == 0) passCut[7] = true;
362
363 if(CleanLeptons->GetEntries() == 2) passCut[8] = true;
364
365 if(maxBtag < 2.1) passCut[9] = true;
366
367 if (finalstateType == 10 || finalstateType == 11){ // mumu/ee
368 if(fabs(dilepton->Mass()-91.1876) > 15.0) passCut[4] = true;
369 if(METdeltaPhilEt > 35) passCut[5] = true;
370 }
371 else if(finalstateType == 12) { // emu
372 passCut[4] = true;
373 if(METdeltaPhilEt > 20) passCut[5] = true;
374 }
375
376 //*********************************************************************************************
377 //Make Selection Histograms. Number of events passing each level of cut
378 //*********************************************************************************************
379 bool passAllCuts = true;
380 for(int c=0; c<nCuts; c++) passAllCuts = passAllCuts & passCut[c];
381
382 //Cut Selection Histograms
383 fHWWSelection->Fill(-1,NNLOWeight->GetVal());
384 if (finalstateType == 10 )
385 fHWWToMuMuSelection->Fill(-1,NNLOWeight->GetVal());
386 else if(finalstateType == 11 )
387 fHWWToEESelection->Fill(-1,NNLOWeight->GetVal());
388 else if(finalstateType == 12 )
389 fHWWToEMuSelection->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 }
409 }
410
411 //*****************************************************************************************
412 //Make Preselection Histograms
413 //*****************************************************************************************
414 fLeptonEta->Fill(CleanLeptons->At(0)->Eta(),NNLOWeight->GetVal());
415 fLeptonEta->Fill(CleanLeptons->At(1)->Eta(),NNLOWeight->GetVal());
416 fLeptonPtMax->Fill(CleanLeptons->At(0)->Pt(),NNLOWeight->GetVal());
417 fLeptonPtMin->Fill(CleanLeptons->At(1)->Pt(),NNLOWeight->GetVal());
418 fMetPtHist->Fill(caloMet->Pt(),NNLOWeight->GetVal());
419 fMetPhiHist->Fill(caloMet->Phi(),NNLOWeight->GetVal());
420 fDeltaPhiLeptons->Fill(deltaPhiLeptons,NNLOWeight->GetVal());
421 fDeltaEtaLeptons->Fill(deltaEtaLeptons,NNLOWeight->GetVal());
422 fDileptonMass->Fill(dilepton->Mass(),NNLOWeight->GetVal());
423
424 //*********************************************************************************************
425 // N-1 Histograms
426 //*********************************************************************************************
427 bool pass;;
428
429 //N Jet Veto
430 pass = true;
431 for (int k=0;k<nCuts;k++) {
432 if (k != 5) {
433 pass = (pass && passCut[k]);
434 }
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) {
444 pass = (pass && passCut[k]);
445 }
446 }
447 if (pass) {
448 fMetPtHist_NMinusOne->Fill(caloMet->Pt(),NNLOWeight->GetVal());
449 }
450
451 // dilepton mass
452 pass = true;
453 for (int k=0;k<nCuts;k++) {
454 if (k != 2 && k != 3)
455 pass = (pass && passCut[k]);
456 }
457 if (pass) {
458 fDileptonMass_NMinusOne->Fill(dilepton->Mass(),NNLOWeight->GetVal());
459 }
460
461 // Lepton Pt Max, Lepton Pt Min, DeltaPhiLeptons
462 pass = true;
463 for (int k=0;k<nCuts;k++) {
464 pass = (pass && passCut[k]);
465 }
466 if (pass) {
467 fLeptonPtMax_NMinusOne->Fill(CleanLeptons->At(0)->Pt(),NNLOWeight->GetVal());
468 fLeptonPtMin_NMinusOne->Fill(CleanLeptons->At(1)->Pt(),NNLOWeight->GetVal());
469 fDeltaPhiLeptons_NMinusOne->Fill(deltaPhiLeptons,NNLOWeight->GetVal());
470 }
471
472 // NSoftMuons
473 pass = true;
474 for (int k=0;k<nCuts;k++) {
475 if (k != 6)
476 pass = (pass && passCut[k]);
477 }
478 if (pass) {
479 fNSoftMuonsHist_NMinusOne->Fill(SoftMuons->GetEntries(),NNLOWeight->GetVal());
480 }
481
482 //*********************************************************************************************
483 //Plots after all Cuts
484 //*********************************************************************************************
485 if (passAllCuts) {
486 fMinDeltaPhiLeptonMet_afterCuts->Fill(minDeltaPhiMetLepton,NNLOWeight->GetVal());
487 fMtLepton1_afterCuts->Fill(mTW[0],NNLOWeight->GetVal());
488 fMtLepton2_afterCuts->Fill(mTW[1],NNLOWeight->GetVal());
489 fMtHiggs_afterCuts->Fill(mtHiggs,NNLOWeight->GetVal());
490 fLeptonPtPlusMet_afterCuts->Fill(CleanLeptons->At(0)->Pt()+CleanLeptons->At(1)->Pt()+caloMet->Pt(),NNLOWeight->GetVal());
491 }
492
493 delete dilepton;
494 delete SoftMuons;
495 for(UInt_t i=0; i<sortedJets.size(); i++) delete sortedJets[i];
496 for(UInt_t i=0; i<sortedJetsAll.size(); i++) delete sortedJetsAll[i];
497 return;
498 }
499
500 //--------------------------------------------------------------------------------------------------
501 void HwwExampleAnalysisMod::SlaveTerminate()
502 {
503
504 // Run finishing code on the computer (slave) that did the analysis. For this
505 // module, we dont do anything here.
506
507 }
508 //--------------------------------------------------------------------------------------------------
509 void HwwExampleAnalysisMod::Terminate()
510 {
511 // Run finishing code on the client computer. For this module, we dont do
512 // anything here.
513
514 }