1 |
loizides |
1.2 |
// $Id: FakeLeptonExampleAnaMod.cc,v 1.1 2009/06/30 10:47:17 loizides Exp $
|
2 |
loizides |
1.1 |
|
3 |
|
|
#include "MitPhysics/FakeMods/interface/FakeLeptonExampleAnaMod.h"
|
4 |
|
|
#include "MitAna/DataUtil/interface/Debug.h"
|
5 |
|
|
#include "MitAna/DataTree/interface/Names.h"
|
6 |
|
|
#include "MitCommon/MathTools/interface/MathUtils.h"
|
7 |
|
|
#include "MitPhysics/Init/interface/ModNames.h"
|
8 |
|
|
#include "MitPhysics/FakeMods/interface/FakeObject.h"
|
9 |
|
|
#include "MitPhysics/FakeMods/interface/FakeEventHeader.h"
|
10 |
|
|
#include <TH1D.h>
|
11 |
|
|
#include <TH2D.h>
|
12 |
|
|
|
13 |
|
|
using namespace mithep;
|
14 |
|
|
|
15 |
|
|
ClassImp(mithep::FakeLeptonExampleAnaMod)
|
16 |
|
|
|
17 |
|
|
//--------------------------------------------------------------------------------------------------
|
18 |
|
|
FakeLeptonExampleAnaMod::FakeLeptonExampleAnaMod(const char *name, const char *title) :
|
19 |
|
|
BaseMod(name,title),
|
20 |
|
|
fUseMCFake(false),
|
21 |
|
|
fPerformFakeMuonMetCorrection(true),
|
22 |
|
|
fSampleName("NotSet"),
|
23 |
|
|
fFakeEventHeaderName(ModNames::gkFakeEventHeadersName),
|
24 |
|
|
fElectronFakeableObjectsName(ModNames::gkElFakeableObjsName),
|
25 |
|
|
fMuonFakeableObjectsName(ModNames::gkMuFakeableObjsName),
|
26 |
|
|
fMCPartBranchName(Names::gkMCPartBrn),
|
27 |
|
|
fGenJetBranchName(Names::gkSC5GenJetBrn),
|
28 |
|
|
fTrackBranchName(Names::gkTrackBrn),
|
29 |
|
|
fMuonBranchName(Names::gkMuonBrn),
|
30 |
|
|
fMetName("NotSet"),
|
31 |
|
|
fCleanJetsName(ModNames::gkCleanJetsName),
|
32 |
|
|
fTriggerObjectsName("NotSet"),
|
33 |
|
|
fParticles(0),
|
34 |
|
|
fGenJets(0),
|
35 |
|
|
fTracks(0),
|
36 |
|
|
fMuons(0),
|
37 |
|
|
fMet(0)
|
38 |
|
|
{
|
39 |
|
|
// Constructor.
|
40 |
|
|
}
|
41 |
|
|
|
42 |
|
|
//--------------------------------------------------------------------------------------------------
|
43 |
|
|
void FakeLeptonExampleAnaMod::Begin()
|
44 |
|
|
{
|
45 |
|
|
// Run startup code on the client machine. For this module, we dont do
|
46 |
|
|
// anything here.
|
47 |
|
|
}
|
48 |
|
|
|
49 |
|
|
//--------------------------------------------------------------------------------------------------
|
50 |
|
|
void FakeLeptonExampleAnaMod::Process()
|
51 |
|
|
{
|
52 |
|
|
// Process entries of the tree.
|
53 |
|
|
LoadBranch(fTrackBranchName);
|
54 |
|
|
LoadBranch(fMuonBranchName);
|
55 |
|
|
|
56 |
|
|
//***********************************************************************************************
|
57 |
|
|
//Import Collections
|
58 |
|
|
//***********************************************************************************************
|
59 |
|
|
|
60 |
|
|
//Obtain all cleaned objects
|
61 |
|
|
ElectronOArr *CleanElectrons = dynamic_cast<ElectronOArr* >
|
62 |
|
|
(FindObjThisEvt(ModNames::gkCleanElectronsName));
|
63 |
|
|
MuonOArr *CleanMuons = dynamic_cast<MuonOArr* >
|
64 |
|
|
(FindObjThisEvt(ModNames::gkCleanMuonsName));
|
65 |
|
|
mithep::ParticleOArr *CleanLeptons = dynamic_cast<mithep::ParticleOArr*>
|
66 |
|
|
(FindObjThisEvt(ModNames::gkMergedLeptonsName));
|
67 |
|
|
|
68 |
|
|
//Get Met
|
69 |
|
|
if (!fMetName.IsNull()) {
|
70 |
|
|
fMet = GetObjThisEvt<MetCol>(fMetName);
|
71 |
|
|
}
|
72 |
|
|
const Met *originalCaloMet = 0;
|
73 |
|
|
if (fMet) {
|
74 |
|
|
originalCaloMet = fMet->At(0);
|
75 |
|
|
} else {
|
76 |
|
|
cout << "Error: Met Collection " << fMetName << " could not be loaded.\n";
|
77 |
|
|
}
|
78 |
|
|
ObjArray<Jet> *OriginalCleanJets = dynamic_cast<ObjArray<Jet>* >
|
79 |
|
|
(FindObjThisEvt(fCleanJetsName.Data()));
|
80 |
|
|
|
81 |
|
|
//Obtain the collection of fake objects
|
82 |
|
|
ElectronCol *ElectronFakeableObjects = 0;
|
83 |
|
|
if(!fElectronFakeableObjectsName.IsNull())
|
84 |
|
|
ElectronFakeableObjects = dynamic_cast<ElectronCol* >
|
85 |
|
|
(FindObjThisEvt(fElectronFakeableObjectsName.Data()));
|
86 |
|
|
MuonCol *MuonFakeableObjects = 0;
|
87 |
|
|
if (!fMuonFakeableObjectsName.IsNull())
|
88 |
|
|
MuonFakeableObjects = dynamic_cast<MuonCol* >
|
89 |
|
|
(FindObjThisEvt(fMuonFakeableObjectsName.Data()));
|
90 |
|
|
ChargedParticleOArr *FakeableObjects = new ChargedParticleOArr;
|
91 |
|
|
if (ElectronFakeableObjects) {
|
92 |
|
|
for (UInt_t i=0; i<ElectronFakeableObjects->GetEntries(); i++)
|
93 |
|
|
FakeableObjects->Add(ElectronFakeableObjects->At(i));
|
94 |
|
|
}
|
95 |
|
|
if (MuonFakeableObjects) {
|
96 |
|
|
for (UInt_t i=0; i<MuonFakeableObjects->GetEntries(); i++)
|
97 |
|
|
FakeableObjects->Add(MuonFakeableObjects->At(i));
|
98 |
|
|
}
|
99 |
|
|
Collection<FakeEventHeader> *FakeEventHeaders = 0;
|
100 |
|
|
if (!fUseMCFake) {
|
101 |
|
|
if (!fFakeEventHeaderName.IsNull()) {
|
102 |
|
|
FakeEventHeaders = dynamic_cast<Collection<FakeEventHeader>* >(FindObjThisEvt(fFakeEventHeaderName.Data()));
|
103 |
|
|
if (!FakeEventHeaders) {
|
104 |
|
|
cout << "Error: FakeEventHeader with name " << fFakeEventHeaderName.Data() << " could not be loaded.\n";
|
105 |
|
|
assert(false);
|
106 |
|
|
}
|
107 |
|
|
} else
|
108 |
|
|
cout << "Error: FakeEventHeaders " << fFakeEventHeaderName.Data() << " could not be loaded.\n";
|
109 |
|
|
}
|
110 |
|
|
|
111 |
|
|
//***********************************************************************************************
|
112 |
|
|
//If we use MC Fakes, then create a new FakeEventHeader containing no fakes with weight = 1
|
113 |
|
|
//This ensures that in the loop over FakeEventHeaders we do the correct thing.
|
114 |
|
|
//***********************************************************************************************
|
115 |
|
|
if (fUseMCFake) {
|
116 |
|
|
ObjArray <FakeEventHeader> *tmpFakeEventHeaders = new ObjArray <FakeEventHeader> ;
|
117 |
|
|
tmpFakeEventHeaders->SetOwner(kTRUE);
|
118 |
|
|
|
119 |
|
|
FakeEventHeader *initialFakeEvent = new FakeEventHeader();
|
120 |
|
|
for (UInt_t j=0;j<OriginalCleanJets->GetEntries();j++)
|
121 |
|
|
initialFakeEvent->AddJets(OriginalCleanJets->At(j));
|
122 |
|
|
|
123 |
|
|
tmpFakeEventHeaders->AddOwned(initialFakeEvent);
|
124 |
|
|
FakeEventHeaders = dynamic_cast<Collection<FakeEventHeader>* > (tmpFakeEventHeaders);
|
125 |
|
|
}
|
126 |
|
|
|
127 |
|
|
//***********************************************************************************************
|
128 |
|
|
//Loop over Fake Event Headers
|
129 |
|
|
//***********************************************************************************************
|
130 |
|
|
for (UInt_t i=0;i<FakeEventHeaders->GetEntries() ; i++) {
|
131 |
|
|
|
132 |
|
|
//Create leptons collection containing real leptons and fake leptons
|
133 |
|
|
ObjArray<Particle> *leptons = NULL;
|
134 |
|
|
if (fUseMCFake) {
|
135 |
|
|
leptons = CleanLeptons;
|
136 |
|
|
} else {
|
137 |
|
|
leptons = new ObjArray<Particle>;
|
138 |
|
|
|
139 |
|
|
for (UInt_t j=0;j<CleanLeptons->GetEntries() ; j++) {
|
140 |
|
|
leptons->Add(CleanLeptons->At(j));
|
141 |
|
|
}
|
142 |
|
|
for (UInt_t j=0;j<FakeEventHeaders->At(i)->FakeObjsSize() ; j++) {
|
143 |
|
|
leptons->Add(FakeEventHeaders->At(i)->FakeObj(j)->FakeParticle());
|
144 |
|
|
}
|
145 |
|
|
}
|
146 |
|
|
//we have to sort leptons
|
147 |
|
|
leptons->Sort();
|
148 |
|
|
|
149 |
|
|
|
150 |
|
|
//Construct the Clean Jet collection.
|
151 |
|
|
ObjArray<Jet> *CleanJets = NULL;
|
152 |
|
|
if (fUseMCFake) {
|
153 |
|
|
CleanJets = OriginalCleanJets;
|
154 |
|
|
} else {
|
155 |
|
|
CleanJets = new ObjArray<Jet>;
|
156 |
|
|
for (UInt_t j=0;j<FakeEventHeaders->At(i)->NJets() ; j++) {
|
157 |
|
|
CleanJets->Add(FakeEventHeaders->At(i)->UnfakedJet(j));
|
158 |
|
|
}
|
159 |
|
|
}
|
160 |
|
|
|
161 |
|
|
//Perform correction for potential fake muons
|
162 |
|
|
//have to add fake muon momentum to originalCaloMet;
|
163 |
|
|
const Met *caloMet = originalCaloMet;
|
164 |
|
|
Double_t FakeMuonMetCorrection_X = 0.0;
|
165 |
|
|
Double_t FakeMuonMetCorrection_Y = 0.0;
|
166 |
|
|
for (UInt_t j=0;j<FakeEventHeaders->At(i)->FakeObjsSize() ; j++) {
|
167 |
|
|
if (FakeEventHeaders->At(i)->FakeObj(j)->ObjType() == kMuon) {
|
168 |
|
|
FakeMuonMetCorrection_X += FakeEventHeaders->At(i)->FakeObj(j)->Px();
|
169 |
|
|
FakeMuonMetCorrection_Y += FakeEventHeaders->At(i)->FakeObj(j)->Py();
|
170 |
|
|
}
|
171 |
|
|
}
|
172 |
|
|
|
173 |
|
|
if (!fUseMCFake && fPerformFakeMuonMetCorrection) {
|
174 |
|
|
caloMet = new Met(originalCaloMet->Px()+FakeMuonMetCorrection_X,
|
175 |
|
|
originalCaloMet->Py()+FakeMuonMetCorrection_Y);
|
176 |
|
|
}
|
177 |
|
|
|
178 |
|
|
//*********************************************************************************************
|
179 |
|
|
//Construct the event weight using fake rate and corrections
|
180 |
|
|
//*********************************************************************************************
|
181 |
|
|
//fake rate has to be corrected by the amount lost when those denominators
|
182 |
|
|
//became fakes in data. If a denominator fakes a lepton in data, it goes in the 2lepton
|
183 |
|
|
//final state, and we don't count it in this prediction. So we have to add it back.
|
184 |
|
|
Double_t eventweight = FakeEventHeaders->At(i)->Weight();
|
185 |
|
|
if (FakeEventHeaders->At(i)->FakeObjsSize() > 0 && FakeEventHeaders->At(i)->Weight() < 1) {
|
186 |
|
|
eventweight = eventweight / (1.0 - FakeEventHeaders->At(i)->Weight());
|
187 |
|
|
}
|
188 |
|
|
|
189 |
|
|
//*********************************************************************************************
|
190 |
|
|
//another correction to account for events lost due to only the fake lepton firing the trigger
|
191 |
|
|
//The numbers need to be changed for your analysis.
|
192 |
|
|
//Given numbers are for the 2 lepton final state.
|
193 |
|
|
//*********************************************************************************************
|
194 |
|
|
if (CleanLeptons->GetEntries() >= 1 && FakeEventHeaders->At(i)->FakeObjsSize() >= 1) {
|
195 |
|
|
if (CleanLeptons->At(0)->ObjType() == kElectron &&
|
196 |
|
|
FakeEventHeaders->At(i)->FakeObj(0)->FakeParticle()->ObjType() == kElectron) {
|
197 |
|
|
eventweight = eventweight * 1.06;
|
198 |
|
|
} else if (CleanLeptons->At(0)->ObjType() == kMuon &&
|
199 |
|
|
FakeEventHeaders->At(i)->FakeObj(0)->FakeParticle()->ObjType() == kMuon) {
|
200 |
|
|
eventweight = eventweight * 1.12;
|
201 |
|
|
} else if (CleanLeptons->At(0)->ObjType() == kElectron &&
|
202 |
|
|
FakeEventHeaders->At(i)->FakeObj(0)->FakeParticle()->ObjType() == kMuon) {
|
203 |
|
|
eventweight = eventweight * 1.17;
|
204 |
|
|
} else if (CleanLeptons->At(0)->ObjType() == kMuon &&
|
205 |
|
|
FakeEventHeaders->At(i)->FakeObj(0)->FakeParticle()->ObjType() == kElectron) {
|
206 |
|
|
eventweight = eventweight * 1.17;
|
207 |
|
|
}
|
208 |
|
|
}
|
209 |
|
|
|
210 |
|
|
//***********************************************************************************************
|
211 |
|
|
//For FR method (fUseMCFake == false)
|
212 |
|
|
//Make analysis specific cuts.
|
213 |
|
|
//For example for 2 lepton final state we require that the event contains
|
214 |
|
|
//one and only one clean lepton with pt > 10 GeV.
|
215 |
|
|
//***********************************************************************************************
|
216 |
|
|
if (!fUseMCFake) {
|
217 |
|
|
if (CleanLeptons->GetEntries() != 1 || CleanLeptons->At(0)->Pt() <= 10.0)
|
218 |
|
|
continue;
|
219 |
|
|
}
|
220 |
|
|
|
221 |
|
|
//*********************************************************************************************
|
222 |
|
|
//Fill some distributions before preselection
|
223 |
|
|
//*********************************************************************************************
|
224 |
|
|
|
225 |
|
|
CompositeParticle *dilepton = NULL;
|
226 |
|
|
|
227 |
|
|
if (leptons->GetEntries()>=2) {
|
228 |
|
|
|
229 |
|
|
dilepton = new CompositeParticle();
|
230 |
|
|
dilepton->AddDaughter(leptons->At(0));
|
231 |
|
|
dilepton->AddDaughter(leptons->At(1));
|
232 |
|
|
|
233 |
|
|
//Dilepton Charge will be filled like this
|
234 |
|
|
// -2: -- , -1: -+, 1: +-, 2:++
|
235 |
|
|
if (dilepton->Charge() == 0) {
|
236 |
|
|
if (leptons->At(0)->Charge() == 1) {
|
237 |
|
|
fDileptonCharge->Fill(1.0,eventweight);
|
238 |
|
|
} else {
|
239 |
|
|
fDileptonCharge->Fill(-1.0,eventweight);
|
240 |
|
|
}
|
241 |
|
|
} else {
|
242 |
|
|
fDileptonCharge->Fill(dilepton->Charge(),eventweight);
|
243 |
|
|
}
|
244 |
|
|
}
|
245 |
|
|
|
246 |
|
|
//*********************************************************************************************
|
247 |
|
|
//Kinematic PreSelection
|
248 |
|
|
//Example given is for the two lepton final state
|
249 |
|
|
//*********************************************************************************************
|
250 |
|
|
//make sure 2nd highest pt lepton has Pt > 10
|
251 |
|
|
if (leptons->GetEntries() < 2 || leptons->At(1)->Pt() <= 10) continue;
|
252 |
|
|
|
253 |
|
|
//make sure the 3rd highest pt lepton has pt <= 10.
|
254 |
|
|
if (leptons->GetEntries() >= 3 && leptons->At(2)->Pt() > 10) continue;
|
255 |
|
|
|
256 |
|
|
//charge of the leptons should be opposite
|
257 |
|
|
if (dilepton->Charge() != 0) continue;
|
258 |
|
|
|
259 |
|
|
//*********************************************************************************************
|
260 |
|
|
//Get nonisolated soft muons
|
261 |
|
|
//*********************************************************************************************
|
262 |
|
|
ObjArray<Muon> *DirtyMuons = new ObjArray<Muon>;
|
263 |
|
|
for (UInt_t m=0; m<fMuons->GetEntries(); ++m) {
|
264 |
|
|
const Muon *mu = fMuons->At(m);
|
265 |
|
|
if(!mu->GlobalTrk()) continue;
|
266 |
|
|
if(mu->Pt() < 5.0) continue;
|
267 |
|
|
|
268 |
|
|
//remove the fake
|
269 |
|
|
bool isFakedMuon = false;
|
270 |
|
|
for (UInt_t f=0;f<FakeEventHeaders->At(i)->FakeObjsSize() ; f++) {
|
271 |
|
|
if (mu->HasTrackerTrk() &&
|
272 |
|
|
(dynamic_cast<const mithep::ChargedParticle*>
|
273 |
|
|
(FakeEventHeaders->At(i)->FakeObj(f)->FakeParticle()))->TrackerTrk() &&
|
274 |
|
|
(dynamic_cast<const mithep::ChargedParticle*>
|
275 |
|
|
(FakeEventHeaders->At(i)->FakeObj(f)->FakeParticle()))->TrackerTrk() ==
|
276 |
|
|
mu->TrackerTrk()
|
277 |
|
|
)
|
278 |
|
|
isFakedMuon = true;
|
279 |
|
|
}
|
280 |
|
|
|
281 |
|
|
//remove clean muons
|
282 |
|
|
bool isCleanMuon = false;
|
283 |
|
|
for (UInt_t j=0; j<CleanMuons->GetEntries(); j++) {
|
284 |
|
|
if(fMuons->At(m) == CleanMuons->At(j)) isCleanMuon = true;
|
285 |
|
|
}
|
286 |
|
|
|
287 |
|
|
if(!isCleanMuon
|
288 |
|
|
&& !(isFakedMuon && !fUseMCFake)
|
289 |
|
|
) DirtyMuons->Add(mu);
|
290 |
|
|
}
|
291 |
|
|
|
292 |
|
|
//*********************************************************************************************
|
293 |
|
|
//Get Clean Tracks excluding the good leptons
|
294 |
|
|
//*********************************************************************************************
|
295 |
|
|
ObjArray<Track> *CleanExtraTracks = new ObjArray<Track>;
|
296 |
|
|
int nTracks = 0;
|
297 |
|
|
|
298 |
|
|
double z0Average = ( (dynamic_cast<const mithep::ChargedParticle*>(leptons->At(0)))->Trk()->Z0()
|
299 |
|
|
+ (dynamic_cast<const mithep::ChargedParticle*>(leptons->At(1)))->Trk()->Z0()) /2 ;
|
300 |
|
|
|
301 |
|
|
for (UInt_t t=0; t<fTracks->GetEntries(); ++t) {
|
302 |
|
|
bool isLepton = false;
|
303 |
|
|
|
304 |
|
|
if (MathUtils::DeltaR(fTracks->At(t)->Phi(),fTracks->At(t)->Eta(),leptons->At(0)->Phi(),
|
305 |
|
|
leptons->At(0)->Eta()) > 0.01 &&
|
306 |
|
|
MathUtils::DeltaR(fTracks->At(t)->Phi(),fTracks->At(t)->Eta(),leptons->At(1)->Phi(),
|
307 |
|
|
leptons->At(1)->Eta()) > 0.01
|
308 |
|
|
) {
|
309 |
|
|
} else {
|
310 |
|
|
isLepton = true;
|
311 |
|
|
}
|
312 |
|
|
|
313 |
|
|
MDB(kAnalysis, 8) {
|
314 |
|
|
cout << "Track " << t << " : " << fTracks->At(t)->Pt() << " " << fTracks->At(t)->Eta()
|
315 |
|
|
<< " " << fTracks->At(t)->Phi() << " islepton=" << isLepton << endl;
|
316 |
|
|
}
|
317 |
|
|
|
318 |
|
|
if ( !isLepton && fTracks->At(t)->Pt() > 3.0
|
319 |
|
|
&& fTracks->At(t)->NHits() >= 8
|
320 |
|
|
&& fabs(z0Average - fTracks->At(t)->Z0()) < 0.5 ) {
|
321 |
|
|
CleanExtraTracks->Add(fTracks->At(t));
|
322 |
|
|
nTracks++;
|
323 |
|
|
}
|
324 |
|
|
}
|
325 |
|
|
|
326 |
|
|
//*********************************************************************************************
|
327 |
|
|
//The code below is an example analysis for the HWW analysis.
|
328 |
|
|
//*********************************************************************************************
|
329 |
|
|
|
330 |
|
|
|
331 |
|
|
//*********************************************************************************************
|
332 |
|
|
//Define Event Variables
|
333 |
|
|
//*********************************************************************************************
|
334 |
|
|
//delta phi between the 2 leptons in degrees
|
335 |
|
|
double deltaPhiLeptons = MathUtils::DeltaPhi(leptons->At(0)->Phi(),
|
336 |
|
|
leptons->At(1)->Phi())* 360.0 / 2 / TMath::Pi();
|
337 |
|
|
|
338 |
|
|
double deltaEtaLeptons = leptons->At(0)->Eta() - leptons->At(1)->Eta();
|
339 |
|
|
|
340 |
|
|
double deltaPhiDileptonMet = MathUtils::DeltaPhi(caloMet->Phi(),
|
341 |
|
|
dilepton->Phi())*360.0 / 2 / TMath::Pi();
|
342 |
|
|
|
343 |
|
|
double mtHiggs = TMath::Sqrt(2.0*dilepton->Pt() * caloMet->Pt()*
|
344 |
|
|
(1.0 - cos(deltaPhiDileptonMet * 2 * TMath::Pi() / 360.0)));
|
345 |
|
|
|
346 |
|
|
//angle between MET and closest lepton
|
347 |
|
|
double deltaPhiMetLepton[2] = {MathUtils::DeltaPhi(caloMet->Phi(), leptons->At(0)->Phi()),
|
348 |
|
|
MathUtils::DeltaPhi(caloMet->Phi(), leptons->At(1)->Phi())};
|
349 |
|
|
|
350 |
|
|
double mTW[2] = {TMath::Sqrt(2.0*leptons->At(0)->Pt()*caloMet->Pt()*
|
351 |
|
|
(1.0 - cos(deltaPhiMetLepton[0]))),
|
352 |
|
|
TMath::Sqrt(2.0*leptons->At(1)->Pt()*caloMet->Pt()*
|
353 |
|
|
(1.0 - cos(deltaPhiMetLepton[1])))};
|
354 |
|
|
|
355 |
|
|
double minDeltaPhiMetLepton = (deltaPhiMetLepton[0] < deltaPhiMetLepton[1])?
|
356 |
|
|
deltaPhiMetLepton[0]:deltaPhiMetLepton[1];
|
357 |
|
|
minDeltaPhiMetLepton = minDeltaPhiMetLepton * 360.0 / 2 / TMath::Pi();
|
358 |
|
|
|
359 |
|
|
//count the number of central Jets for vetoing
|
360 |
|
|
int nCentralJets = 0;
|
361 |
|
|
for (UInt_t j=0; j<CleanJets->GetEntries(); j++) {
|
362 |
|
|
if (fabs(CleanJets->At(j)->Eta()) < 2.5)
|
363 |
|
|
nCentralJets++;
|
364 |
|
|
}
|
365 |
|
|
|
366 |
|
|
//Lepton Type
|
367 |
|
|
int finalstateType = -1;
|
368 |
|
|
if (leptons->At(0)->ObjType() == kMuon && leptons->At(1)->ObjType() == kMuon ){ // mumu
|
369 |
|
|
finalstateType = 10;
|
370 |
|
|
} else if(leptons->At(0)->ObjType() == kElectron && leptons->At(1)->ObjType() == kElectron ){ // ee
|
371 |
|
|
finalstateType = 11;
|
372 |
|
|
} else if((leptons->At(0)->ObjType() == kElectron && leptons->At(1)->ObjType() == kMuon) ||
|
373 |
|
|
(leptons->At(0)->ObjType() == kMuon && leptons->At(1)->ObjType() == kElectron)) {
|
374 |
|
|
finalstateType = 12;
|
375 |
|
|
} else {
|
376 |
|
|
cerr << "Error: finalstate lepton type not supported\n";
|
377 |
|
|
}
|
378 |
|
|
|
379 |
|
|
//*********************************************************************************************
|
380 |
|
|
//Define Cuts
|
381 |
|
|
//*********************************************************************************************
|
382 |
|
|
const int nCuts = 9;
|
383 |
|
|
bool passCut[nCuts] = {false, false, false, false,
|
384 |
|
|
false, false, false, false, false};
|
385 |
|
|
|
386 |
|
|
if(leptons->At(0)->Pt() > 20.0 &&
|
387 |
|
|
leptons->At(1)->Pt() > 10.0 &&
|
388 |
|
|
caloMet->Pt() > 30.0 &&
|
389 |
|
|
dilepton->Mass() > 12.0
|
390 |
|
|
) passCut[0] = true;
|
391 |
|
|
//above cuts are for preselction to be fed into TMVA
|
392 |
|
|
|
393 |
|
|
if(nCentralJets < 1) passCut[1] = true;
|
394 |
|
|
|
395 |
|
|
if (finalstateType == 10){ // mumu
|
396 |
|
|
if(caloMet->Pt() > 50.0 &&
|
397 |
|
|
caloMet->Pt() < 200.0) passCut[2] = true;
|
398 |
|
|
if(deltaPhiLeptons < 45.0) passCut[3] = true;
|
399 |
|
|
if(dilepton->Mass() < 50.0) passCut[4] = true;
|
400 |
|
|
if(leptons->At(0)->Pt() > 35.0 &&
|
401 |
|
|
leptons->At(0)->Pt() < 55.0) passCut[5] = true;
|
402 |
|
|
if(leptons->At(1)->Pt() > 25.0) passCut[6] = true;
|
403 |
|
|
}
|
404 |
|
|
else if(finalstateType == 11 ){ // ee
|
405 |
|
|
if(caloMet->Pt() > 51.0 &&
|
406 |
|
|
caloMet->Pt() < 200.0) passCut[2] = true;
|
407 |
|
|
if(deltaPhiLeptons < 45.0) passCut[3] = true;
|
408 |
|
|
if(dilepton->Mass() < 40.0) passCut[4] = true;
|
409 |
|
|
if(leptons->At(0)->Pt() > 25.0 &&
|
410 |
|
|
leptons->At(0)->Pt() < 49.0) passCut[5] = true;
|
411 |
|
|
if(leptons->At(1)->Pt() > 25.0) passCut[6] = true;
|
412 |
|
|
}
|
413 |
|
|
else if(finalstateType == 12) { //emu
|
414 |
|
|
if(caloMet->Pt() > 45.0 &&
|
415 |
|
|
caloMet->Pt() < 105.0) passCut[2] = true;
|
416 |
|
|
if(deltaPhiLeptons < 70.0) passCut[3] = true;
|
417 |
|
|
if(dilepton->Mass() < 45.0) passCut[4] = true;
|
418 |
|
|
if(leptons->At(0)->Pt() > 25.0 &&
|
419 |
|
|
leptons->At(0)->Pt() < 50.0) passCut[5] = true;
|
420 |
|
|
if(leptons->At(1)->Pt() > 25.0) passCut[6] = true;
|
421 |
|
|
}
|
422 |
|
|
|
423 |
|
|
if (DirtyMuons->GetEntries() < 1) passCut[7] = true;
|
424 |
|
|
if (CleanExtraTracks->GetEntries() < 4) passCut[8] = true;
|
425 |
|
|
|
426 |
|
|
//*********************************************************************************************
|
427 |
|
|
//Final Decision
|
428 |
|
|
//*********************************************************************************************
|
429 |
|
|
bool passAllCuts = true;
|
430 |
|
|
for(int c=0; c<nCuts; c++) passAllCuts = passAllCuts & passCut[c];
|
431 |
|
|
|
432 |
|
|
//*****************************************************************************************
|
433 |
|
|
//Histograms after no cuts
|
434 |
|
|
//*****************************************************************************************
|
435 |
|
|
fLeptonPtMax->Fill(leptons->At(0)->Pt(),eventweight);
|
436 |
|
|
fLeptonPtMin->Fill(leptons->At(1)->Pt(),eventweight);
|
437 |
|
|
fMetPtHist->Fill(caloMet->Pt(),eventweight);
|
438 |
|
|
fDeltaPhiLeptons->Fill(deltaPhiLeptons,eventweight);
|
439 |
|
|
fDeltaEtaLeptons->Fill(deltaEtaLeptons,eventweight);
|
440 |
|
|
fDileptonMass->Fill(dilepton->Mass(),eventweight);
|
441 |
|
|
|
442 |
|
|
//*********************************************************************************************
|
443 |
|
|
// N-1 Histograms
|
444 |
|
|
//*********************************************************************************************
|
445 |
|
|
|
446 |
|
|
//N Jet Veto
|
447 |
|
|
bool pass = true;
|
448 |
|
|
for (int k=0;k<nCuts;k++) {
|
449 |
|
|
if (k != 1) {
|
450 |
|
|
pass = (pass && passCut[k]);
|
451 |
|
|
}
|
452 |
|
|
}
|
453 |
|
|
if (pass) {
|
454 |
|
|
fNCentralJets_NMinusOne->Fill(nCentralJets,eventweight);
|
455 |
|
|
}
|
456 |
|
|
|
457 |
|
|
//Met Cut
|
458 |
|
|
pass = true;
|
459 |
|
|
for (int k=0;k<nCuts;k++) {
|
460 |
|
|
if (k != 2) {
|
461 |
|
|
pass = (pass && passCut[k]);
|
462 |
|
|
}
|
463 |
|
|
}
|
464 |
|
|
if (pass) {
|
465 |
|
|
fMetPtHist_NMinusOne->Fill(caloMet->Pt(),eventweight);
|
466 |
|
|
}
|
467 |
|
|
|
468 |
|
|
//DeltaPhiLeptons
|
469 |
|
|
pass = true;
|
470 |
|
|
for (int k=0;k<nCuts;k++) {
|
471 |
|
|
if (k != 3) {
|
472 |
|
|
pass = (pass && passCut[k]);
|
473 |
|
|
}
|
474 |
|
|
}
|
475 |
|
|
if (pass) {
|
476 |
|
|
fDeltaPhiLeptons_NMinusOne->Fill(deltaPhiLeptons,eventweight);
|
477 |
|
|
}
|
478 |
|
|
|
479 |
|
|
//dilepton mass
|
480 |
|
|
pass = true;
|
481 |
|
|
for (int k=0;k<nCuts;k++) {
|
482 |
|
|
if (k != 4)
|
483 |
|
|
pass = (pass && passCut[k]);
|
484 |
|
|
}
|
485 |
|
|
if (pass) {
|
486 |
|
|
fDileptonMass_NMinusOne->Fill(dilepton->Mass(),eventweight);
|
487 |
|
|
}
|
488 |
|
|
|
489 |
|
|
//Lepton Pt Max
|
490 |
|
|
pass = true;
|
491 |
|
|
for (int k=0;k<nCuts;k++) {
|
492 |
|
|
if (k != 5) {
|
493 |
|
|
pass = (pass && passCut[k]);
|
494 |
|
|
}
|
495 |
|
|
}
|
496 |
|
|
if (pass) {
|
497 |
|
|
fLeptonPtMax_NMinusOne->Fill(leptons->At(0)->Pt(),eventweight);
|
498 |
|
|
}
|
499 |
|
|
|
500 |
|
|
//Lepton Pt Min
|
501 |
|
|
pass = true;
|
502 |
|
|
for (int k=0;k<nCuts;k++) {
|
503 |
|
|
if (k != 6) {
|
504 |
|
|
pass = (pass && passCut[k]);
|
505 |
|
|
}
|
506 |
|
|
}
|
507 |
|
|
if (pass) {
|
508 |
|
|
fLeptonPtMin_NMinusOne->Fill(leptons->At(1)->Pt(),eventweight);
|
509 |
|
|
}
|
510 |
|
|
|
511 |
|
|
//NDirtyMuons
|
512 |
|
|
pass = true;
|
513 |
|
|
for (int k=0;k<nCuts;k++) {
|
514 |
|
|
if (k != 7)
|
515 |
|
|
pass = (pass && passCut[k]);
|
516 |
|
|
}
|
517 |
|
|
if (pass) {
|
518 |
|
|
fNDirtyMuonsHist_NMinusOne->Fill(DirtyMuons->GetEntries(),eventweight);
|
519 |
|
|
}
|
520 |
|
|
|
521 |
|
|
//NCleanExtraTracks
|
522 |
|
|
pass = true;
|
523 |
|
|
for (int k=0;k<nCuts;k++) {
|
524 |
|
|
if (k != 8)
|
525 |
|
|
pass = (pass && passCut[k]);
|
526 |
|
|
}
|
527 |
|
|
if (pass) {
|
528 |
|
|
fNCleanExtraTracksHist_NMinusOne->Fill(CleanExtraTracks->GetEntries(),
|
529 |
|
|
eventweight);
|
530 |
|
|
}
|
531 |
|
|
|
532 |
|
|
//*********************************************************************************************
|
533 |
|
|
//Plots after all Cuts
|
534 |
|
|
//*********************************************************************************************
|
535 |
|
|
if (passAllCuts) {
|
536 |
|
|
fMinDeltaPhiLeptonMet_afterCuts->Fill(minDeltaPhiMetLepton,eventweight);
|
537 |
|
|
fMtLepton1_afterCuts->Fill(mTW[0],eventweight);
|
538 |
|
|
fMtLepton2_afterCuts->Fill(mTW[1],eventweight);
|
539 |
|
|
fMtHiggs_afterCuts->Fill(mtHiggs,eventweight);
|
540 |
|
|
}
|
541 |
|
|
|
542 |
|
|
if (!fUseMCFake) {
|
543 |
|
|
delete leptons;
|
544 |
|
|
delete CleanJets;
|
545 |
|
|
delete caloMet;
|
546 |
|
|
}
|
547 |
|
|
delete dilepton;
|
548 |
|
|
delete DirtyMuons;
|
549 |
|
|
delete CleanExtraTracks;
|
550 |
|
|
}
|
551 |
|
|
|
552 |
|
|
delete FakeableObjects;
|
553 |
|
|
if (fUseMCFake) {
|
554 |
|
|
delete FakeEventHeaders;
|
555 |
|
|
}
|
556 |
|
|
return;
|
557 |
|
|
}
|
558 |
|
|
|
559 |
|
|
//--------------------------------------------------------------------------------------------------
|
560 |
|
|
void FakeLeptonExampleAnaMod::SlaveBegin()
|
561 |
|
|
{
|
562 |
|
|
// Run startup code on the computer (slave) doing the actual analysis. Here,
|
563 |
|
|
// we typically initialize histograms and other analysis objects and request
|
564 |
|
|
// branches. For this module, we request a branch of the MitTree.
|
565 |
|
|
|
566 |
|
|
ReqBranch(fTrackBranchName, fTracks);
|
567 |
|
|
ReqBranch(fMuonBranchName, fMuons);
|
568 |
|
|
|
569 |
|
|
//Create your histograms here
|
570 |
|
|
|
571 |
|
|
|
572 |
|
|
//***********************************************************************************************
|
573 |
|
|
// Before preselection
|
574 |
|
|
//***********************************************************************************************
|
575 |
|
|
AddTH1(fDileptonCharge, "hDileptonCharge", ";DileptonCharge;Number of Events", 5, -2.5, 2.5);
|
576 |
|
|
TAxis *xa = fDileptonCharge->GetXaxis();
|
577 |
|
|
xa->SetBinLabel(1,"--");
|
578 |
|
|
xa->SetBinLabel(2,"-+");
|
579 |
|
|
xa->SetBinLabel(4,"+-");
|
580 |
|
|
xa->SetBinLabel(5,"++");
|
581 |
|
|
|
582 |
|
|
//***********************************************************************************************
|
583 |
|
|
// Histograms after preselection
|
584 |
|
|
//***********************************************************************************************
|
585 |
|
|
|
586 |
|
|
AddTH1(fLeptonPtMax ,"hLeptonPtMax",";Lepton P_t Max;Number of Events",150,0.,150.);
|
587 |
|
|
AddTH1(fLeptonPtMin ,"hLeptonPtMin",";Lepton P_t Min;Number of Events",150,0.,150.);
|
588 |
|
|
AddTH1(fMetPtHist ,"hMetPtHist",";Met;Number of Events",150,0.,300.);
|
589 |
|
|
AddTH1(fDeltaPhiLeptons ,"hDeltaPhiLeptons",";#Delta#phi_{ll};Number of Events",90,0,180);
|
590 |
|
|
AddTH1(fDeltaEtaLeptons ,"hDeltaEtaLeptons",";#Delta#eta_{ll};Number of Events",100,-50.,5.0);
|
591 |
|
|
AddTH1(fDileptonMass ,"hDileptonMass",";Mass_{ll};Number of Events",150,0.,300.);
|
592 |
|
|
|
593 |
|
|
//***********************************************************************************************
|
594 |
|
|
// N-1 Histograms
|
595 |
|
|
//***********************************************************************************************
|
596 |
|
|
//All events
|
597 |
|
|
AddTH1(fLeptonPtMax_NMinusOne ,"hLeptonPtMax_NMinusOne",
|
598 |
|
|
";Lepton P_t Max;Number of Events",150,0.,150.);
|
599 |
|
|
AddTH1(fLeptonPtMin_NMinusOne ,"hLeptonPtMin_NMinusOne",
|
600 |
|
|
";Lepton P_t Min;Number of Events",150,0.,150.);
|
601 |
|
|
AddTH1(fMetPtHist_NMinusOne ,"hMetPtHist_NMinusOne",
|
602 |
|
|
";Met;Number of Events",150,0.,300.);
|
603 |
|
|
AddTH1(fMetPhiHist_NMinusOne ,"hMetPhiHist_NMinusOne",
|
604 |
|
|
";#phi;Number of Events",28,-3.5,3.5);
|
605 |
|
|
AddTH1(fMETdeltaPhilEtHist_NMinusOne ,"hMETdeltaPhilEtHist_NMinusOne",
|
606 |
|
|
";METdeltaPhilEtHist;Number of Events",150,0.,300.);
|
607 |
|
|
AddTH1(fNCentralJets_NMinusOne ,"hNCentralJets_NMinusOne",
|
608 |
|
|
";Number of Central Jets;Number of Events",6,-0.5,5.5);
|
609 |
|
|
AddTH1(fNDirtyMuonsHist_NMinusOne ,"hNDirtyMuonsHist_NMinusOne",
|
610 |
|
|
";Number of Dirty Muons; Number of Events",6,-0.5,5.5);
|
611 |
|
|
AddTH1(fNCleanExtraTracksHist_NMinusOne ,"hNCleanExtraTracksHist_NMinusOne",
|
612 |
|
|
";Number of Clean Extra Tracks; Number of Events",
|
613 |
|
|
15,-0.5,14.5);
|
614 |
|
|
AddTH1(fDeltaPhiLeptons_NMinusOne ,"hDeltaPhiLeptons_NMinusOne",
|
615 |
|
|
";#Delta#phi_{ll};Number of Events",90,0,180);
|
616 |
|
|
AddTH1(fDeltaEtaLeptons_NMinusOne ,"hDeltaEtaLeptons_NMinusOne",
|
617 |
|
|
";#Delta#eta_{ll};Number of Events",100,-5.0,5.0);
|
618 |
|
|
AddTH1(fDileptonMass_NMinusOne ,"hDileptonMass_NMinusOne",
|
619 |
|
|
";Mass_{ll};Number of Events",150,0.,300.);
|
620 |
|
|
AddTH1(fMinDeltaPhiLeptonMet_NMinusOne ,"hMinDeltaPhiLeptonMet_NMinusOne",
|
621 |
|
|
";Min #Delta#phi_{l,Met};Number of Events",90,0.,180);
|
622 |
|
|
|
623 |
|
|
|
624 |
|
|
//***********************************************************************************************
|
625 |
|
|
// After all cuts Histograms
|
626 |
|
|
//***********************************************************************************************
|
627 |
|
|
|
628 |
|
|
AddTH1(fMinDeltaPhiLeptonMet_afterCuts ,"hMinDeltaPhiLeptonMet_afterCuts",
|
629 |
|
|
";Min #Delta#phi_{l,Met};Number of Events",90,0.,180);
|
630 |
|
|
AddTH1(fMtLepton1_afterCuts ,"hMtLepton1_afterCuts",
|
631 |
|
|
";M_t (Lepton1,Met);Number of Events",100,0.,200.);
|
632 |
|
|
AddTH1(fMtLepton2_afterCuts ,"hMtLepton2_afterCuts",
|
633 |
|
|
";M_t (Lepton2,Met);Number of Events",100,0.,200.);
|
634 |
|
|
AddTH1(fMtHiggs_afterCuts ,"hMtHiggs_afterCuts",
|
635 |
|
|
";M_t (l1+l2+Met);Number of Events",150,0.,300.);
|
636 |
|
|
|
637 |
|
|
}
|
638 |
|
|
|
639 |
|
|
//--------------------------------------------------------------------------------------------------
|
640 |
|
|
void FakeLeptonExampleAnaMod::SlaveTerminate()
|
641 |
|
|
{
|
642 |
|
|
// Run finishing code on the computer (slave) that did the analysis. For this
|
643 |
|
|
// module, we dont do anything here.
|
644 |
|
|
}
|
645 |
|
|
|
646 |
|
|
//--------------------------------------------------------------------------------------------------
|
647 |
|
|
void FakeLeptonExampleAnaMod::Terminate()
|
648 |
|
|
{
|
649 |
|
|
// Run finishing code on the client computer. For this module, we dont do
|
650 |
|
|
// anything here.
|
651 |
|
|
}
|