Fork me on GitHub

source: git/modules/BTagging.cc@ 84edab9

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 84edab9 was 1180bc1, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

fixed BTagging without LHEF

  • 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
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 "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
53BTagging::BTagging() :
54 fItJetInputArray(0)
55{
56}
57
58//------------------------------------------------------------------------------
59
60BTagging::~BTagging()
61{
62}
63
64//------------------------------------------------------------------------------
65
66void BTagging::Init()
67{
68 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
69 ExRootConfParam param;
70 DelphesFormula *formula;
71 Int_t i, size;
72
73 fBitNumber = GetInt("BitNumber", 0);
74
75 // read efficiency formulas
76 param = GetParam("EfficiencyFormula");
77 size = param.GetSize();
78
79 fEfficiencyMap.clear();
80 for(i = 0; i < size/2; ++i)
81 {
82 formula = new DelphesFormula;
83 formula->Compile(param[i*2 + 1].GetString());
84
85 fEfficiencyMap[param[i*2].GetInt()] = formula;
86 }
87
88 // set default efficiency formula
89 itEfficiencyMap = fEfficiencyMap.find(0);
90 if(itEfficiencyMap == fEfficiencyMap.end())
91 {
92 formula = new DelphesFormula;
93 formula->Compile("0.0");
94
95 fEfficiencyMap[0] = formula;
96 }
97
98 // import input array(s)
99
100 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
101 fItJetInputArray = fJetInputArray->MakeIterator();
102}
103
104//------------------------------------------------------------------------------
105
106void BTagging::Finish()
107{
108 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
109 DelphesFormula *formula;
110
111 if(fItJetInputArray) delete fItJetInputArray;
112
113 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap)
114 {
115 formula = itEfficiencyMap->second;
116 if(formula) delete formula;
117 }
118}
119
120//------------------------------------------------------------------------------
121
122void BTagging::Process()
123{
124 Candidate *jet;
125 Double_t pt, eta, phi, e;
126 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap;
127 DelphesFormula *formula;
128
129 // loop over all input jets
130 fItJetInputArray->Reset();
131 while((jet = static_cast<Candidate*>(fItJetInputArray->Next())))
132 {
133 const TLorentzVector &jetMomentum = jet->Momentum;
134 eta = jetMomentum.Eta();
135 phi = jetMomentum.Phi();
136 pt = jetMomentum.Pt();
137 e = jetMomentum.E();
138
139 // find an efficiency formula
140 itEfficiencyMap = fEfficiencyMap.find(jet->Flavor);
141 if(itEfficiencyMap == fEfficiencyMap.end())
142 {
143 itEfficiencyMap = fEfficiencyMap.find(0);
144 }
145 formula = itEfficiencyMap->second;
146
147 // apply an efficiency formula
148 jet->BTag |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber;
149
150 // find an efficiency formula for algo flavor definition
151 itEfficiencyMap = fEfficiencyMap.find(jet->FlavorAlgo);
152 if(itEfficiencyMap == fEfficiencyMap.end())
153 {
154 itEfficiencyMap = fEfficiencyMap.find(0);
155 }
156 formula = itEfficiencyMap->second;
157
158 // apply an efficiency formula
159 jet->BTagAlgo |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber;
160
161 // find an efficiency formula for phys flavor definition
162 itEfficiencyMap = fEfficiencyMap.find(jet->FlavorPhys);
163 if(itEfficiencyMap == fEfficiencyMap.end())
164 {
165 itEfficiencyMap = fEfficiencyMap.find(0);
166 }
167 formula = itEfficiencyMap->second;
168
169 // apply an efficiency formula
170 jet->BTagPhys |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber;
171 }
172}
173
174//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.