Fork me on GitHub

source: git/modules/TrackCountingBTagging.cc@ 723b77d

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 723b77d was 122e1e5, checked in by Pavel Demin <pavel.demin@…>, 9 years ago

replace Sqrt with Hypot

  • Property mode set to 100644
File size: 3.9 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 // import input array(s)
75
76 fTrackInputArray = ImportArray(GetString("TrackInputArray", "Calorimeter/eflowTracks"));
77 fItTrackInputArray = fTrackInputArray->MakeIterator();
78
79 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
80 fItJetInputArray = fJetInputArray->MakeIterator();
81}
82
83//------------------------------------------------------------------------------
84
85void TrackCountingBTagging::Finish()
86{
87 if(fItTrackInputArray) delete fItTrackInputArray;
88 if(fItJetInputArray) delete fItJetInputArray;
89}
90
91//------------------------------------------------------------------------------
92
93void TrackCountingBTagging::Process()
94{
95 Candidate *jet, *track;
96
97 Double_t jpx, jpy;
98 Double_t dr, tpx, tpy, tpt;
99 Double_t xd, yd, dxy, ddxy, ip, sip;
100
101 Int_t sign;
102
103 Int_t count;
104
105 // loop over all input jets
106 fItJetInputArray->Reset();
107 while((jet = static_cast<Candidate*>(fItJetInputArray->Next())))
108 {
109 const TLorentzVector &jetMomentum = jet->Momentum;
110 jpx = jetMomentum.Px();
111 jpy = jetMomentum.Py();
112
113 // loop over all input tracks
114 fItTrackInputArray->Reset();
115 count = 0;
116 while((track = static_cast<Candidate*>(fItTrackInputArray->Next())))
117 {
118 const TLorentzVector &trkMomentum = track->Momentum;
119
120 dr = jetMomentum.DeltaR(trkMomentum);
121
122 tpt = trkMomentum.Pt();
123 tpx = trkMomentum.Px();
124 tpy = trkMomentum.Py();
125
126 xd = track->Xd;
127 yd = track->Yd;
128 dxy = TMath::Hypot(xd, yd);
129 ddxy = track->SDxy;
130
131 if(tpt < fPtMin) continue;
132 if(dr > fDeltaR) continue;
133 if(dxy > fIPmax) continue;
134
135 sign = (jpx*xd + jpy*yd > 0.0) ? 1 : -1;
136
137 ip = sign*dxy;
138 sip = ip / TMath::Abs(ddxy);
139
140 if(sip > fSigMin) count++;
141 }
142
143 // set BTag flag to true if count >= Ntracks
144 jet->BTag |= (count >= fNtracks) << fBitNumber;
145 }
146}
147
148//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.