1 | /** \class VertexFinder
|
---|
2 | *
|
---|
3 | * Cluster vertices from tracks
|
---|
4 | *
|
---|
5 | * \authors A. Hart, M. Selvaggi
|
---|
6 | *
|
---|
7 | */
|
---|
8 |
|
---|
9 |
|
---|
10 | #include "modules/VertexFinder.h"
|
---|
11 | #include "classes/DelphesClasses.h"
|
---|
12 | #include "classes/DelphesFactory.h"
|
---|
13 | #include "classes/DelphesFormula.h"
|
---|
14 | #include "classes/DelphesPileUpReader.h"
|
---|
15 |
|
---|
16 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
17 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
18 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
19 |
|
---|
20 | #include "TMath.h"
|
---|
21 | #include "TString.h"
|
---|
22 | #include "TFormula.h"
|
---|
23 | #include "TRandom3.h"
|
---|
24 | #include "TObjArray.h"
|
---|
25 | #include "TDatabasePDG.h"
|
---|
26 | #include "TLorentzVector.h"
|
---|
27 | #include "TMatrixT.h"
|
---|
28 | #include "TVector3.h"
|
---|
29 |
|
---|
30 | static const Double_t mm = 1.;
|
---|
31 | static const Double_t m = 1000.*mm;
|
---|
32 | static const Double_t ns = 1.;
|
---|
33 | static const Double_t s = 1.e+9 *ns;
|
---|
34 | static const Double_t c_light = 2.99792458e+8 * m/s;
|
---|
35 |
|
---|
36 | //------------------------------------------------------------------------------
|
---|
37 |
|
---|
38 | VertexFinder::VertexFinder() :
|
---|
39 | fSigma(0), fMinPT(0), fMaxEta(0), fSeedMinPT(0), fMinNDF(0), fGrowSeeds(0)
|
---|
40 | {
|
---|
41 | }
|
---|
42 |
|
---|
43 | //------------------------------------------------------------------------------
|
---|
44 |
|
---|
45 | VertexFinder::~VertexFinder()
|
---|
46 | {
|
---|
47 | }
|
---|
48 |
|
---|
49 | //------------------------------------------------------------------------------
|
---|
50 |
|
---|
51 | void VertexFinder::Init()
|
---|
52 | {
|
---|
53 | fSigma = GetDouble("Sigma", 3.0);
|
---|
54 | fMinPT = GetDouble("MinPT", 0.1);
|
---|
55 | fMaxEta = GetDouble("MaxEta", 10.0);
|
---|
56 | fSeedMinPT = GetDouble("SeedMinPT", 5.0);
|
---|
57 | fMinNDF = GetInt("MinNDF", 4);
|
---|
58 | fGrowSeeds = GetInt("GrowSeeds", 1);
|
---|
59 |
|
---|
60 | fInputArray = ImportArray(GetString("InputArray", "TrackSmearing/tracks"));
|
---|
61 | fItInputArray = fInputArray->MakeIterator();
|
---|
62 |
|
---|
63 | fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
|
---|
64 | fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
|
---|
65 | }
|
---|
66 |
|
---|
67 | //------------------------------------------------------------------------------
|
---|
68 |
|
---|
69 | void VertexFinder::Finish()
|
---|
70 | {
|
---|
71 | if(fItInputArray) delete fItInputArray;
|
---|
72 | }
|
---|
73 |
|
---|
74 | //------------------------------------------------------------------------------
|
---|
75 | //
|
---|
76 | Bool_t VertexFinder::secondAscending (pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
|
---|
77 | {
|
---|
78 | return (pair0.second < pair1.second);
|
---|
79 | }
|
---|
80 |
|
---|
81 | Bool_t VertexFinder::secondDescending (pair<UInt_t, Double_t> pair0, pair<UInt_t, Double_t> pair1)
|
---|
82 | {
|
---|
83 | return (pair0.second > pair1.second);
|
---|
84 | }
|
---|
85 |
|
---|
86 | void VertexFinder::Process()
|
---|
87 | {
|
---|
88 | Candidate *candidate;
|
---|
89 |
|
---|
90 | ////////////////////////////////////////////////////////////////////////////////
|
---|
91 | // Clear the track and cluster maps before starting
|
---|
92 | ////////////////////////////////////////////////////////////////////////////////
|
---|
93 | trackIDToDouble.clear ();
|
---|
94 | trackIDToInt.clear ();
|
---|
95 | trackIDToBool.clear ();
|
---|
96 | clusterIDToDouble.clear ();
|
---|
97 | clusterIDToInt.clear ();
|
---|
98 | clusterIDToBool.clear ();
|
---|
99 | trackPT.clear ();
|
---|
100 | clusterSumPT2.clear ();
|
---|
101 | ////////////////////////////////////////////////////////////////////////////////
|
---|
102 |
|
---|
103 | ////////////////////////////////////////////////////////////////////////////////
|
---|
104 | // Create the initial cluster seeds
|
---|
105 | ////////////////////////////////////////////////////////////////////////////////
|
---|
106 | createSeeds ();
|
---|
107 | ////////////////////////////////////////////////////////////////////////////////
|
---|
108 |
|
---|
109 | ////////////////////////////////////////////////////////////////////////////////
|
---|
110 | // In order of descending seed pt, grow each cluster. If a cluster ends up with
|
---|
111 | // fewer than MinNDF tracks, release the tracks for other clusters to claim.
|
---|
112 | ////////////////////////////////////////////////////////////////////////////////
|
---|
113 | sort (clusterSumPT2.begin (), clusterSumPT2.end (), secondDescending);
|
---|
114 | for (vector<pair<UInt_t, Double_t> >::const_iterator cluster = clusterSumPT2.begin (); cluster != clusterSumPT2.end (); cluster++)
|
---|
115 | {
|
---|
116 | // Skip the cluster if it no longer has any tracks
|
---|
117 | if (!clusterIDToInt.at (cluster->first).at ("ndf"))
|
---|
118 | continue;
|
---|
119 |
|
---|
120 | // Grow the cluster if GrowSeeds is true
|
---|
121 | if (fGrowSeeds)
|
---|
122 | growCluster (cluster->first);
|
---|
123 |
|
---|
124 | // If the cluster still has fewer than MinNDF tracks, release the tracks;
|
---|
125 | // otherwise, mark the seed track as claimed
|
---|
126 |
|
---|
127 | if ((Int_t) clusterIDToInt.at (cluster->first).at ("ndf") < fMinNDF)
|
---|
128 | {
|
---|
129 | for (map<UInt_t, map<string, Int_t> >::iterator track = trackIDToInt.begin (); track != trackIDToInt.end (); track++)
|
---|
130 | {
|
---|
131 | if (track->second.at ("clusterIndex") != (Int_t) cluster->first)
|
---|
132 | continue;
|
---|
133 | track->second["clusterIndex"] = -1;
|
---|
134 | trackIDToBool[track->first]["claimed"] = false;
|
---|
135 | }
|
---|
136 | }
|
---|
137 | else
|
---|
138 | trackIDToBool[clusterIDToInt.at (cluster->first).at ("seed")]["claimed"] = true;
|
---|
139 | }
|
---|
140 | ////////////////////////////////////////////////////////////////////////////////
|
---|
141 |
|
---|
142 | ////////////////////////////////////////////////////////////////////////////////
|
---|
143 | // Add tracks to the output array after updating their ClusterIndex.
|
---|
144 | ////////////////////////////////////////////////////////////////////////////////
|
---|
145 | fItInputArray->Reset ();
|
---|
146 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
147 | {
|
---|
148 | if (candidate->Momentum.Pt () < fMinPT || fabs (candidate->Momentum.Eta ()) > fMaxEta)
|
---|
149 | continue;
|
---|
150 | candidate->ClusterIndex = trackIDToInt.at (candidate->GetUniqueID ()).at ("clusterIndex");
|
---|
151 | fOutputArray->Add(candidate);
|
---|
152 | }
|
---|
153 | ////////////////////////////////////////////////////////////////////////////////
|
---|
154 |
|
---|
155 | ////////////////////////////////////////////////////////////////////////////////
|
---|
156 | // Add clusters with at least MinNDF tracks to the output array in order of
|
---|
157 | // descending sum(pt**2).
|
---|
158 | ////////////////////////////////////////////////////////////////////////////////
|
---|
159 | clusterSumPT2.clear ();
|
---|
160 | for (map<UInt_t, map<string, Int_t> >::const_iterator cluster = clusterIDToInt.begin (); cluster != clusterIDToInt.end (); cluster++)
|
---|
161 | {
|
---|
162 |
|
---|
163 | if (cluster->second.at ("ndf") < fMinNDF)
|
---|
164 | continue;
|
---|
165 | clusterSumPT2.push_back (make_pair (cluster->first, clusterIDToDouble.at (cluster->first).at ("sumPT2")));
|
---|
166 | }
|
---|
167 | sort (clusterSumPT2.begin (), clusterSumPT2.end (), secondDescending);
|
---|
168 |
|
---|
169 | for (vector<pair<UInt_t, Double_t> >::const_iterator cluster = clusterSumPT2.begin (); cluster != clusterSumPT2.end (); cluster++)
|
---|
170 | {
|
---|
171 | DelphesFactory *factory = GetFactory();
|
---|
172 | candidate = factory->NewCandidate();
|
---|
173 |
|
---|
174 | candidate->ClusterIndex = cluster->first;
|
---|
175 | candidate->ClusterNDF = clusterIDToInt.at (cluster->first).at ("ndf");
|
---|
176 | candidate->ClusterSigma = fSigma;
|
---|
177 | candidate->SumPT2 = cluster->second;
|
---|
178 | candidate->Position.SetXYZT(0.0, 0.0, clusterIDToDouble.at (cluster->first).at ("z"), 0.0);
|
---|
179 | candidate->PositionError.SetXYZT(0.0, 0.0, clusterIDToDouble.at (cluster->first).at ("ez"), 0.0);
|
---|
180 |
|
---|
181 | fVertexOutputArray->Add(candidate);
|
---|
182 | }
|
---|
183 | ////////////////////////////////////////////////////////////////////////////////
|
---|
184 | }
|
---|
185 |
|
---|
186 | void
|
---|
187 | VertexFinder::createSeeds ()
|
---|
188 | {
|
---|
189 | Candidate *candidate;
|
---|
190 | UInt_t clusterIndex = 0, maxSeeds = 0;
|
---|
191 |
|
---|
192 | ////////////////////////////////////////////////////////////////////////////////
|
---|
193 | // Loop over all tracks, initializing some variables.
|
---|
194 | ////////////////////////////////////////////////////////////////////////////////
|
---|
195 | fItInputArray->Reset();
|
---|
196 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
197 | {
|
---|
198 | if (candidate->Momentum.Pt () < fMinPT || fabs (candidate->Momentum.Eta ()) > fMaxEta)
|
---|
199 | continue;
|
---|
200 |
|
---|
201 | trackIDToDouble[candidate->GetUniqueID ()]["pt"] = candidate->Momentum.Pt ();
|
---|
202 | trackIDToDouble[candidate->GetUniqueID ()]["ept"] = candidate->ErrorPT ? candidate->ErrorPT : 1.0e-15;;
|
---|
203 | trackIDToDouble[candidate->GetUniqueID ()]["eta"] = candidate->Momentum.Eta ();
|
---|
204 |
|
---|
205 | trackIDToDouble[candidate->GetUniqueID ()]["z"] = candidate->DZ;
|
---|
206 | trackIDToDouble[candidate->GetUniqueID ()]["ez"] = candidate->ErrorDZ ? candidate->ErrorDZ : 1.0e-15;
|
---|
207 |
|
---|
208 | trackIDToInt[candidate->GetUniqueID ()]["clusterIndex"] = -1;
|
---|
209 | trackIDToInt[candidate->GetUniqueID ()]["interactionIndex"] = candidate->IsPU;
|
---|
210 |
|
---|
211 | trackIDToBool[candidate->GetUniqueID ()]["claimed"] = false;
|
---|
212 |
|
---|
213 | trackPT.push_back (make_pair (candidate->GetUniqueID (), candidate->Momentum.Pt ()));
|
---|
214 | }
|
---|
215 | ////////////////////////////////////////////////////////////////////////////////
|
---|
216 |
|
---|
217 | ////////////////////////////////////////////////////////////////////////////////
|
---|
218 | // Sort tracks by pt and leave only the SeedMinPT highest pt ones in the
|
---|
219 | // trackPT vector.
|
---|
220 | ////////////////////////////////////////////////////////////////////////////////
|
---|
221 | sort (trackPT.begin (), trackPT.end (), secondDescending);
|
---|
222 | for (vector<pair<UInt_t, Double_t> >::const_iterator track = trackPT.begin (); track != trackPT.end (); track++, maxSeeds++)
|
---|
223 | {
|
---|
224 | if (track->second < fSeedMinPT)
|
---|
225 | break;
|
---|
226 | }
|
---|
227 | // If there are no tracks with pt above MinSeedPT, create just one seed from
|
---|
228 | // the highest pt track.
|
---|
229 | if (!maxSeeds)
|
---|
230 | maxSeeds++;
|
---|
231 | if (trackPT.size () > maxSeeds)
|
---|
232 | {
|
---|
233 | trackPT.erase (trackPT.begin () + maxSeeds, trackPT.end ());
|
---|
234 | }
|
---|
235 | ////////////////////////////////////////////////////////////////////////////////
|
---|
236 |
|
---|
237 | ////////////////////////////////////////////////////////////////////////////////
|
---|
238 | // Create the seeds from the SeedMinPT highest pt tracks.
|
---|
239 | ////////////////////////////////////////////////////////////////////////////////
|
---|
240 | for (vector<pair<UInt_t, Double_t> >::const_iterator track = trackPT.begin (); track != trackPT.end (); track++, clusterIndex++)
|
---|
241 | {
|
---|
242 | addTrackToCluster (track->first, clusterIndex);
|
---|
243 | clusterSumPT2.push_back (make_pair (clusterIndex, track->second * track->second));
|
---|
244 | }
|
---|
245 | ////////////////////////////////////////////////////////////////////////////////
|
---|
246 | }
|
---|
247 |
|
---|
248 | void
|
---|
249 | VertexFinder::growCluster (const UInt_t clusterIndex)
|
---|
250 | {
|
---|
251 | Bool_t done = false;
|
---|
252 | UInt_t nearestID;
|
---|
253 | Int_t oldClusterIndex;
|
---|
254 | Double_t nearestDistance;
|
---|
255 | vector<UInt_t> nearTracks;
|
---|
256 | nearTracks.clear ();
|
---|
257 |
|
---|
258 | ////////////////////////////////////////////////////////////////////////////////
|
---|
259 | // Grow the cluster until there are no more tracks within Sigma standard
|
---|
260 | // deviations of the cluster.
|
---|
261 | ////////////////////////////////////////////////////////////////////////////////
|
---|
262 | while (!done)
|
---|
263 | {
|
---|
264 | done = true;
|
---|
265 | nearestID = 0;
|
---|
266 | nearestDistance = -1.0;
|
---|
267 |
|
---|
268 | // These two loops are for finding the nearest track to the cluster. The
|
---|
269 | // first time, the ID of each track within 10*Sigma of the cluster is
|
---|
270 | // saved in the nearTracks vector; subsequently, to save time, only the
|
---|
271 | // tracks in this vector are checked.
|
---|
272 | if (!nearTracks.size ())
|
---|
273 | {
|
---|
274 |
|
---|
275 | for (map<UInt_t, map<string, Double_t> >::const_iterator track = trackIDToDouble.begin (); track != trackIDToDouble.end (); track++)
|
---|
276 | {
|
---|
277 | if (trackIDToBool.at (track->first).at ("claimed") || trackIDToInt.at (track->first).at ("clusterIndex") == (Int_t) clusterIndex)
|
---|
278 | continue;
|
---|
279 |
|
---|
280 | Double_t distance = fabs (clusterIDToDouble.at (clusterIndex).at ("z") - track->second.at ("z")) / hypot (clusterIDToDouble.at (clusterIndex).at ("ez"), track->second.at ("ez"));
|
---|
281 | if (nearestDistance < 0.0 || distance < nearestDistance)
|
---|
282 | {
|
---|
283 | nearestID = track->first;
|
---|
284 | nearestDistance = distance;
|
---|
285 | }
|
---|
286 | if (distance < 10.0 * fSigma)
|
---|
287 | nearTracks.push_back (track->first);
|
---|
288 | }
|
---|
289 | }
|
---|
290 |
|
---|
291 | else
|
---|
292 | {
|
---|
293 | for (vector<UInt_t>::const_iterator track = nearTracks.begin (); track != nearTracks.end (); track++)
|
---|
294 | {
|
---|
295 | if (trackIDToBool.at (*track).at ("claimed") || trackIDToInt.at (*track).at ("clusterIndex") == (Int_t) clusterIndex)
|
---|
296 | continue;
|
---|
297 | Double_t distance = fabs (clusterIDToDouble.at (clusterIndex).at ("z") - trackIDToDouble.at (*track).at ("z")) / hypot (clusterIDToDouble.at (clusterIndex).at ("ez"), trackIDToDouble.at (*track).at ("ez"));
|
---|
298 | if (nearestDistance < 0.0 || distance < nearestDistance)
|
---|
299 | {
|
---|
300 | nearestID = *track;
|
---|
301 | nearestDistance = distance;
|
---|
302 | }
|
---|
303 | }
|
---|
304 | }
|
---|
305 |
|
---|
306 | // If no tracks within Sigma of the cluster were found, stop growing.
|
---|
307 | done = nearestDistance > fSigma || nearestDistance < 0.0;
|
---|
308 | if (done)
|
---|
309 | {
|
---|
310 | continue;
|
---|
311 | }
|
---|
312 |
|
---|
313 | // Add the nearest track within Sigma to the cluster. If it already
|
---|
314 | // belonged to another cluster, remove it from that cluster first.
|
---|
315 | if (nearestDistance < fSigma)
|
---|
316 | {
|
---|
317 | oldClusterIndex = trackIDToInt.at (nearestID).at ("clusterIndex");
|
---|
318 | if (oldClusterIndex >= 0)
|
---|
319 | removeTrackFromCluster (nearestID, oldClusterIndex);
|
---|
320 |
|
---|
321 | trackIDToBool[nearestID]["claimed"] = true;
|
---|
322 | addTrackToCluster (nearestID, clusterIndex);
|
---|
323 | }
|
---|
324 | }
|
---|
325 | ////////////////////////////////////////////////////////////////////////////////
|
---|
326 | }
|
---|
327 |
|
---|
328 | Double_t
|
---|
329 | VertexFinder::weight (const UInt_t trackID)
|
---|
330 | {
|
---|
331 | return ((trackIDToDouble.at (trackID).at ("pt") / (trackIDToDouble.at (trackID).at ("ept") * trackIDToDouble.at (trackID).at ("ez"))) * (trackIDToDouble.at (trackID).at ("pt") / (trackIDToDouble.at (trackID).at ("ept") * trackIDToDouble.at (trackID).at ("ez"))));
|
---|
332 | }
|
---|
333 |
|
---|
334 | void
|
---|
335 | VertexFinder::removeTrackFromCluster (const UInt_t trackID, const UInt_t clusterID)
|
---|
336 | {
|
---|
337 | Double_t wz = weight (trackID);
|
---|
338 |
|
---|
339 | trackIDToInt[trackID]["clusterIndex"] = -1;
|
---|
340 | clusterIDToInt[clusterID]["ndf"]--;
|
---|
341 |
|
---|
342 | clusterIDToDouble[clusterID]["sumZ"] -= wz * trackIDToDouble.at (trackID).at ("z");
|
---|
343 | clusterIDToDouble[clusterID]["errorSumZ"] -= wz * trackIDToDouble.at (trackID).at ("ez") * trackIDToDouble.at (trackID).at ("ez");
|
---|
344 | clusterIDToDouble[clusterID]["sumOfWeightsZ"] -= wz;
|
---|
345 | clusterIDToDouble[clusterID]["z"] = clusterIDToDouble.at (clusterID).at ("sumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ");
|
---|
346 | clusterIDToDouble[clusterID]["ez"] = sqrt ((1.0 / clusterIDToInt.at (clusterID).at ("ndf")) * (clusterIDToDouble.at (clusterID).at ("errorSumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ")));
|
---|
347 | clusterIDToDouble[clusterID]["sumPT2"] -= trackIDToDouble.at (trackID).at ("pt") * trackIDToDouble.at (trackID).at ("pt");
|
---|
348 | }
|
---|
349 |
|
---|
350 | void
|
---|
351 | VertexFinder::addTrackToCluster (const UInt_t trackID, const UInt_t clusterID)
|
---|
352 | {
|
---|
353 | Double_t wz = weight (trackID);
|
---|
354 |
|
---|
355 | if (!clusterIDToInt.count (clusterID))
|
---|
356 | {
|
---|
357 | clusterIDToInt[clusterID]["ndf"] = 0;
|
---|
358 | clusterIDToInt[clusterID]["seed"] = trackID;
|
---|
359 | clusterIDToDouble[clusterID]["sumZ"] = 0.0;
|
---|
360 | clusterIDToDouble[clusterID]["errorSumZ"] = 0.0;
|
---|
361 | clusterIDToDouble[clusterID]["sumOfWeightsZ"] = 0.0;
|
---|
362 | clusterIDToDouble[clusterID]["sumPT2"] = 0.0;
|
---|
363 | }
|
---|
364 |
|
---|
365 | trackIDToInt[trackID]["clusterIndex"] = clusterID;
|
---|
366 | clusterIDToInt[clusterID]["ndf"]++;
|
---|
367 |
|
---|
368 | clusterIDToDouble[clusterID]["sumZ"] += wz * trackIDToDouble.at (trackID).at ("z");
|
---|
369 | clusterIDToDouble[clusterID]["errorSumZ"] += wz * trackIDToDouble.at (trackID).at ("ez") * trackIDToDouble.at (trackID).at ("ez");
|
---|
370 | clusterIDToDouble[clusterID]["sumOfWeightsZ"] += wz;
|
---|
371 | clusterIDToDouble[clusterID]["z"] = clusterIDToDouble.at (clusterID).at ("sumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ");
|
---|
372 | clusterIDToDouble[clusterID]["ez"] = sqrt ((1.0 / clusterIDToInt.at (clusterID).at ("ndf")) * (clusterIDToDouble.at (clusterID).at ("errorSumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ")));
|
---|
373 | clusterIDToDouble[clusterID]["sumPT2"] += trackIDToDouble.at (trackID).at ("pt") * trackIDToDouble.at (trackID).at ("pt");
|
---|
374 | }
|
---|
375 |
|
---|
376 | //------------------------------------------------------------------------------
|
---|