Fork me on GitHub

source: git/modules/TrackCountingBTagging.cc@ a190d94

ImprovedOutputFile Timing dual_readout llp
Last change on this file since a190d94 was 1fa50c2, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix GPLv3 header

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