Fork me on GitHub

source: git/modules/TauTagging.cc@ 64a5c3f

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 64a5c3f was 7e227ae, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

added BitNumber to TauTagging

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