ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/FakeMods/src/FakeLeptonExampleAnaMod.cc
Revision: 1.1
Committed: Tue Jun 30 10:47:17 2009 UTC (15 years, 10 months ago) by loizides
Content type: text/plain
Branch: MAIN
Log Message:
Added FakeMods.

File Contents

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