Fork me on GitHub

source: git/modules/BTagging.cc@ 51cabe8

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