1 | /** \class VertexSorter
|
---|
2 | *
|
---|
3 | *
|
---|
4 | * Sorts vertices according to different criteria
|
---|
5 | *
|
---|
6 | * \authors A. Hart, M. Selvaggi
|
---|
7 | *
|
---|
8 | *
|
---|
9 | */
|
---|
10 |
|
---|
11 | #include "modules/VertexSorter.h"
|
---|
12 |
|
---|
13 | #include "classes/DelphesClasses.h"
|
---|
14 | #include "classes/DelphesFactory.h"
|
---|
15 | #include "classes/DelphesFormula.h"
|
---|
16 | #include "classes/DelphesPileUpReader.h"
|
---|
17 |
|
---|
18 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
19 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
20 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
21 |
|
---|
22 | #include "TMath.h"
|
---|
23 | #include "TString.h"
|
---|
24 | #include "TFormula.h"
|
---|
25 | #include "TRandom3.h"
|
---|
26 | #include "TObjArray.h"
|
---|
27 | #include "TDatabasePDG.h"
|
---|
28 | #include "TLorentzVector.h"
|
---|
29 | #include "TMatrixT.h"
|
---|
30 | #include "TVector3.h"
|
---|
31 |
|
---|
32 | #include <utility>
|
---|
33 | #include <algorithm>
|
---|
34 | #include <stdexcept>
|
---|
35 | #include <iostream>
|
---|
36 | #include <vector>
|
---|
37 | #include <map>
|
---|
38 | #include <string>
|
---|
39 |
|
---|
40 | using namespace std;
|
---|
41 |
|
---|
42 | static const Double_t mm = 1.;
|
---|
43 | static const Double_t m = 1000.*mm;
|
---|
44 | static const Double_t ns = 1.;
|
---|
45 | static const Double_t s = 1.e+9 *ns;
|
---|
46 | static const Double_t c_light = 2.99792458e+8 * m/s;
|
---|
47 |
|
---|
48 | //------------------------------------------------------------------------------
|
---|
49 |
|
---|
50 | VertexSorter::VertexSorter() :
|
---|
51 | fInputArray(NULL), fTrackInputArray(NULL), fItTrackInputArray(NULL), fJetInputArray(NULL), fItJetInputArray(NULL), fOutputArray(NULL)
|
---|
52 | {
|
---|
53 | }
|
---|
54 |
|
---|
55 | //------------------------------------------------------------------------------
|
---|
56 |
|
---|
57 | VertexSorter::~VertexSorter()
|
---|
58 | {
|
---|
59 | }
|
---|
60 |
|
---|
61 | //------------------------------------------------------------------------------
|
---|
62 |
|
---|
63 | void VertexSorter::Init()
|
---|
64 | {
|
---|
65 | fInputArray = ImportArray(GetString("InputArray", "VertexFinder/vertices"));
|
---|
66 |
|
---|
67 | fTrackInputArray = ImportArray(GetString("TrackInputArray", "VertexFinder/tracks"));
|
---|
68 | fItTrackInputArray = fTrackInputArray->MakeIterator();
|
---|
69 |
|
---|
70 | if (string (GetString("JetInputArray", "")) != "")
|
---|
71 | {
|
---|
72 | fJetInputArray = ImportArray(GetString("JetInputArray", ""));
|
---|
73 | fItJetInputArray = fJetInputArray->MakeIterator();
|
---|
74 | }
|
---|
75 |
|
---|
76 | // import beamspot
|
---|
77 | try
|
---|
78 | {
|
---|
79 | fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
|
---|
80 | }
|
---|
81 | catch(runtime_error &e)
|
---|
82 | {
|
---|
83 | fBeamSpotInputArray = 0;
|
---|
84 | }
|
---|
85 |
|
---|
86 | fOutputArray = ExportArray(GetString("OutputArray", "clusters"));
|
---|
87 |
|
---|
88 | fMethod = GetString("Method", "BTV");
|
---|
89 | }
|
---|
90 |
|
---|
91 | //------------------------------------------------------------------------------
|
---|
92 |
|
---|
93 | void VertexSorter::Finish()
|
---|
94 | {
|
---|
95 | if(fItTrackInputArray) delete fItTrackInputArray;
|
---|
96 | if(fItJetInputArray) delete fItJetInputArray;
|
---|
97 | }
|
---|
98 |
|
---|
99 | //------------------------------------------------------------------------------
|
---|
100 |
|
---|
101 | static Bool_t secondDescending (pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
|
---|
102 | {
|
---|
103 | return (pair0.second > pair1.second);
|
---|
104 | }
|
---|
105 |
|
---|
106 | static Bool_t secondAscending (pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
|
---|
107 | {
|
---|
108 | return (pair0.second < pair1.second);
|
---|
109 | }
|
---|
110 |
|
---|
111 | //------------------------------------------------------------------------------
|
---|
112 |
|
---|
113 | void VertexSorter::Process()
|
---|
114 | {
|
---|
115 | Candidate *candidate, *jetCandidate, *beamSpotCandidate;
|
---|
116 | map<Int_t, UInt_t> clusterIDToIndex;
|
---|
117 | map<Int_t, Double_t> clusterIDToSumPT2;
|
---|
118 | vector<pair<Int_t, Double_t> > sortedClusterIDs;
|
---|
119 |
|
---|
120 | for (Int_t iCluster = 0; iCluster < fInputArray->GetEntries (); iCluster++)
|
---|
121 | {
|
---|
122 | const Candidate &cluster = *((Candidate *) fInputArray->At (iCluster));
|
---|
123 | clusterIDToIndex[cluster.ClusterIndex] = iCluster;
|
---|
124 | clusterIDToSumPT2[cluster.ClusterIndex] = 0.0;
|
---|
125 | }
|
---|
126 |
|
---|
127 | if(fMethod == "BTV")
|
---|
128 | {
|
---|
129 | if (!fJetInputArray)
|
---|
130 | {
|
---|
131 | cout << "BTV PV sorting selected, but no jet collection given!" << endl;
|
---|
132 | throw 0;
|
---|
133 | }
|
---|
134 |
|
---|
135 | fItTrackInputArray->Reset();
|
---|
136 | while((candidate = static_cast<Candidate*>(fItTrackInputArray->Next())))
|
---|
137 | {
|
---|
138 | if (candidate->Momentum.Pt () < 1.0)
|
---|
139 | continue;
|
---|
140 | if (candidate->ClusterIndex < 0)
|
---|
141 | continue;
|
---|
142 | TLorentzVector p (candidate->Momentum.Px (), candidate->Momentum.Py (), candidate->Momentum.Pz (), candidate->Momentum.E ());
|
---|
143 | Bool_t isInJet = false;
|
---|
144 |
|
---|
145 | fItJetInputArray->Reset();
|
---|
146 | while((jetCandidate = static_cast<Candidate*>(fItJetInputArray->Next())))
|
---|
147 | {
|
---|
148 | if (jetCandidate->Momentum.Pt () < 30.0)
|
---|
149 | continue;
|
---|
150 | TLorentzVector q (jetCandidate->Momentum.Px (), jetCandidate->Momentum.Py (), jetCandidate->Momentum.Pz (), jetCandidate->Momentum.E ());
|
---|
151 |
|
---|
152 | if (p.DeltaR (q) > 0.4)
|
---|
153 | continue;
|
---|
154 | isInJet = true;
|
---|
155 | break;
|
---|
156 | }
|
---|
157 | if (!isInJet)
|
---|
158 | continue;
|
---|
159 |
|
---|
160 | clusterIDToSumPT2.at (candidate->ClusterIndex) += candidate->Momentum.Pt () * candidate->Momentum.Pt ();
|
---|
161 | }
|
---|
162 |
|
---|
163 | for (const auto &clusterID : clusterIDToSumPT2)
|
---|
164 | sortedClusterIDs.push_back (make_pair (clusterID.first, clusterID.second));
|
---|
165 | sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondDescending);
|
---|
166 | }
|
---|
167 | else if (fMethod == "GenClosest")
|
---|
168 | {
|
---|
169 | if (!fBeamSpotInputArray)
|
---|
170 | {
|
---|
171 | cout << "GenClosest PV sorting selected, but no beamspot collection given!" << endl;
|
---|
172 | throw 0;
|
---|
173 | }
|
---|
174 | if (!fBeamSpotInputArray->GetSize ())
|
---|
175 | {
|
---|
176 | cout << "Beamspot collection is empty!" << endl;
|
---|
177 | throw 0;
|
---|
178 | }
|
---|
179 |
|
---|
180 | beamSpotCandidate = (Candidate *) fBeamSpotInputArray->At (0);
|
---|
181 | for (Int_t iCluster = 0; iCluster < fInputArray->GetEntries (); iCluster++)
|
---|
182 | {
|
---|
183 | const Candidate &cluster = *((Candidate *) fInputArray->At (iCluster));
|
---|
184 | sortedClusterIDs.push_back (make_pair (cluster.ClusterIndex, fabs (cluster.Position.Z () - beamSpotCandidate->Position.Z ())));
|
---|
185 | }
|
---|
186 | sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondAscending);
|
---|
187 | }
|
---|
188 | else if (fMethod == "GenBest")
|
---|
189 | {
|
---|
190 | fItTrackInputArray->Reset();
|
---|
191 | while((candidate = static_cast<Candidate*>(fItTrackInputArray->Next())))
|
---|
192 | {
|
---|
193 | if (candidate->IsPU)
|
---|
194 | continue;
|
---|
195 | for (const auto &clusterID : clusterIDToIndex)
|
---|
196 | {
|
---|
197 | if (candidate->ClusterIndex != clusterID.first)
|
---|
198 | continue;
|
---|
199 | clusterIDToSumPT2.at (clusterID.first) += candidate->Momentum.Pt () * candidate->Momentum.Pt ();
|
---|
200 | }
|
---|
201 | }
|
---|
202 |
|
---|
203 | for (const auto &clusterID : clusterIDToSumPT2)
|
---|
204 | sortedClusterIDs.push_back (make_pair (clusterID.first, clusterID.second));
|
---|
205 | sort (sortedClusterIDs.begin (), sortedClusterIDs.end (), secondDescending);
|
---|
206 | }
|
---|
207 | else
|
---|
208 | {
|
---|
209 | cout << "\"" << fMethod << "\" is not a valid sorting method!" << endl;
|
---|
210 | cout << "Valid methods are:" << endl;
|
---|
211 | cout << " BTV" << endl;
|
---|
212 | cout << " GenClosest" << endl;
|
---|
213 | cout << " GenBest" << endl;
|
---|
214 | throw 0;
|
---|
215 | }
|
---|
216 |
|
---|
217 | for (const auto &clusterID : sortedClusterIDs)
|
---|
218 | {
|
---|
219 | Candidate *cluster = (Candidate *) fInputArray->At (clusterIDToIndex.at (clusterID.first));
|
---|
220 | if (fMethod == "BTV")
|
---|
221 | cluster->BTVSumPT2 = clusterID.second;
|
---|
222 | else if (fMethod == "GenClosest")
|
---|
223 | cluster->GenDeltaZ = clusterID.second;
|
---|
224 | else if (fMethod == "GenBest")
|
---|
225 | cluster->GenSumPT2 = clusterID.second;
|
---|
226 | fOutputArray->Add (cluster);
|
---|
227 | }
|
---|
228 | }
|
---|
229 |
|
---|
230 | //------------------------------------------------------------------------------
|
---|