Fork me on GitHub

source: git/modules/TrackCountingBTagging.cc@ 9e3f2fb

Timing
Last change on this file since 9e3f2fb 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.4 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 TrackCountingBTagging
20 *
21 * b-tagging algorithm based on counting tracks with large impact parameter
22 *
23 * \author M. Selvaggi - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/TrackCountingBTagging.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "TFormula.h"
34#include "TLorentzVector.h"
35#include "TMath.h"
36#include "TObjArray.h"
37#include "TRandom3.h"
38#include "TString.h"
39
40#include <algorithm>
41#include <iostream>
42#include <sstream>
43#include <stdexcept>
44
45using namespace std;
46
47//------------------------------------------------------------------------------
48
49TrackCountingBTagging::TrackCountingBTagging() :
50 fItTrackInputArray(0), fItJetInputArray(0)
51{
52}
53
54//------------------------------------------------------------------------------
55
56TrackCountingBTagging::~TrackCountingBTagging()
57{
58}
59
60//------------------------------------------------------------------------------
61
62void TrackCountingBTagging::Init()
63{
64 fBitNumber = GetInt("BitNumber", 0);
65
66 fPtMin = GetDouble("TrackPtMin", 1.0);
67 fDeltaR = GetDouble("DeltaR", 0.3);
68 fIPmax = GetDouble("TrackIPMax", 2.0);
69
70 fSigMin = GetDouble("SigMin", 6.5);
71 fNtracks = GetInt("Ntracks", 3);
72
73 fUse3D = GetBool("Use3D", false);
74
75 // import input array(s)
76
77 fTrackInputArray = ImportArray(GetString("TrackInputArray", "Calorimeter/eflowTracks"));
78 fItTrackInputArray = fTrackInputArray->MakeIterator();
79
80 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
81 fItJetInputArray = fJetInputArray->MakeIterator();
82}
83
84//------------------------------------------------------------------------------
85
86void TrackCountingBTagging::Finish()
87{
88 if(fItTrackInputArray) delete fItTrackInputArray;
89 if(fItJetInputArray) delete fItJetInputArray;
90}
91
92//------------------------------------------------------------------------------
93
94void TrackCountingBTagging::Process()
95{
96 Candidate *jet, *track;
97
98 Double_t jpx, jpy, jpz;
99 Double_t dr, tpt;
100 Double_t xd, yd, zd, d0, dd0, dz, ddz, sip;
101
102 Int_t sign;
103
104 Int_t count;
105
106 // loop over all input jets
107 fItJetInputArray->Reset();
108 while((jet = static_cast<Candidate *>(fItJetInputArray->Next())))
109 {
110 const TLorentzVector &jetMomentum = jet->Momentum;
111 jpx = jetMomentum.Px();
112 jpy = jetMomentum.Py();
113 jpz = jetMomentum.Pz();
114
115 // loop over all input tracks
116 fItTrackInputArray->Reset();
117 count = 0;
118 // stop once we have enough tracks
119 while((track = static_cast<Candidate *>(fItTrackInputArray->Next())) and count < fNtracks)
120 {
121 const TLorentzVector &trkMomentum = track->Momentum;
122 tpt = trkMomentum.Pt();
123 if(tpt < fPtMin) continue;
124
125 d0 = TMath::Abs(track->D0);
126 if(d0 > fIPmax) continue;
127
128 dr = jetMomentum.DeltaR(trkMomentum);
129 if(dr > fDeltaR) continue;
130
131 xd = track->Xd;
132 yd = track->Yd;
133 zd = track->Zd;
134 dd0 = TMath::Abs(track->ErrorD0);
135 dz = TMath::Abs(track->DZ);
136 ddz = TMath::Abs(track->ErrorDZ);
137
138 if(fUse3D)
139 {
140 sign = (jpx * xd + jpy * yd + jpz * zd > 0.0) ? 1 : -1;
141 //add transverse and longitudinal significances in quadrature
142 sip = sign * TMath::Sqrt(TMath::Power(d0 / dd0, 2) + TMath::Power(dz / ddz, 2));
143 }
144 else
145 {
146 sign = (jpx * xd + jpy * yd > 0.0) ? 1 : -1;
147 sip = sign * d0 / TMath::Abs(dd0);
148 }
149
150 if(sip > fSigMin) count++;
151 }
152
153 // set BTag flag to true if count >= Ntracks
154 jet->BTag |= (count >= fNtracks) << fBitNumber;
155 }
156}
157
158//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.