Fork me on GitHub

source: git/modules/TauTagging.cc@ 400597a

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 400597a was 7c0fcd5, checked in by Pavel Demin <demin@…>, 10 years ago

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