Fork me on GitHub

source: git/modules/TauTagging.cc@ b78adf8

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b78adf8 was cab38f6, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

remove svn tags and fix formatting

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