Fork me on GitHub

source: git/modules/BTagging.cc@ a85a257

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

fix EOL characters in GPLv3 header

  • Property mode set to 100644
File size: 6.0 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 BTagging
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/BTagging.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 BTaggingPartonClassifier : public ExRootClassifier
62{
63public:
64
65 BTaggingPartonClassifier() {}
66
67 Int_t GetCategory(TObject *object);
68
69 Double_t fEtaMax, fPTMin;
70};
71
72//------------------------------------------------------------------------------
73
74Int_t BTaggingPartonClassifier::GetCategory(TObject *object)
75{
76 Candidate *parton = static_cast<Candidate*>(object);
77 const TLorentzVector &momentum = parton->Momentum;
78 Int_t pdgCode;
79
80 if(momentum.Pt() <= fPTMin || TMath::Abs(momentum.Eta()) > fEtaMax) return -1;
81
82 pdgCode = TMath::Abs(parton->PID);
83 if(pdgCode != 21 && pdgCode > 5) return -1;
84
85 return 0;
86}
87
88//------------------------------------------------------------------------------
89
90BTagging::BTagging() :
91 fClassifier(0), fFilter(0),
92 fItPartonInputArray(0), fItJetInputArray(0)
93{
94 fClassifier = new BTaggingPartonClassifier;
95}
96
97//------------------------------------------------------------------------------
98
99BTagging::~BTagging()
100{
101 if(fClassifier) delete fClassifier;
102}
103
104//------------------------------------------------------------------------------
105
106void BTagging::Init()
107{
108 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
109 ExRootConfParam param;
110 DelphesFormula *formula;
111 Int_t i, size;
112
113 fBitNumber = GetInt("BitNumber", 0);
114
115 fDeltaR = GetDouble("DeltaR", 0.5);
116
117 fClassifier->fPTMin = GetDouble("PartonPTMin", 1.0);
118 fClassifier->fEtaMax = GetDouble("PartonEtaMax", 2.5);
119
120 // read efficiency formulas
121 param = GetParam("EfficiencyFormula");
122 size = param.GetSize();
123
124 fEfficiencyMap.clear();
125 for(i = 0; i < size/2; ++i)
126 {
127 formula = new DelphesFormula;
128 formula->Compile(param[i*2 + 1].GetString());
129
130 fEfficiencyMap[param[i*2].GetInt()] = formula;
131 }
132
133 // set default efficiency formula
134 itEfficiencyMap = fEfficiencyMap.find(0);
135 if(itEfficiencyMap == fEfficiencyMap.end())
136 {
137 formula = new DelphesFormula;
138 formula->Compile("0.0");
139
140 fEfficiencyMap[0] = formula;
141 }
142
143 // import input array(s)
144
145 fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));
146 fItPartonInputArray = fPartonInputArray->MakeIterator();
147
148 fFilter = new ExRootFilter(fPartonInputArray);
149
150 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
151 fItJetInputArray = fJetInputArray->MakeIterator();
152}
153
154//------------------------------------------------------------------------------
155
156void BTagging::Finish()
157{
158 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
159 DelphesFormula *formula;
160
161 if(fFilter) delete fFilter;
162 if(fItJetInputArray) delete fItJetInputArray;
163 if(fItPartonInputArray) delete fItPartonInputArray;
164
165 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
166 {
167 formula = itEfficiencyMap->second;
168 if(formula) delete formula;
169 }
170}
171
172//------------------------------------------------------------------------------
173
174void BTagging::Process()
175{
176 Candidate *jet, *parton;
177 Double_t pt, eta, phi;
178 TObjArray *partonArray;
179 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
180 DelphesFormula *formula;
181 Int_t pdgCode, pdgCodeMax;
182
183 // select quark and gluons
184 fFilter->Reset();
185 partonArray = fFilter->GetSubArray(fClassifier, 0);
186
187 if(partonArray == 0) return;
188
189 TIter itPartonArray(partonArray);
190
191 // loop over all input jets
192 fItJetInputArray->Reset();
193 while((jet = static_cast<Candidate*>(fItJetInputArray->Next())))
194 {
195 const TLorentzVector &jetMomentum = jet->Momentum;
196 pdgCodeMax = -1;
197 eta = jetMomentum.Eta();
198 phi = jetMomentum.Phi();
199 pt = jetMomentum.Pt();
200
201 // loop over all input partons
202 itPartonArray.Reset();
203 while((parton = static_cast<Candidate*>(itPartonArray.Next())))
204 {
205 pdgCode = TMath::Abs(parton->PID);
206 if(pdgCode == 21) pdgCode = 0;
207 if(jetMomentum.DeltaR(parton->Momentum) <= fDeltaR)
208 {
209 if(pdgCodeMax < pdgCode) pdgCodeMax = pdgCode;
210 }
211 }
212 if(pdgCodeMax == 0) pdgCodeMax = 21;
213 if(pdgCodeMax == -1) pdgCodeMax = 0;
214
215 // find an efficency formula
216 itEfficiencyMap = fEfficiencyMap.find(pdgCodeMax);
217 if(itEfficiencyMap == fEfficiencyMap.end())
218 {
219 itEfficiencyMap = fEfficiencyMap.find(0);
220 }
221 formula = itEfficiencyMap->second;
222
223 // apply an efficency formula
224 jet->BTag |= (gRandom->Uniform() <= formula->Eval(pt, eta)) << fBitNumber;
225 }
226}
227
228//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.