Fork me on GitHub

source: git/modules/BTagging.cc@ b4786d3

Last change on this file since b4786d3 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

  • Property mode set to 100644
File size: 4.8 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/** \class BTagging
20 *
21 * Determines origin of jet,
22 * applies b-tagging efficiency (miss identification rate) formulas
23 * and sets b-tagging flags
24 *
25 * \author P. Demin - UCL, Louvain-la-Neuve
26 *
27 */
28
29#include "modules/BTagging.h"
30
31#include "classes/DelphesClasses.h"
32#include "classes/DelphesFactory.h"
33#include "classes/DelphesFormula.h"
34
35#include "TDatabasePDG.h"
36#include "TFormula.h"
37#include "TLorentzVector.h"
38#include "TMath.h"
39#include "TObjArray.h"
40#include "TRandom3.h"
41#include "TString.h"
42
43#include <algorithm>
44#include <iostream>
45#include <sstream>
46#include <stdexcept>
47
48using namespace std;
49
50//------------------------------------------------------------------------------
51
52BTagging::BTagging() :
53 fItJetInputArray(0)
54{
55}
56
57//------------------------------------------------------------------------------
58
59BTagging::~BTagging()
60{
61}
62
63//------------------------------------------------------------------------------
64
65void BTagging::Init()
66{
67 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
68 ExRootConfParam param;
69 DelphesFormula *formula;
70 Int_t i, size;
71
72 fBitNumber = GetInt("BitNumber", 0);
73
74 // read efficiency formulas
75 param = GetParam("EfficiencyFormula");
76 size = param.GetSize();
77
78 fEfficiencyMap.clear();
79 for(i = 0; i < size / 2; ++i)
80 {
81 formula = new DelphesFormula;
82 formula->Compile(param[i * 2 + 1].GetString());
83
84 fEfficiencyMap[param[i * 2].GetInt()] = formula;
85 }
86
87 // set default efficiency formula
88 itEfficiencyMap = fEfficiencyMap.find(0);
89 if(itEfficiencyMap == fEfficiencyMap.end())
90 {
91 formula = new DelphesFormula;
92 formula->Compile("0.0");
93
94 fEfficiencyMap[0] = formula;
95 }
96
97 // import input array(s)
98
99 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
100 fItJetInputArray = fJetInputArray->MakeIterator();
101}
102
103//------------------------------------------------------------------------------
104
105void BTagging::Finish()
106{
107 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
108 DelphesFormula *formula;
109
110 if(fItJetInputArray) delete fItJetInputArray;
111
112 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
113 {
114 formula = itEfficiencyMap->second;
115 if(formula) delete formula;
116 }
117}
118
119//------------------------------------------------------------------------------
120
121void BTagging::Process()
122{
123 Candidate *jet;
124 Double_t pt, eta, phi, e;
125 map<Int_t, DelphesFormula *>::iterator itEfficiencyMap;
126 DelphesFormula *formula;
127
128 // loop over all input jets
129 fItJetInputArray->Reset();
130 while((jet = static_cast<Candidate *>(fItJetInputArray->Next())))
131 {
132 const TLorentzVector &jetMomentum = jet->Momentum;
133 eta = jetMomentum.Eta();
134 phi = jetMomentum.Phi();
135 pt = jetMomentum.Pt();
136 e = jetMomentum.E();
137
138 // find an efficiency formula
139 itEfficiencyMap = fEfficiencyMap.find(jet->Flavor);
140 if(itEfficiencyMap == fEfficiencyMap.end())
141 {
142 itEfficiencyMap = fEfficiencyMap.find(0);
143 }
144 formula = itEfficiencyMap->second;
145
146 // apply an efficiency formula
147 jet->BTag |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber;
148
149 // find an efficiency formula for algo flavor definition
150 itEfficiencyMap = fEfficiencyMap.find(jet->FlavorAlgo);
151 if(itEfficiencyMap == fEfficiencyMap.end())
152 {
153 itEfficiencyMap = fEfficiencyMap.find(0);
154 }
155 formula = itEfficiencyMap->second;
156
157 // apply an efficiency formula
158 jet->BTagAlgo |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber;
159
160 // find an efficiency formula for phys flavor definition
161 itEfficiencyMap = fEfficiencyMap.find(jet->FlavorPhys);
162 if(itEfficiencyMap == fEfficiencyMap.end())
163 {
164 itEfficiencyMap = fEfficiencyMap.find(0);
165 }
166 formula = itEfficiencyMap->second;
167
168 // apply an efficiency formula
169 jet->BTagPhys |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber;
170 }
171}
172
173//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.