Fork me on GitHub

source: git/modules/TrackCountingBTagging.cc@ 529fe78

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 529fe78 was e4c3fef, checked in by pavel <pavel@…>, 11 years ago

several minor corrections

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