Fork me on GitHub

source: git/modules/TauTagging.cc@ 224c929

Last change on this file since 224c929 was 9e08f27, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

added e/mu mistag rate for taus (need to update partonArray in other than FWLite readers)

  • Property mode set to 100644
File size: 8.1 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
19/** \class TauTagging
20 *
21 * Determines origin of jet,
22 * applies b-tagging efficiency (miss identification rate) formulas
23 * and sets b-tagging flags
24 *
25 * \author P. Demin - UCL, Louvain-la-Neuve
26 *
27 */
28
29#include "modules/TauTagging.h"
30
31#include "classes/DelphesClasses.h"
32#include "classes/DelphesFactory.h"
33#include "classes/DelphesFormula.h"
34
35#include "TDatabasePDG.h"
36#include "TFormula.h"
37#include "TLorentzVector.h"
38#include "TMath.h"
39#include "TObjArray.h"
40#include "TRandom3.h"
41#include "TString.h"
42
43#include <algorithm>
44#include <iostream>
45#include <sstream>
46#include <stdexcept>
47
48using namespace std;
49
50//------------------------------------------------------------------------------
51TauTaggingPartonClassifier::TauTaggingPartonClassifier(const TObjArray *array) :
52 fParticleInputArray(array)
53{
54}
55
56//------------------------------------------------------------------------------
57
58Int_t TauTaggingPartonClassifier::GetCategory(TObject *object)
59{
60 Candidate *tau = static_cast<Candidate *>(object);
61 Candidate *daughter1 = 0;
62 Candidate *daughter2 = 0;
63
64 const TLorentzVector &momentum = tau->Momentum;
65 Int_t pdgCode, i, j;
66
67 pdgCode = TMath::Abs(tau->PID);
68 if(pdgCode != 15) return -1;
69
70 if(momentum.Pt() <= fPTMin || TMath::Abs(momentum.Eta()) > fEtaMax) return -1;
71
72 if(tau->D1 < 0) return -1;
73
74 if(tau->D2 < tau->D1) return -1;
75
76 if(tau->D1 >= fParticleInputArray->GetEntriesFast() || tau->D2 >= fParticleInputArray->GetEntriesFast())
77 {
78 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
79 }
80
81 for(i = tau->D1; i <= tau->D2; ++i)
82 {
83 daughter1 = static_cast<Candidate *>(fParticleInputArray->At(i));
84 pdgCode = TMath::Abs(daughter1->PID);
85 //if(pdgCode == 11 || pdgCode == 13 || pdgCode == 15)
86 // return -1;
87 if(pdgCode == 24)
88 {
89 if(daughter1->D1 < 0) return -1;
90 for(j = daughter1->D1; j <= daughter1->D2; ++j)
91 {
92 daughter2 = static_cast<Candidate *>(fParticleInputArray->At(j));
93 pdgCode = TMath::Abs(daughter2->PID);
94 if(pdgCode == 11 || pdgCode == 13) return -1;
95 }
96 }
97 }
98 return 0;
99}
100
101//------------------------------------------------------------------------------
102
103TauTagging::TauTagging() :
104 fClassifier(0), fFilter(0),
105 fItPartonInputArray(0), fItJetInputArray(0)
106{
107}
108
109//------------------------------------------------------------------------------
110
111TauTagging::~TauTagging()
112{
113}
114
115//------------------------------------------------------------------------------
116
117void TauTagging::Init()
118{
119 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
120 ExRootConfParam param;
121 DelphesFormula *formula;
122 Int_t i, size;
123
124 fBitNumber = GetInt("BitNumber", 0);
125
126 fDeltaR = GetDouble("DeltaR", 0.5);
127
128 // read efficiency formulas
129 param = GetParam("EfficiencyFormula");
130 size = param.GetSize();
131
132 fEfficiencyMap.clear();
133 for(i = 0; i < size / 2; ++i)
134 {
135 formula = new DelphesFormula;
136 formula->Compile(param[i * 2 + 1].GetString());
137
138 fEfficiencyMap[param[i * 2].GetInt()] = formula;
139 }
140
141 // set default efficiency formula
142 itEfficiencyMap = fEfficiencyMap.find(0);
143 if(itEfficiencyMap == fEfficiencyMap.end())
144 {
145 formula = new DelphesFormula;
146 formula->Compile("0.0");
147
148 fEfficiencyMap[0] = formula;
149 }
150
151 // import input array(s)
152
153 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "Delphes/allParticles"));
154
155 fClassifier = new TauTaggingPartonClassifier(fParticleInputArray);
156 fClassifier->fPTMin = GetDouble("TauPTMin", 1.0);
157 fClassifier->fEtaMax = GetDouble("TauEtaMax", 2.5);
158
159 fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
160 fItPartonInputArray = fPartonInputArray->MakeIterator();
161
162 fFilter = new ExRootFilter(fPartonInputArray);
163
164 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
165 fItJetInputArray = fJetInputArray->MakeIterator();
166}
167
168//------------------------------------------------------------------------------
169
170void TauTagging::Finish()
171{
172 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
173 DelphesFormula *formula;
174
175 if(fFilter) delete fFilter;
176 if(fClassifier) delete fClassifier;
177 if(fItJetInputArray) delete fItJetInputArray;
178 if(fItPartonInputArray) delete fItPartonInputArray;
179
180 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
181 {
182 formula = itEfficiencyMap->second;
183 if(formula) delete formula;
184 }
185}
186
187//------------------------------------------------------------------------------
188
189void TauTagging::Process()
190{
191 Candidate *jet, *tau, *daughter, *part;
192 TLorentzVector tauMomentum;
193 Double_t pt, eta, phi, e, eff;
194 TObjArray *tauArray;
195 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
196 DelphesFormula *formula;
197 Int_t pdgCode, charge, i;
198
199 // select taus
200 fFilter->Reset();
201 tauArray = fFilter->GetSubArray(fClassifier, 0);
202
203 // loop over all input jets
204 fItJetInputArray->Reset();
205
206 while((jet = static_cast<Candidate *>(fItJetInputArray->Next())))
207 {
208
209 const TLorentzVector &jetMomentum = jet->Momentum;
210 pdgCode = 0;
211 charge = gRandom->Uniform() > 0.5 ? 1 : -1;
212 eta = jetMomentum.Eta();
213 phi = jetMomentum.Phi();
214 pt = jetMomentum.Pt();
215 e = jetMomentum.E();
216
217 // loop over all input taus
218 if(tauArray)
219 {
220 TIter itTauArray(tauArray);
221 while((tau = static_cast<Candidate *>(itTauArray.Next())))
222 {
223 if(tau->D1 < 0) continue;
224
225 if(tau->D1 >= fParticleInputArray->GetEntriesFast() || tau->D2 >= fParticleInputArray->GetEntriesFast())
226 {
227 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
228 }
229
230 tauMomentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
231
232 for(i = tau->D1; i <= tau->D2; ++i)
233 {
234 daughter = static_cast<Candidate *>(fParticleInputArray->At(i));
235 if(TMath::Abs(daughter->PID) == 16) continue;
236 tauMomentum += daughter->Momentum;
237 }
238
239 if(jetMomentum.DeltaR(tauMomentum) <= fDeltaR)
240 {
241 pdgCode = 15;
242 charge = tau->Charge;
243 }
244 }
245 }
246
247 // fake electrons and muons
248
249 if (pdgCode == 0)
250 {
251
252 Double_t drMin = fDeltaR;
253 fItPartonInputArray->Reset();
254 while((part = static_cast<Candidate *>(fItPartonInputArray->Next())))
255 {
256 if(TMath::Abs(part->PID) == 11 || TMath::Abs(part->PID) == 13)
257 {
258 tauMomentum = part->Momentum;
259 if (tauMomentum.Pt() < fClassifier->fPTMin) continue;
260 if (TMath::Abs(tauMomentum.Eta()) > fClassifier->fEtaMax) continue;
261
262 Double_t dr = jetMomentum.DeltaR(tauMomentum);
263 if( dr < drMin)
264 {
265 drMin = dr;
266 pdgCode = TMath::Abs(part->PID);
267 charge = part->Charge;
268 }
269 }
270 }
271 }
272
273 // find an efficency formula
274 itEfficiencyMap = fEfficiencyMap.find(pdgCode);
275 if(itEfficiencyMap == fEfficiencyMap.end())
276 {
277 itEfficiencyMap = fEfficiencyMap.find(0);
278 }
279 formula = itEfficiencyMap->second;
280
281 // apply an efficency formula
282 eff = formula->Eval(pt, eta, phi, e);
283 jet->TauTag |= (gRandom->Uniform() <= eff) << fBitNumber;
284 jet->TauWeight = eff;
285
286 // set tau charge
287 jet->Charge = charge;
288 }
289}
290
291//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.