Fork me on GitHub

source: git/modules/TrackCountingTauTagging.cc@ 9867b49

Timing
Last change on this file since 9867b49 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

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