Fork me on GitHub

source: git/modules/VertexFinder.cc@ 641cb3d

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 641cb3d was 641cb3d, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

removed couts

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