Fork me on GitHub

source: git/modules/TrackCountingBTagging.cc@ 5a1f24a

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5a1f24a was 9458496b, checked in by Kevin Pedro <kpedro88@…>, 7 years ago

optimize track counting algorithm

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