1 |
// $Id: GenFakesMod.cc,v 1.6 2009/08/11 11:19:40 phedex Exp $
|
2 |
|
3 |
#include "MitPhysics/FakeMods/interface/GenFakesMod.h"
|
4 |
#include "MitCommon/MathTools/interface/MathUtils.h"
|
5 |
#include "MitAna/DataUtil/interface/Debug.h"
|
6 |
#include "MitAna/DataTree/interface/ElectronCol.h"
|
7 |
#include "MitAna/DataTree/interface/MuonCol.h"
|
8 |
#include "MitAna/DataTree/interface/PhotonCol.h"
|
9 |
#include "MitAna/DataTree/interface/JetCol.h"
|
10 |
#include "MitAna/DataTree/interface/ParticleCol.h"
|
11 |
#include "MitAna/DataTree/interface/MCParticleCol.h"
|
12 |
#include "MitPhysics/Init/interface/ModNames.h"
|
13 |
#include "MitPhysics/FakeMods/interface/FakeEventHeader.h"
|
14 |
#include "MitPhysics/FakeMods/interface/FakeRate.h"
|
15 |
|
16 |
using namespace mithep;
|
17 |
|
18 |
ClassImp(mithep::GenFakesMod)
|
19 |
|
20 |
//--------------------------------------------------------------------------------------------------
|
21 |
GenFakesMod::GenFakesMod(const char *name, const char *title) :
|
22 |
BaseMod(name,title),
|
23 |
fElectronFRFilename("InputRequired"),
|
24 |
fMuonFRFilename("InputRequired"),
|
25 |
fUse2DFakeRate(false),
|
26 |
fUseFitFunction(false),
|
27 |
fElectronFRFunctionName("InputRequired"),
|
28 |
fMuonFRFunctionName("InputRequired"),
|
29 |
fElectronFRHistName("InputRequired"),
|
30 |
fMuonFRHistName("InputRequired"),
|
31 |
fCleanElectronsName(ModNames::gkCleanElectronsName),
|
32 |
fCleanMuonsName(ModNames::gkCleanMuonsName),
|
33 |
fCleanPhotonsName(ModNames::gkCleanPhotonsName),
|
34 |
fCleanJetsName(ModNames::gkCleanJetsName),
|
35 |
fMCLeptonsName(ModNames::gkMCLeptonsName),
|
36 |
fMCTausName(ModNames::gkMCTausName),
|
37 |
fElFakeableObjsName(ModNames::gkElFakeableObjsName),
|
38 |
fMuFakeableObjsName(ModNames::gkMuFakeableObjsName),
|
39 |
fFakeEventHeadersName(ModNames::gkFakeEventHeadersName)
|
40 |
{
|
41 |
// Constructor.
|
42 |
}
|
43 |
|
44 |
//--------------------------------------------------------------------------------------------------
|
45 |
void GenFakesMod::LoadFakeRate()
|
46 |
{
|
47 |
//Load FakeRate Probabilities.
|
48 |
fFakeRate = new FakeRate(fElectronFRFilename,fMuonFRFilename,fElectronFRFunctionName,
|
49 |
fMuonFRFunctionName,fElectronFRHistName,
|
50 |
fMuonFRHistName,fUse2DFakeRate, fUseFitFunction );
|
51 |
}
|
52 |
|
53 |
//--------------------------------------------------------------------------------------------------
|
54 |
void GenFakesMod::Process()
|
55 |
{
|
56 |
// Process entries of the tree.
|
57 |
|
58 |
// get input Fakeable object collections
|
59 |
const ElectronCol *ElFakeableObjs = 0;
|
60 |
if (!fElFakeableObjsName.IsNull())
|
61 |
ElFakeableObjs = GetObjThisEvt<ElectronCol>(fElFakeableObjsName);
|
62 |
const MuonCol *MuFakeableObjs = 0;
|
63 |
if (!fMuFakeableObjsName.IsNull())
|
64 |
MuFakeableObjs = GetObjThisEvt<MuonCol>(fMuFakeableObjsName);
|
65 |
|
66 |
// get input clean object collections
|
67 |
const ElectronCol *CleanElectrons = 0;
|
68 |
if (!fCleanElectronsName.IsNull())
|
69 |
CleanElectrons = GetObjThisEvt<ElectronCol>(fCleanElectronsName);
|
70 |
const MuonCol *CleanMuons = 0;
|
71 |
if (!fCleanMuonsName.IsNull())
|
72 |
CleanMuons = GetObjThisEvt<MuonCol>(fCleanMuonsName);
|
73 |
const PhotonCol *CleanPhotons = 0;
|
74 |
if (!fCleanPhotonsName.IsNull())
|
75 |
CleanPhotons = GetObjThisEvt<PhotonCol>(fCleanPhotonsName);
|
76 |
const JetCol *CleanJets = 0;
|
77 |
if (!fCleanJetsName.IsNull())
|
78 |
CleanJets = GetObjThisEvt<JetCol>(fCleanJetsName);
|
79 |
mithep::ParticleOArr *CleanLeptons = dynamic_cast<mithep::ParticleOArr*>
|
80 |
(FindObjThisEvt(ModNames::gkMergedLeptonsName));
|
81 |
|
82 |
//get monte carlo collections
|
83 |
const MCParticleCol *GenLeptons = 0;
|
84 |
if (!fMCLeptonsName.IsNull())
|
85 |
GenLeptons = GetObjThisEvt<MCParticleCol>(fMCLeptonsName);
|
86 |
const MCParticleCol *GenTaus = 0;
|
87 |
if (!fMCTausName.IsNull())
|
88 |
GenTaus = GetObjThisEvt<MCParticleCol>(fMCTausName);
|
89 |
ObjArray<MCParticle> *GenLeptonsAndTaus = new ObjArray<MCParticle>;
|
90 |
if (GenLeptons) {
|
91 |
for (UInt_t i=0; i<GenLeptons->GetEntries(); i++)
|
92 |
GenLeptonsAndTaus->Add(GenLeptons->At(i));
|
93 |
}
|
94 |
if (GenTaus) {
|
95 |
for (UInt_t i=0; i<GenTaus->GetEntries(); i++)
|
96 |
GenLeptonsAndTaus->Add(GenTaus->At(i));
|
97 |
}
|
98 |
|
99 |
// create final output collection
|
100 |
ObjArray <FakeEventHeader> *FakeEventHeaders = new ObjArray <FakeEventHeader> ;
|
101 |
FakeEventHeaders->SetOwner(kTRUE);
|
102 |
|
103 |
//initialize with one fake event containing no fake objects and all jets.
|
104 |
FakeEventHeader *initialFakeEvent = new FakeEventHeader();
|
105 |
for (UInt_t j=0;j<CleanJets->GetEntries();j++)
|
106 |
initialFakeEvent->AddJets(CleanJets->At(j));
|
107 |
|
108 |
FakeEventHeaders->AddOwned(initialFakeEvent);
|
109 |
|
110 |
// *****************************************************************************************
|
111 |
// Fake into Muons
|
112 |
// Loop through all Muon Fakeable objects and consider the fake possibility.
|
113 |
// *****************************************************************************************
|
114 |
for (UInt_t n = 0; n < MuFakeableObjs->GetEntries(); ++n) {
|
115 |
|
116 |
//make temporary fake event headers array
|
117 |
ObjArray <FakeEventHeader> *tempFakeEventHeaders = new ObjArray <FakeEventHeader> ;
|
118 |
tempFakeEventHeaders->SetOwner(kTRUE);
|
119 |
|
120 |
//loop over all fake events generated so far - and perform an additional fake if necessary
|
121 |
for (UInt_t i=0; i<FakeEventHeaders->GetEntries();i++) {
|
122 |
|
123 |
// *****************************************************************************************
|
124 |
// Determine if the fakeable object was a clean lepton
|
125 |
// *****************************************************************************************
|
126 |
Bool_t isCleanLepton = false;
|
127 |
for (UInt_t j = 0; j < CleanLeptons->GetEntries() ; j++) {
|
128 |
Double_t deltaR = MathUtils::DeltaR(MuFakeableObjs->At(n)->Phi(),
|
129 |
MuFakeableObjs->At(n)->Eta(),
|
130 |
CleanLeptons->At(j)->Phi(), CleanLeptons->At(j)->Eta());
|
131 |
|
132 |
if (deltaR < 0.3) {
|
133 |
isCleanLepton = true;
|
134 |
break;
|
135 |
}
|
136 |
}
|
137 |
|
138 |
// *****************************************************************************************
|
139 |
// Determine if the fakeable object was a real muon from Monte Carlo
|
140 |
// *****************************************************************************************
|
141 |
Bool_t isGenMuon = false;
|
142 |
for (UInt_t l=0; l<GenLeptonsAndTaus->GetEntries(); l++) {
|
143 |
if (MathUtils::DeltaR(MuFakeableObjs->At(n)->Phi(), MuFakeableObjs->At(n)->Eta(),
|
144 |
GenLeptonsAndTaus->At(l)->Phi(),
|
145 |
GenLeptonsAndTaus->At(l)->Eta()) < 0.3) {
|
146 |
isGenMuon = true;
|
147 |
break;
|
148 |
}
|
149 |
}
|
150 |
|
151 |
//this is used to determine the weight of the unfaked event.
|
152 |
double totalCumulativeFakeProbability = 0.0;
|
153 |
|
154 |
// *****************************************************************************************
|
155 |
// Perform Muon Fake
|
156 |
// *****************************************************************************************
|
157 |
|
158 |
//match fake to one of the jets
|
159 |
int fakeToJetMatch = -1; //index of the jet that matches to the fake
|
160 |
double minDR = 5000;
|
161 |
for (UInt_t jj=0;jj<FakeEventHeaders->At(i)->NJets();jj++) {
|
162 |
Double_t deltaR = MathUtils::DeltaR(FakeEventHeaders->At(i)->UnfakedJet(jj)->Mom(),
|
163 |
MuFakeableObjs->At(n)->Mom());
|
164 |
if (deltaR < minDR) {
|
165 |
minDR = deltaR;
|
166 |
fakeToJetMatch = jj;
|
167 |
}
|
168 |
}
|
169 |
if (!(minDR < 0.5)) {
|
170 |
fakeToJetMatch = -1;
|
171 |
}
|
172 |
|
173 |
//Obtain the muon FakeRate
|
174 |
Double_t muonFakeProb = 0.0;
|
175 |
Double_t muonFakeProbStatErrorLow = 0.0;
|
176 |
Double_t muonFakeProbStatErrorHigh = 0.0;
|
177 |
Double_t muonFakeProbSysErrorLow = 0.0;
|
178 |
Double_t muonFakeProbSysErrorHigh = 0.0;
|
179 |
Double_t muonFakeProbErrorLow = 0.0;
|
180 |
Double_t muonFakeProbErrorHigh = 0.0;
|
181 |
if(fFakeRate) {
|
182 |
muonFakeProb = fFakeRate->MuonFakeRate(MuFakeableObjs->At(n)->Et(),
|
183 |
MuFakeableObjs->At(n)->Eta(),
|
184 |
MuFakeableObjs->At(n)->Phi());
|
185 |
muonFakeProbStatErrorLow = fFakeRate->MuonFakeRateStatErrorLow(MuFakeableObjs->At(n)->Et(),
|
186 |
MuFakeableObjs->At(n)->Eta(),
|
187 |
MuFakeableObjs->At(n)->Phi());
|
188 |
muonFakeProbStatErrorHigh =
|
189 |
fFakeRate->MuonFakeRateStatErrorHigh(MuFakeableObjs->At(n)->Et(),
|
190 |
MuFakeableObjs->At(n)->Eta(),
|
191 |
MuFakeableObjs->At(n)->Phi());
|
192 |
muonFakeProbSysErrorLow = fFakeRate->MuonFakeRateSysErrorLow(MuFakeableObjs->At(n)->Et(),
|
193 |
MuFakeableObjs->At(n)->Eta(),
|
194 |
MuFakeableObjs->At(n)->Phi());
|
195 |
muonFakeProbSysErrorHigh = fFakeRate->MuonFakeRateSysErrorHigh(MuFakeableObjs->At(n)->Et(),
|
196 |
MuFakeableObjs->At(n)->Eta(),
|
197 |
MuFakeableObjs->At(n)->Phi());
|
198 |
muonFakeProbErrorLow = TMath::Sqrt(muonFakeProbStatErrorLow*muonFakeProbStatErrorLow +
|
199 |
muonFakeProbSysErrorLow*muonFakeProbSysErrorLow);
|
200 |
muonFakeProbErrorHigh = TMath::Sqrt(muonFakeProbStatErrorHigh*muonFakeProbStatErrorHigh +
|
201 |
muonFakeProbSysErrorHigh*muonFakeProbSysErrorHigh);
|
202 |
|
203 |
} else {
|
204 |
Fatal("Process()","Error: fFakeRate is a NULL pointer.");
|
205 |
}
|
206 |
|
207 |
//only fake into a muon if the fakeable object did not match to a clean lepton
|
208 |
if (!isCleanLepton) {
|
209 |
|
210 |
//create new fake event header
|
211 |
FakeEventHeader *fakeMuonEvent = new FakeEventHeader();
|
212 |
fakeMuonEvent->SetWeight(FakeEventHeaders->At(i)->Weight() * muonFakeProb);
|
213 |
Double_t weightLowError = 0;
|
214 |
Double_t weightHighError = 0;
|
215 |
if (muonFakeProb > 0) {
|
216 |
weightLowError = FakeEventHeaders->At(i)->Weight()*muonFakeProb*
|
217 |
TMath::Sqrt((FakeEventHeaders->At(i)->WeightLowError()/
|
218 |
FakeEventHeaders->At(i)->Weight())*
|
219 |
(FakeEventHeaders->At(i)->WeightLowError()/
|
220 |
FakeEventHeaders->At(i)->Weight()) +
|
221 |
(muonFakeProbErrorLow/muonFakeProb)*
|
222 |
(muonFakeProbErrorLow/muonFakeProb));
|
223 |
weightHighError = FakeEventHeaders->At(i)->Weight()*muonFakeProb*
|
224 |
TMath::Sqrt((FakeEventHeaders->At(i)->WeightHighError()/
|
225 |
FakeEventHeaders->At(i)->Weight())*
|
226 |
(FakeEventHeaders->At(i)->WeightHighError()/
|
227 |
FakeEventHeaders->At(i)->Weight()) +
|
228 |
(muonFakeProbErrorHigh/muonFakeProb)*
|
229 |
(muonFakeProbErrorHigh/muonFakeProb));
|
230 |
}
|
231 |
fakeMuonEvent->SetWeightLowError(weightLowError);
|
232 |
fakeMuonEvent->SetWeightHighError(weightHighError);
|
233 |
|
234 |
//add all previous fakes
|
235 |
for (UInt_t f=0;f<FakeEventHeaders->At(i)->FakeObjsSize();f++) {
|
236 |
fakeMuonEvent->AddFakeObject(FakeEventHeaders->At(i)->FakeObj(f));
|
237 |
}
|
238 |
//add new fake
|
239 |
fakeMuonEvent->AddFakeObject(MuFakeableObjs->At(n),kMuon,isCleanLepton,isGenMuon);
|
240 |
|
241 |
//add all previous jets except the one matching to the new fake
|
242 |
for (UInt_t jj=0;jj<FakeEventHeaders->At(i)->NJets();jj++) {
|
243 |
if ((int)jj != fakeToJetMatch)
|
244 |
fakeMuonEvent->AddJets(FakeEventHeaders->At(i)->UnfakedJet(jj));
|
245 |
}
|
246 |
|
247 |
//add fake event to the temporary fake event header array
|
248 |
tempFakeEventHeaders->AddOwned(fakeMuonEvent);
|
249 |
|
250 |
//increase cumulative fake probability
|
251 |
totalCumulativeFakeProbability += muonFakeProb;
|
252 |
}
|
253 |
|
254 |
// *****************************************************************************************
|
255 |
// Do not fake into Muon
|
256 |
// *****************************************************************************************
|
257 |
//create new fake event header
|
258 |
FakeEventHeader *notFakeEvent = new FakeEventHeader();
|
259 |
notFakeEvent->SetWeight(FakeEventHeaders->At(i)->Weight() *
|
260 |
(1-totalCumulativeFakeProbability));
|
261 |
//add previous fakes
|
262 |
for (UInt_t f=0;f<FakeEventHeaders->At(i)->FakeObjsSize();f++) {
|
263 |
notFakeEvent->AddFakeObject(FakeEventHeaders->At(i)->FakeObj(f));
|
264 |
}
|
265 |
//add previous jets
|
266 |
for (UInt_t jj=0;jj<FakeEventHeaders->At(i)->NJets();jj++) {
|
267 |
notFakeEvent->AddJets(FakeEventHeaders->At(i)->UnfakedJet(jj));
|
268 |
}
|
269 |
tempFakeEventHeaders->AddOwned(notFakeEvent);
|
270 |
|
271 |
} //loop over all current fake event headers
|
272 |
|
273 |
//replace current fake event headers with the new temporary ones.
|
274 |
delete FakeEventHeaders;
|
275 |
FakeEventHeaders = tempFakeEventHeaders;
|
276 |
} //loop over all muon fakeable objects
|
277 |
|
278 |
|
279 |
// *****************************************************************************************
|
280 |
// Fake into Electrons
|
281 |
// Loop through all Electron Fakeable objects and consider the fake possibility.
|
282 |
// *****************************************************************************************
|
283 |
for (UInt_t n = 0; n < ElFakeableObjs->GetEntries();n++) {
|
284 |
|
285 |
//make temporary fake event headers array
|
286 |
ObjArray <FakeEventHeader> *tempFakeEventHeaders = new ObjArray <FakeEventHeader> ;
|
287 |
tempFakeEventHeaders->SetOwner(kTRUE);
|
288 |
|
289 |
//loop over all fake events generated so far - and perform an additional fake if necessary
|
290 |
for (UInt_t i=0; i<FakeEventHeaders->GetEntries();i++) {
|
291 |
|
292 |
// *****************************************************************************************
|
293 |
// Determine if the fakeable object was a clean lepton
|
294 |
// *****************************************************************************************
|
295 |
Bool_t isCleanLepton = false;
|
296 |
for (UInt_t j = 0; j < CleanLeptons->GetEntries() ; j++) {
|
297 |
Double_t deltaR = MathUtils::DeltaR(ElFakeableObjs->At(n)->Phi(),
|
298 |
ElFakeableObjs->At(n)->Eta(),
|
299 |
CleanLeptons->At(j)->Phi(), CleanLeptons->At(j)->Eta());
|
300 |
|
301 |
if (deltaR < 0.3) {
|
302 |
isCleanLepton = true;
|
303 |
break;
|
304 |
}
|
305 |
}
|
306 |
|
307 |
// *****************************************************************************************
|
308 |
// Determine if the fakeable object was a real electron from Monte Carlo
|
309 |
// *****************************************************************************************
|
310 |
Bool_t isGenElectron = false;
|
311 |
for (UInt_t l=0; l<GenLeptonsAndTaus->GetEntries(); l++) {
|
312 |
if (MathUtils::DeltaR(ElFakeableObjs->At(n)->Phi(),
|
313 |
ElFakeableObjs->At(n)->Eta(),
|
314 |
GenLeptonsAndTaus->At(l)->Phi(),
|
315 |
GenLeptonsAndTaus->At(l)->Eta()) < 0.3) {
|
316 |
isGenElectron = true;
|
317 |
}
|
318 |
}
|
319 |
|
320 |
|
321 |
//this is used to determine the weight of the unfaked event.
|
322 |
double totalCumulativeFakeProbability = 0.0;
|
323 |
|
324 |
// *****************************************************************************************
|
325 |
// Determine if the current electron fakeable object already corresponds to one of the
|
326 |
// fake muons already in the FakeEventHeader, determined based on deltaR proximity.
|
327 |
// If the current electron fakeable object corresponds to one of the fake muon, then
|
328 |
// we do not allow it to fake an electron, since one denominator cannot fake two leptons.
|
329 |
// *****************************************************************************************
|
330 |
Bool_t alreadyFaked = false;
|
331 |
for (UInt_t f = 0; f < FakeEventHeaders->At(i)->FakeObjsSize() ; f++) {
|
332 |
double deltaR = MathUtils::DeltaR(FakeEventHeaders->At(i)->FakeObj(f)->Mom(),
|
333 |
ElFakeableObjs->At(n)->Mom());
|
334 |
if (deltaR < 0.3)
|
335 |
alreadyFaked = true;
|
336 |
}
|
337 |
if (!alreadyFaked) {
|
338 |
|
339 |
// *****************************************************************************************
|
340 |
// Perform Electron Fake
|
341 |
// *****************************************************************************************
|
342 |
|
343 |
//match fake to one of the jets
|
344 |
int fakeToJetMatch = -1; //index of the jet that matches to the fake
|
345 |
double minDR = 5000;
|
346 |
for (UInt_t jj=0;jj<FakeEventHeaders->At(i)->NJets();jj++) {
|
347 |
Double_t deltaR = MathUtils::DeltaR(FakeEventHeaders->At(i)->UnfakedJet(jj)->Mom(),
|
348 |
ElFakeableObjs->At(n)->Mom());
|
349 |
if (deltaR < minDR) {
|
350 |
minDR = deltaR;
|
351 |
fakeToJetMatch = jj;
|
352 |
}
|
353 |
}
|
354 |
if (!(minDR < 0.5)) {
|
355 |
fakeToJetMatch = -1;
|
356 |
}
|
357 |
|
358 |
//Obtain the electron FakeRate
|
359 |
Double_t electronFakeProb = 0.0;
|
360 |
Double_t electronFakeProbStatErrorLow = 0.0;
|
361 |
Double_t electronFakeProbStatErrorHigh = 0.0;
|
362 |
Double_t electronFakeProbSysErrorLow = 0.0;
|
363 |
Double_t electronFakeProbSysErrorHigh = 0.0;
|
364 |
Double_t electronFakeProbErrorLow = 0.0;
|
365 |
Double_t electronFakeProbErrorHigh = 0.0;
|
366 |
if(fFakeRate) {
|
367 |
electronFakeProb = fFakeRate->ElectronFakeRate(ElFakeableObjs->At(n)->Et(),
|
368 |
ElFakeableObjs->At(n)->Eta(),
|
369 |
ElFakeableObjs->At(n)->Phi());
|
370 |
electronFakeProbStatErrorLow =
|
371 |
fFakeRate->ElectronFakeRateStatErrorLow(ElFakeableObjs->At(n)->Et(),
|
372 |
ElFakeableObjs->At(n)->Eta(),
|
373 |
ElFakeableObjs->At(n)->Phi());
|
374 |
electronFakeProbStatErrorHigh =
|
375 |
fFakeRate->ElectronFakeRateStatErrorHigh(ElFakeableObjs->At(n)->Et(),
|
376 |
ElFakeableObjs->At(n)->Eta(),
|
377 |
ElFakeableObjs->At(n)->Phi());
|
378 |
electronFakeProbSysErrorLow =
|
379 |
fFakeRate->ElectronFakeRateSysErrorLow(ElFakeableObjs->At(n)->Et(),
|
380 |
ElFakeableObjs->At(n)->Eta(),
|
381 |
ElFakeableObjs->At(n)->Phi());
|
382 |
electronFakeProbSysErrorHigh =
|
383 |
fFakeRate->ElectronFakeRateSysErrorHigh(ElFakeableObjs->At(n)->Et(),
|
384 |
ElFakeableObjs->At(n)->Eta(),
|
385 |
ElFakeableObjs->At(n)->Phi());
|
386 |
electronFakeProbErrorLow =
|
387 |
TMath::Sqrt(electronFakeProbStatErrorLow*electronFakeProbStatErrorLow +
|
388 |
electronFakeProbSysErrorLow*electronFakeProbSysErrorLow);
|
389 |
electronFakeProbErrorHigh =
|
390 |
TMath::Sqrt(electronFakeProbStatErrorHigh*electronFakeProbStatErrorHigh +
|
391 |
electronFakeProbSysErrorHigh*electronFakeProbSysErrorHigh);
|
392 |
} else {
|
393 |
Fatal("Process()","Error: fFakeRate is a NULL pointer.");
|
394 |
}
|
395 |
|
396 |
//only fake into an electron if the fakeable object did not match to a clean lepton
|
397 |
if (!isCleanLepton) {
|
398 |
//create new fake event header
|
399 |
FakeEventHeader *fakeElectronEvent = new FakeEventHeader();
|
400 |
fakeElectronEvent->SetWeight(FakeEventHeaders->At(i)->Weight() * electronFakeProb);
|
401 |
Double_t weightLowError = 0;
|
402 |
Double_t weightHighError = 0;
|
403 |
if (electronFakeProb) {
|
404 |
weightLowError = FakeEventHeaders->At(i)->Weight()*electronFakeProb*
|
405 |
TMath::Sqrt((FakeEventHeaders->At(i)->WeightLowError()/
|
406 |
FakeEventHeaders->At(i)->Weight())*
|
407 |
(FakeEventHeaders->At(i)->WeightLowError()/
|
408 |
FakeEventHeaders->At(i)->Weight()) +
|
409 |
(electronFakeProbErrorLow/electronFakeProb)*
|
410 |
(electronFakeProbErrorLow/electronFakeProb));
|
411 |
weightHighError = FakeEventHeaders->At(i)->Weight()*electronFakeProb*
|
412 |
TMath::Sqrt((FakeEventHeaders->At(i)->WeightHighError()/
|
413 |
FakeEventHeaders->At(i)->Weight())*
|
414 |
(FakeEventHeaders->At(i)->WeightHighError()/
|
415 |
FakeEventHeaders->At(i)->Weight()) +
|
416 |
(electronFakeProbErrorHigh/electronFakeProb)*
|
417 |
(electronFakeProbErrorHigh/electronFakeProb));
|
418 |
}
|
419 |
fakeElectronEvent->SetWeightLowError(weightLowError);
|
420 |
fakeElectronEvent->SetWeightHighError(weightHighError);
|
421 |
|
422 |
//add previous fakes
|
423 |
for (UInt_t f=0;f<FakeEventHeaders->At(i)->FakeObjsSize();f++) {
|
424 |
fakeElectronEvent->AddFakeObject(FakeEventHeaders->At(i)->FakeObj(f));
|
425 |
}
|
426 |
//add the new fake
|
427 |
fakeElectronEvent->AddFakeObject(ElFakeableObjs->At(n),
|
428 |
kElectron,isCleanLepton,isGenElectron);
|
429 |
//add previous jets that do not match to the new fake
|
430 |
for (UInt_t jj=0;jj<FakeEventHeaders->At(i)->NJets();jj++) {
|
431 |
if ((int)jj != fakeToJetMatch)
|
432 |
fakeElectronEvent->AddJets(FakeEventHeaders->At(i)->UnfakedJet(jj));
|
433 |
}
|
434 |
|
435 |
//add fake event to the temporary fake event header array
|
436 |
tempFakeEventHeaders->AddOwned(fakeElectronEvent);
|
437 |
//increase cumulative fake probability
|
438 |
totalCumulativeFakeProbability += electronFakeProb;
|
439 |
}
|
440 |
}
|
441 |
|
442 |
// *****************************************************************************************
|
443 |
// Do not fake into anything
|
444 |
// *****************************************************************************************
|
445 |
//create new fake event header
|
446 |
FakeEventHeader *notFakeEvent = new FakeEventHeader();
|
447 |
notFakeEvent->SetWeight(FakeEventHeaders->At(i)->Weight() *
|
448 |
(1-totalCumulativeFakeProbability));
|
449 |
//add previous fakes
|
450 |
for (UInt_t f=0;f<FakeEventHeaders->At(i)->FakeObjsSize();f++) {
|
451 |
notFakeEvent->AddFakeObject(FakeEventHeaders->At(i)->FakeObj(f));
|
452 |
}
|
453 |
//add previous jets
|
454 |
for (UInt_t jj=0;jj<FakeEventHeaders->At(i)->NJets();jj++) {
|
455 |
notFakeEvent->AddJets(FakeEventHeaders->At(i)->UnfakedJet(jj));
|
456 |
}
|
457 |
tempFakeEventHeaders->AddOwned(notFakeEvent);
|
458 |
|
459 |
} //for all current fake event headers
|
460 |
|
461 |
//replace current fake event headers with the new temporary ones.
|
462 |
delete FakeEventHeaders;
|
463 |
FakeEventHeaders = tempFakeEventHeaders;
|
464 |
} //loop over all fakeable objects
|
465 |
|
466 |
// export FakeEventHeaders for other modules to use
|
467 |
FakeEventHeaders->SetName(fFakeEventHeadersName);
|
468 |
AddObjThisEvt(FakeEventHeaders);
|
469 |
|
470 |
//delete temporary collections
|
471 |
delete GenLeptonsAndTaus;
|
472 |
}
|