ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitPhysics/FakeMods/src/FakeLeptonExampleAnaMod.cc
Revision: 1.4
Committed: Mon Aug 10 16:07:26 2009 UTC (15 years, 8 months ago) by phedex
Content type: text/plain
Branch: MAIN
Changes since 1.3: +1 -3 lines
Log Message:
Use TH2DAsymError class to store asymmetric statistical and systematic errors. Do comment and style cleanup.

File Contents

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