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