Fork me on GitHub

source: git/modules/TauTagging.cc@ fa33983

ImprovedOutputFile Timing dual_readout llp
Last change on this file since fa33983 was 50d7259, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

corrected small bug in TauTaggingPartonClassifier responsible for auto-vetoing taus with explicit W daughter, moved TauTaggingPartonClassifier definition from TauTagging.cc to TauTagging.h

  • Property mode set to 100644
File size: 7.1 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[d7d2da3]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);
[50d7259]63 Candidate *daughter1 = 0;
64 Candidate *daughter2 = 0;
65
[d7d2da3]66 const TLorentzVector &momentum = tau->Momentum;
[50d7259]67 Int_t pdgCode, i, j;
[d7d2da3]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
[6c565cc]76 if(tau->D1 >= fParticleInputArray->GetEntriesFast() ||
77 tau->D2 >= fParticleInputArray->GetEntriesFast())
78 {
79 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
80 }
81
[d7d2da3]82 for(i = tau->D1; i <= tau->D2; ++i)
83 {
[50d7259]84 daughter1 = static_cast<Candidate *>(fParticleInputArray->At(i));
85 pdgCode = TMath::Abs(daughter1->PID);
86 if(pdgCode == 11 || pdgCode == 13 || pdgCode == 15) 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 }
[d7d2da3]98 }
99
100 return 0;
101}
102
103//------------------------------------------------------------------------------
104
105TauTagging::TauTagging() :
106 fClassifier(0), fFilter(0),
107 fItPartonInputArray(0), fItJetInputArray(0)
108{
109}
110
111//------------------------------------------------------------------------------
112
113TauTagging::~TauTagging()
114{
115}
116
117//------------------------------------------------------------------------------
118
119void TauTagging::Init()
120{
121 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
122 ExRootConfParam param;
123 DelphesFormula *formula;
124 Int_t i, size;
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{
[f229d8a]191 Candidate *jet, *tau, *daughter;
192 TLorentzVector tauMomentum;
[d7d2da3]193 Double_t pt, eta, phi;
194 TObjArray *tauArray;
195 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
196 DelphesFormula *formula;
[f229d8a]197 Int_t pdgCode, charge, i;
[d7d2da3]198
199 // select taus
200 fFilter->Reset();
201 tauArray = fFilter->GetSubArray(fClassifier, 0);
202
203 if(tauArray == 0) return;
204
205 TIter itTauArray(tauArray);
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
218 // loop over all input taus
219 itTauArray.Reset();
220 while((tau = static_cast<Candidate *>(itTauArray.Next())))
221 {
[7eccbe7]222 if(tau->D1 < 0) continue;
223
224 if(tau->D1 >= fParticleInputArray->GetEntriesFast() ||
225 tau->D2 >= fParticleInputArray->GetEntriesFast())
226 {
227 throw runtime_error("tau's daughter index is greater than the ParticleInputArray size");
228 }
229
[f229d8a]230 tauMomentum.SetPxPyPzE(0.0, 0.0, 0.0, 0.0);
[7eccbe7]231
[f229d8a]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 }
[7eccbe7]238
[f229d8a]239 if(jetMomentum.DeltaR(tauMomentum) <= fDeltaR)
[d7d2da3]240 {
241 pdgCode = 15;
242 charge = tau->Charge;
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 jet->TauTag = gRandom->Uniform() <= formula->Eval(pt, eta);
255 // set tau charge
256 jet->Charge = charge;
257 }
258}
259
260//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.