ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/FakeMods/src/FakeLeptonExampleAnaMod.cc
Revision: 1.3
Committed: Mon Jul 13 11:27:13 2009 UTC (15 years, 9 months ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_010
Changes since 1.2: +11 -5 lines
Log Message:
Fix compiler warnings... Still this code needs further cleanup.

File Contents

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