ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/SelMods/src/HwwExampleAnalysisMod.cc
Revision: 1.22
Committed: Sat Apr 7 11:46:43 2012 UTC (13 years, 1 month ago) by ceballos
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_029c, Mit_029b, Mit_029a, Mit_028a, Mit_028, Mit_027, Mit_027a, HEAD
Changes since 1.21: +7 -10 lines
Log Message:
removing cleanjets

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