source: trunk/modules/MadGraphIsolatedLeptonFinder.cc@ 15

Last change on this file since 15 was 2, checked in by Pavel Demin, 16 years ago

first commit

File size: 4.3 KB
RevLine 
[2]1
2#include "modules/MadGraphIsolatedLeptonFinder.h"
3
4#include "modules/MadGraphParticleClassifier.h"
5
6#include "ExRootAnalysis/ExRootClasses.h"
7#include "ExRootAnalysis/ExRootFilter.h"
8
9#include "ExRootAnalysis/ExRootFactory.h"
10#include "ExRootAnalysis/ExRootCandidate.h"
11
12#include "TString.h"
13#include "TClonesArray.h"
14#include "TLorentzVector.h"
15
16#include <map>
17#include <set>
18#include <deque>
19
20#include <iostream>
21
22using namespace std;
23
24//------------------------------------------------------------------------------
25
26Double_t MadGraphIsolatedLeptonFinder::GetMinDeltaR(ExRootGenParticle *lepton)
27{
28 Double_t distMin = 1.0e6;
29 Double_t pt1, pt2, dist;
30 TLorentzVector vector1, vector2;
31 ExRootGenParticle *particle;
32
33 vector1.SetPxPyPzE(lepton->Px, lepton->Py, lepton->Pz, lepton->E);
34 pt1 = vector1.Pt();
35
36 fItParticle->Reset();
37 while(particle = static_cast<ExRootGenParticle*>(fItParticle->Next()))
38 {
39 if(particle->Status != 1) continue;
40
41 vector2.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
42 pt2 = vector2.Pt();
43
44 dist = (pt1 == 0.0 || pt2 == 0.0 ? 1.0e6 : vector1.DeltaR(vector2));
45
46 if(pt2 > fMinPT && dist < distMin)
47 {
48 distMin = dist;
49 }
50 }
51
52 return distMin;
53}
54
55//------------------------------------------------------------------------------
56
57MadGraphIsolatedLeptonFinder::MadGraphIsolatedLeptonFinder() :
58 fItParticle(0)
59{
60}
61
62//------------------------------------------------------------------------------
63
64MadGraphIsolatedLeptonFinder::~MadGraphIsolatedLeptonFinder()
65{
66 if(fItParticle) delete fItParticle;
67}
68
69//------------------------------------------------------------------------------
70
71void MadGraphIsolatedLeptonFinder::Init()
72{
73 TString className;
74 ExRootConfParam param, classParticles;
75
76 Int_t i, j, pid, sizeParam, sizeParticles;
77
78 // import ROOT tree branch
79
80 fBranchParticle = UseBranch("GenParticle");
81
82 fItParticle = fBranchParticle->MakeIterator();
83
84 // create classifier and filter
85
86 fMinPT = GetDouble("MinPT", 1.0);
87 fMinDeltaR = GetDouble("MinDR", 0.1);
88
89 fClassifier = new MadGraphParticleClassifier();
90 fFilter = new ExRootFilter(fBranchParticle);
91
92 // read particle classes from configuration file and setup classifier
93
94 param = GetParam("ClassParticles");
95 sizeParam = param.GetSize();
96
97 for(i = 0; i < sizeParam/2; ++i)
98 {
99 className = param[i*2].GetString();
100 classParticles = param[i*2 + 1];
101 sizeParticles = classParticles.GetSize();
102
103 for(j = 0; j < sizeParticles; ++j)
104 {
105 pid = classParticles[j].GetInt();
106 fClassifier->InsertClassPID(className, pid);
107 }
108 }
109
110 fClassifier->SetExtendable(kFALSE);
111
112 // create output arrays
113
114 fOutputArray = ExportArray("candidates");
115}
116
117//------------------------------------------------------------------------------
118
119void MadGraphIsolatedLeptonFinder::Finish()
120{
121 if(fFilter) delete fFilter;
122 if(fClassifier) delete fClassifier;
123}
124
125//------------------------------------------------------------------------------
126
127void MadGraphIsolatedLeptonFinder::Process()
128{
129 TObjArray *array = 0;
130 ExRootGenParticle *particle = 0;
131 ExRootCandidate *candidate = 0;
132 ExRootFactory *factory = GetFactory();
133
134 Int_t category;
135 TString name;
136 TLorentzVector momentum;
137
138 fFilter->Reset();
139
140 // make filter classify particles and fill all subarrays
141 // at this point classifier creates additional/missing classes
142 fFilter->GetSubArray(fClassifier, 0);
143
144 // loop over all classes and export class names and classified particles
145 for(category = 0; category < fClassifier->GetMaxCategories(); ++category)
146 {
147 array = fFilter->GetSubArray(fClassifier, category);
148 name = fClassifier->GetCategoryClassName(category);
149
150 if(array == 0) continue;
151
152 // sort particles by PT
153 ExRootGenParticle::fgCompare = ExRootComparePT<ExRootGenParticle>::Instance();
154 array->Sort();
155
156 TIter itArray(array);
157
158 while(particle = static_cast<ExRootGenParticle*>(itArray.Next()))
159 {
160 if(GetMinDeltaR(particle) < fMinDeltaR) continue;
161
162 momentum.SetPxPyPzE(particle->Px, particle->Py, particle->Pz, particle->E);
163
164 candidate = factory->NewCandidate();
165
166 candidate->SetP4(momentum);
167 candidate->SetName(name);
168
169 fOutputArray->Add(candidate);
170 }
171 }
172}
173
Note: See TracBrowser for help on using the repository browser.