Fork me on GitHub

source: git/modules/TrackCountingTauTagging.cc@ e5ea42e

Last change on this file since e5ea42e was 7e227ae, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

added BitNumber to TauTagging

  • Property mode set to 100644
File size: 8.1 KB
Line 
1
2/** \class TrackCountingTauTagging
3 *
4 * Determines origin of jet,
5 * applies b-tagging efficiency (miss identification rate) formulas
6 * and sets b-tagging flags
7 *
8 * $Date: 2013-07-12 00:22:27 +0200 (Fri, 12 Jul 2013) $
9 * $Revision: 1217 $
10 *
11 *
12 * \author P. Demin - UCL, Louvain-la-Neuve
13 *
14 */
15
16#include "modules/TrackCountingTauTagging.h"
17
18#include "classes/DelphesClasses.h"
19#include "classes/DelphesFactory.h"
20#include "classes/DelphesFormula.h"
21
22#include "ExRootAnalysis/ExRootResult.h"
23#include "ExRootAnalysis/ExRootFilter.h"
24#include "ExRootAnalysis/ExRootClassifier.h"
25
26#include "TMath.h"
27#include "TString.h"
28#include "TFormula.h"
29#include "TRandom3.h"
30#include "TObjArray.h"
31#include "TDatabasePDG.h"
32#include "TLorentzVector.h"
33
34#include <algorithm>
35#include <stdexcept>
36#include <iostream>
37#include <sstream>
38
39using namespace std;
40
41//------------------------------------------------------------------------------
42
43class TrackCountingTauTaggingPartonClassifier : public ExRootClassifier
44{
45public:
46
47 TrackCountingTauTaggingPartonClassifier(const TObjArray *array);
48
49 Int_t GetCategory(TObject *object);
50
51 Double_t fEtaMax, fPTMin;
52
53 const TObjArray *fParticleInputArray;
54};
55
56//------------------------------------------------------------------------------
57TrackCountingTauTaggingPartonClassifier::TrackCountingTauTaggingPartonClassifier(const TObjArray *array) :
58 fParticleInputArray(array)
59{
60}
61
62//------------------------------------------------------------------------------
63
64Int_t TrackCountingTauTaggingPartonClassifier::GetCategory(TObject *object)
65{
66 Candidate *tau = static_cast<Candidate *>(object);
67 Candidate *daughter1 = 0;
68 Candidate *daughter2 = 0;
69
70 const TLorentzVector &momentum = tau->Momentum;
71 Int_t pdgCode, i, j;
72
73 pdgCode = TMath::Abs(tau->PID);
74 if(pdgCode != 15) return -1;
75
76 if(momentum.Pt() <= fPTMin || TMath::Abs(momentum.Eta()) > fEtaMax) return -1;
77
78 if(tau->D1 < 0) return -1;
79
80 if(tau->D2 < tau->D1) return -1;
81
82 if(tau->D1 >= fParticleInputArray->GetEntriesFast() ||
83 tau->D2 >= fParticleInputArray->GetEntriesFast())
84 {
85 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
86 }
87
88 for(i = tau->D1; i <= tau->D2; ++i)
89 {
90 daughter1 = static_cast<Candidate *>(fParticleInputArray->At(i));
91 pdgCode = TMath::Abs(daughter1->PID);
92 if(pdgCode == 11 || pdgCode == 13 || pdgCode == 15) return -1;
93 else if(pdgCode == 24)
94 {
95 if(daughter1->D1 < 0) return -1;
96 for(j = daughter1->D1; j <= daughter1->D2; ++j)
97 {
98 daughter2 = static_cast<Candidate*>(fParticleInputArray->At(j));
99 pdgCode = TMath::Abs(daughter2->PID);
100 if(pdgCode == 11 || pdgCode == 13) return -1;
101 }
102
103 }
104 }
105
106 return 0;
107}
108
109//------------------------------------------------------------------------------
110
111TrackCountingTauTagging::TrackCountingTauTagging() :
112 fClassifier(0), fFilter(0),
113 fItPartonInputArray(0), fItTrackInputArray(0), fItJetInputArray(0)
114{
115}
116
117//------------------------------------------------------------------------------
118
119TrackCountingTauTagging::~TrackCountingTauTagging()
120{
121}
122
123//------------------------------------------------------------------------------
124
125void TrackCountingTauTagging::Init()
126{
127 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
128 ExRootConfParam param;
129 DelphesFormula *formula;
130 Int_t i, size;
131
132 fBitNumber = GetInt("BitNumber", 0);
133
134 fDeltaR = GetDouble("DeltaR", 0.5);
135 fDeltaRTrack = GetDouble("DeltaRTrack", 0.2);
136 fTrackPTMin = GetDouble("TrackPTMin", 1.0);
137
138 // read efficiency formulas
139 param = GetParam("EfficiencyFormula");
140 size = param.GetSize();
141
142 fEfficiencyMap.clear();
143 for(i = 0; i < size/2; ++i)
144 {
145 formula = new DelphesFormula;
146 formula->Compile(param[i*2 + 1].GetString());
147
148 fEfficiencyMap[param[i*2].GetInt()] = formula;
149 }
150
151 // set default efficiency formula
152 itEfficiencyMap = fEfficiencyMap.find(0);
153 if(itEfficiencyMap == fEfficiencyMap.end())
154 {
155 formula = new DelphesFormula;
156 formula->Compile("0.0");
157
158 fEfficiencyMap[0] = formula;
159 }
160
161 // import input array(s)
162
163 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "Delphes/allParticles"));
164
165 fClassifier = new TrackCountingTauTaggingPartonClassifier(fParticleInputArray);
166 fClassifier->fPTMin = GetDouble("TauPTMin", 1.0);
167 fClassifier->fEtaMax = GetDouble("TauEtaMax", 2.5);
168
169 fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
170 fItPartonInputArray = fPartonInputArray->MakeIterator();
171
172 fTrackInputArray = ImportArray(GetString("TrackInputArray", "TrackMerger/tracks"));
173 fItTrackInputArray = fTrackInputArray->MakeIterator();
174
175 fFilter = new ExRootFilter(fPartonInputArray);
176
177 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
178 fItJetInputArray = fJetInputArray->MakeIterator();
179}
180
181//------------------------------------------------------------------------------
182
183void TrackCountingTauTagging::Finish()
184{
185 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
186 DelphesFormula *formula;
187
188 if(fFilter) delete fFilter;
189 if(fClassifier) delete fClassifier;
190 if(fItJetInputArray) delete fItJetInputArray;
191 if(fItTrackInputArray) delete fItTrackInputArray;
192 if(fItPartonInputArray) delete fItPartonInputArray;
193
194 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
195 {
196 formula = itEfficiencyMap->second;
197 if(formula) delete formula;
198 }
199}
200
201//------------------------------------------------------------------------------
202
203void TrackCountingTauTagging::Process()
204{
205 Candidate *jet, *tau, *track, *daughter;
206 TLorentzVector tauMomentum;
207 Double_t pt, eta, phi, e;
208 TObjArray *tauArray;
209 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
210 DelphesFormula *formula;
211 Int_t pdgCode, charge, i, identifier;
212
213 // select taus
214 fFilter->Reset();
215 tauArray = fFilter->GetSubArray(fClassifier, 0);
216
217 if(tauArray == 0) return;
218
219 TIter itTauArray(tauArray);
220
221 // loop over all input jets
222 fItJetInputArray->Reset();
223 while((jet = static_cast<Candidate *>(fItJetInputArray->Next())))
224 {
225 identifier = 0;
226 const TLorentzVector &jetMomentum = jet->Momentum;
227 pdgCode = 0;
228 charge = 0;
229 eta = jetMomentum.Eta();
230 phi = jetMomentum.Phi();
231 pt = jetMomentum.Pt();
232 e = jetMomentum.E();
233
234
235// loop over all input tracks
236 fItTrackInputArray->Reset();
237 while((track = static_cast<Candidate *>(fItTrackInputArray->Next())))
238 {
239 if((track->Momentum).Pt() < fTrackPTMin) continue;
240 if(jetMomentum.DeltaR(track->Momentum) <= fDeltaRTrack) {
241 identifier -= 1;
242 charge += track->Charge;
243 }
244 }
245
246 // loop over all input taus
247 itTauArray.Reset();
248 bool matchedTau = false;
249 while((tau = static_cast<Candidate *>(itTauArray.Next())))
250 {
251 if(tau->D1 < 0) continue;
252
253 if(tau->D1 >= fParticleInputArray->GetEntriesFast() ||
254 tau->D2 >= fParticleInputArray->GetEntriesFast())
255 {
256 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
257 }
258
259 tauMomentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
260
261 for(i = tau->D1; i <= tau->D2; ++i)
262 {
263 daughter = static_cast<Candidate *>(fParticleInputArray->At(i));
264 if(TMath::Abs(daughter->PID) == 16) continue;
265 tauMomentum += daughter->Momentum;
266 }
267
268 if(jetMomentum.DeltaR(tauMomentum) <= fDeltaR)
269 {
270 matchedTau = true;
271 pdgCode = 15;
272 }
273 }
274 if(matchedTau)
275 identifier *= -1;
276 // find an efficency formula
277 // If the identifier is larger than 2, set it to 2 (multiprong requires at least 2 tracks)
278 if (identifier > 2)
279 identifier = 2;
280 else if (identifier < -2)
281 identifier = -2;
282
283
284 itEfficiencyMap = fEfficiencyMap.find(identifier);
285 if(itEfficiencyMap == fEfficiencyMap.end())
286 {
287 itEfficiencyMap = fEfficiencyMap.find(0);
288 }
289 formula = itEfficiencyMap->second;
290
291 // apply an efficency formula
292
293 // apply an efficency formula
294 jet->TauTag |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber;
295
296
297 // set tau charge
298 jet->Charge = charge;
299 }
300}
301
302//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.