Fork me on GitHub

source: git/modules/TrackCountingBTagging.cc@ 2c8865f

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

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

  • Property mode set to 100644
File size: 4.4 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[a0431dc]19/** \class TrackCountingBTagging
20 *
[e4c3fef]21 * b-tagging algorithm based on counting tracks with large impact parameter
[a0431dc]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"
[341014c]35#include "TMath.h"
36#include "TObjArray.h"
37#include "TRandom3.h"
38#include "TString.h"
[a0431dc]39
[e4c3fef]40#include <algorithm>
[a0431dc]41#include <iostream>
42#include <sstream>
[341014c]43#include <stdexcept>
[a0431dc]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);
[e4c3fef]65
[6cdc544]66 fPtMin = GetDouble("TrackPtMin", 1.0);
67 fDeltaR = GetDouble("DeltaR", 0.3);
68 fIPmax = GetDouble("TrackIPMax", 2.0);
[e4c3fef]69
[6cdc544]70 fSigMin = GetDouble("SigMin", 6.5);
71 fNtracks = GetInt("Ntracks", 3);
[e4c3fef]72
[e40b9cf]73 fUse3D = GetBool("Use3D", false);
74
[a0431dc]75 // import input array(s)
76
77 fTrackInputArray = ImportArray(GetString("TrackInputArray", "Calorimeter/eflowTracks"));
78 fItTrackInputArray = fTrackInputArray->MakeIterator();
[e4c3fef]79
[a0431dc]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;
[e4c3fef]97
[e40b9cf]98 Double_t jpx, jpy, jpz;
[1dad056]99 Double_t dr, tpt;
[e40b9cf]100 Double_t xd, yd, zd, d0, dd0, dz, ddz, sip;
[e4c3fef]101
[a0431dc]102 Int_t sign;
[e4c3fef]103
[a0431dc]104 Int_t count;
[e4c3fef]105
[a0431dc]106 // loop over all input jets
107 fItJetInputArray->Reset();
[341014c]108 while((jet = static_cast<Candidate *>(fItJetInputArray->Next())))
[a0431dc]109 {
110 const TLorentzVector &jetMomentum = jet->Momentum;
111 jpx = jetMomentum.Px();
112 jpy = jetMomentum.Py();
[e40b9cf]113 jpz = jetMomentum.Pz();
[e4c3fef]114
[a0431dc]115 // loop over all input tracks
116 fItTrackInputArray->Reset();
117 count = 0;
[9458496b]118 // stop once we have enough tracks
[341014c]119 while((track = static_cast<Candidate *>(fItTrackInputArray->Next())) and count < fNtracks)
[a0431dc]120 {
[e4c3fef]121 const TLorentzVector &trkMomentum = track->Momentum;
[9458496b]122 tpt = trkMomentum.Pt();
123 if(tpt < fPtMin) continue;
[341014c]124
[9458496b]125 d0 = TMath::Abs(track->D0);
126 if(d0 > fIPmax) continue;
127
[e4c3fef]128 dr = jetMomentum.DeltaR(trkMomentum);
[9458496b]129 if(dr > fDeltaR) continue;
130
[6cdc544]131 xd = track->Xd;
132 yd = track->Yd;
[e40b9cf]133 zd = track->Zd;
134 dd0 = TMath::Abs(track->ErrorD0);
135 dz = TMath::Abs(track->DZ);
136 ddz = TMath::Abs(track->ErrorDZ);
[e4c3fef]137
[341014c]138 if(fUse3D)
139 {
140 sign = (jpx * xd + jpy * yd + jpz * zd > 0.0) ? 1 : -1;
[9458496b]141 //add transverse and longitudinal significances in quadrature
[341014c]142 sip = sign * TMath::Sqrt(TMath::Power(d0 / dd0, 2) + TMath::Power(dz / ddz, 2));
[e40b9cf]143 }
[341014c]144 else
145 {
146 sign = (jpx * xd + jpy * yd > 0.0) ? 1 : -1;
[e40b9cf]147 sip = sign * d0 / TMath::Abs(dd0);
148 }
[e4c3fef]149
150 if(sip > fSigMin) count++;
[a0431dc]151 }
[e4c3fef]152
[a0431dc]153 // set BTag flag to true if count >= Ntracks
[e4c3fef]154 jet->BTag |= (count >= fNtracks) << fBitNumber;
[a0431dc]155 }
156}
157
158//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.