Fork me on GitHub

source: git/modules/TauTagging.cc@ 4b3df0d

Last change on this file since 4b3df0d 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: 7.3 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 else 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
99 return 0;
100}
101
102//------------------------------------------------------------------------------
103
104TauTagging::TauTagging() :
105 fClassifier(0), fFilter(0),
106 fItPartonInputArray(0), fItJetInputArray(0)
107{
108}
109
110//------------------------------------------------------------------------------
111
112TauTagging::~TauTagging()
113{
114}
115
116//------------------------------------------------------------------------------
117
118void TauTagging::Init()
119{
120 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
121 ExRootConfParam param;
122 DelphesFormula *formula;
123 Int_t i, size;
124
125 fBitNumber = GetInt("BitNumber", 0);
126
127 fDeltaR = GetDouble("DeltaR", 0.5);
128
129 // read efficiency formulas
130 param = GetParam("EfficiencyFormula");
131 size = param.GetSize();
132
133 fEfficiencyMap.clear();
134 for(i = 0; i < size / 2; ++i)
135 {
136 formula = new DelphesFormula;
137 formula->Compile(param[i * 2 + 1].GetString());
138
139 fEfficiencyMap[param[i * 2].GetInt()] = formula;
140 }
141
142 // set default efficiency formula
143 itEfficiencyMap = fEfficiencyMap.find(0);
144 if(itEfficiencyMap == fEfficiencyMap.end())
145 {
146 formula = new DelphesFormula;
147 formula->Compile("0.0");
148
149 fEfficiencyMap[0] = formula;
150 }
151
152 // import input array(s)
153
154 fParticleInputArray = ImportArray(GetString("ParticleInputArray", "Delphes/allParticles"));
155
156 fClassifier = new TauTaggingPartonClassifier(fParticleInputArray);
157 fClassifier->fPTMin = GetDouble("TauPTMin", 1.0);
158 fClassifier->fEtaMax = GetDouble("TauEtaMax", 2.5);
159
160 fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
161 fItPartonInputArray = fPartonInputArray->MakeIterator();
162
163 fFilter = new ExRootFilter(fPartonInputArray);
164
165 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
166 fItJetInputArray = fJetInputArray->MakeIterator();
167}
168
169//------------------------------------------------------------------------------
170
171void TauTagging::Finish()
172{
173 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
174 DelphesFormula *formula;
175
176 if(fFilter) delete fFilter;
177 if(fClassifier) delete fClassifier;
178 if(fItJetInputArray) delete fItJetInputArray;
179 if(fItPartonInputArray) delete fItPartonInputArray;
180
181 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
182 {
183 formula = itEfficiencyMap->second;
184 if(formula) delete formula;
185 }
186}
187
188//------------------------------------------------------------------------------
189
190void TauTagging::Process()
191{
192 Candidate *jet, *tau, *daughter;
193 TLorentzVector tauMomentum;
194 Double_t pt, eta, phi, e, eff;
195 TObjArray *tauArray;
196 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
197 DelphesFormula *formula;
198 Int_t pdgCode, charge, i;
199
200 // select taus
201 fFilter->Reset();
202 tauArray = fFilter->GetSubArray(fClassifier, 0);
203
204 // loop over all input jets
205 fItJetInputArray->Reset();
206 while((jet = static_cast<Candidate *>(fItJetInputArray->Next())))
207 {
208 const TLorentzVector &jetMomentum = jet->Momentum;
209 pdgCode = 0;
210 charge = gRandom->Uniform() > 0.5 ? 1 : -1;
211 eta = jetMomentum.Eta();
212 phi = jetMomentum.Phi();
213 pt = jetMomentum.Pt();
214 e = jetMomentum.E();
215
216 // loop over all input taus
217 if(tauArray)
218 {
219 TIter itTauArray(tauArray);
220 while((tau = static_cast<Candidate *>(itTauArray.Next())))
221 {
222 if(tau->D1 < 0) continue;
223
224 if(tau->D1 >= fParticleInputArray->GetEntriesFast() || tau->D2 >= fParticleInputArray->GetEntriesFast())
225 {
226 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
227 }
228
229 tauMomentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
230
231 for(i = tau->D1; i <= tau->D2; ++i)
232 {
233 daughter = static_cast<Candidate *>(fParticleInputArray->At(i));
234 if(TMath::Abs(daughter->PID) == 16) continue;
235 tauMomentum += daughter->Momentum;
236 }
237
238 if(jetMomentum.DeltaR(tauMomentum) <= fDeltaR)
239 {
240 pdgCode = 15;
241 charge = tau->Charge;
242 }
243 }
244 }
245 // find an efficency formula
246 itEfficiencyMap = fEfficiencyMap.find(pdgCode);
247 if(itEfficiencyMap == fEfficiencyMap.end())
248 {
249 itEfficiencyMap = fEfficiencyMap.find(0);
250 }
251 formula = itEfficiencyMap->second;
252
253 // apply an efficency formula
254 eff = formula->Eval(pt, eta, phi, e);
255 jet->TauTag |= (gRandom->Uniform() <= eff) << fBitNumber;
256 jet->TauWeight = eff;
257
258 // set tau charge
259 jet->Charge = charge;
260 }
261}
262
263//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.