Fork me on GitHub

source: git/modules/VertexFinder.cc@ b0fa9af

ImprovedOutputFile Timing dual_readout llp
Last change on this file since b0fa9af was 0e2f49b, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

added classes for Vertex Finding

  • Property mode set to 100644
File size: 15.4 KB
Line 
1/** \class VertexFinder
2 *
3 * Merges particles from pile-up sample into event
4 *
5 *
6 * $Date: 2013-02-12 15:13:59 +0100 (Tue, 12 Feb 2013) $
7 * $Revision: 907 $
8 *
9 *
10 * \author M. Selvaggi - UCL, Louvain-la-Neuve
11 *
12 */
13
14#include "modules/VertexFinder.h"
15
16//#include "CLHEP/Units/GlobalSystemOfUnits.h"
17//#include "CLHEP/Units/GlobalPhysicalConstants.h"
18
19#include "classes/DelphesClasses.h"
20#include "classes/DelphesFactory.h"
21#include "classes/DelphesFormula.h"
22#include "classes/DelphesPileUpReader.h"
23
24#include "ExRootAnalysis/ExRootResult.h"
25#include "ExRootAnalysis/ExRootFilter.h"
26#include "ExRootAnalysis/ExRootClassifier.h"
27
28#include "TMath.h"
29#include "TString.h"
30#include "TFormula.h"
31#include "TRandom3.h"
32#include "TObjArray.h"
33#include "TDatabasePDG.h"
34#include "TLorentzVector.h"
35#include "TMatrixT.h"
36#include "TVector3.h"
37
38static const double mm = 1.;
39static const double m = 1000.*mm;
40static const double ns = 1.;
41static const double s = 1.e+9 *ns;
42static const double c_light = 2.99792458e+8 * m/s;
43
44//------------------------------------------------------------------------------
45
46VertexFinder::VertexFinder() :
47 fSigma(0), fMinPT(0), fMaxEta(0), fSeedMinPT(0), fMinNDF(0), fGrowSeeds(0)
48{
49}
50
51//------------------------------------------------------------------------------
52
53VertexFinder::~VertexFinder()
54{
55}
56
57//------------------------------------------------------------------------------
58
59void VertexFinder::Init()
60{
61 fSigma = GetDouble("Sigma", 3.0);
62 fMinPT = GetDouble("MinPT", 0.1);
63 fMaxEta = GetDouble("MaxEta", 10.0);
64 fSeedMinPT = GetDouble("SeedMinPT", 5.0);
65 fMinNDF = GetInt("MinNDF", 4);
66 fGrowSeeds = GetInt("GrowSeeds", 1);
67
68 fInputArray = ImportArray(GetString("InputArray", "TrackSmearing/tracks"));
69 fItInputArray = fInputArray->MakeIterator();
70
71 fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
72 fClusterOutputArray = ExportArray(GetString("ClusterOutputArray", "clusters"));
73}
74
75//------------------------------------------------------------------------------
76
77void VertexFinder::Finish()
78{
79 if(fItInputArray) delete fItInputArray;
80}
81
82//------------------------------------------------------------------------------
83//
84bool VertexFinder::secondAscending (pair<unsigned, double> pair0, pair<unsigned, double> pair1)
85{
86 return (pair0.second < pair1.second);
87}
88
89bool VertexFinder::secondDescending (pair<unsigned, double> pair0, pair<unsigned, double> pair1)
90{
91 return (pair0.second > pair1.second);
92}
93
94void VertexFinder::Process()
95{
96 Candidate *candidate;
97
98////////////////////////////////////////////////////////////////////////////////
99// Clear the track and cluster maps before starting
100////////////////////////////////////////////////////////////////////////////////
101 trackIDToDouble.clear ();
102 trackIDToInt.clear ();
103 trackIDToBool.clear ();
104 clusterIDToDouble.clear ();
105 clusterIDToInt.clear ();
106 clusterIDToBool.clear ();
107 trackPT.clear ();
108 clusterSumPT2.clear ();
109////////////////////////////////////////////////////////////////////////////////
110
111////////////////////////////////////////////////////////////////////////////////
112// Create the initial cluster seeds
113////////////////////////////////////////////////////////////////////////////////
114 createSeeds ();
115////////////////////////////////////////////////////////////////////////////////
116
117////////////////////////////////////////////////////////////////////////////////
118// In order of descending seed pt, grow each cluster. If a cluster ends up with
119// fewer than MinNDF tracks, release the tracks for other clusters to claim.
120////////////////////////////////////////////////////////////////////////////////
121 sort (clusterSumPT2.begin (), clusterSumPT2.end (), secondDescending);
122 for (vector<pair<unsigned, double> >::const_iterator cluster = clusterSumPT2.begin (); cluster != clusterSumPT2.end (); cluster++)
123 {
124 // Skip the cluster if it no longer has any tracks
125 if (!clusterIDToInt.at (cluster->first).at ("ndf"))
126 continue;
127 // Grow the cluster if GrowSeeds is true
128 if (fGrowSeeds)
129 growCluster (cluster->first);
130
131 // If the cluster still has fewer than MinNDF tracks, release the tracks;
132 // otherwise, mark the seed track as claimed
133 if ((int) clusterIDToInt.at (cluster->first).at ("ndf") < fMinNDF)
134 {
135 for (map<unsigned, map<string, int> >::iterator track = trackIDToInt.begin (); track != trackIDToInt.end (); track++)
136 {
137 if (track->second.at ("clusterIndex") != (int) cluster->first)
138 continue;
139 track->second["clusterIndex"] = -1;
140 trackIDToBool[track->first]["claimed"] = false;
141 }
142 }
143 else
144 trackIDToBool[clusterIDToInt.at (cluster->first).at ("seed")]["claimed"] = true;
145 }
146////////////////////////////////////////////////////////////////////////////////
147
148////////////////////////////////////////////////////////////////////////////////
149// Add tracks to the output array after updating their ClusterIndex.
150////////////////////////////////////////////////////////////////////////////////
151 fItInputArray->Reset ();
152 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
153 {
154 if (candidate->Momentum.Pt () < fMinPT || fabs (candidate->Momentum.Eta ()) > fMaxEta)
155 continue;
156 candidate->ClusterIndex = trackIDToInt.at (candidate->GetUniqueID ()).at ("clusterIndex");
157 fOutputArray->Add(candidate);
158 }
159////////////////////////////////////////////////////////////////////////////////
160
161////////////////////////////////////////////////////////////////////////////////
162// Add clusters with at least MinNDF tracks to the output array in order of
163// descending sum(pt**2).
164////////////////////////////////////////////////////////////////////////////////
165 clusterSumPT2.clear ();
166 for (map<unsigned, map<string, int> >::const_iterator cluster = clusterIDToInt.begin (); cluster != clusterIDToInt.end (); cluster++)
167 {
168 if (cluster->second.at ("ndf") < fMinNDF)
169 continue;
170 clusterSumPT2.push_back (make_pair (cluster->first, clusterIDToDouble.at (cluster->first).at ("sumPT2")));
171 }
172 sort (clusterSumPT2.begin (), clusterSumPT2.end (), secondDescending);
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 fClusterOutputArray->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 for (map<unsigned, map<string, double> >::const_iterator track = trackIDToDouble.begin (); track != trackIDToDouble.end (); track++)
279 {
280 if (trackIDToBool.at (track->first).at ("claimed") || trackIDToInt.at (track->first).at ("clusterIndex") == (int) clusterIndex)
281 continue;
282 double distance = fabs (clusterIDToDouble.at (clusterIndex).at ("z") - track->second.at ("z")) / hypot (clusterIDToDouble.at (clusterIndex).at ("ez"), track->second.at ("ez"));
283 if (nearestDistance < 0.0 || distance < nearestDistance)
284 {
285 nearestID = track->first;
286 nearestDistance = distance;
287 }
288 if (distance < 10.0 * fSigma)
289 nearTracks.push_back (track->first);
290 }
291 }
292 else
293 {
294 for (vector<unsigned>::const_iterator track = nearTracks.begin (); track != nearTracks.end (); track++)
295 {
296 if (trackIDToBool.at (*track).at ("claimed") || trackIDToInt.at (*track).at ("clusterIndex") == (int) clusterIndex)
297 continue;
298 double distance = fabs (clusterIDToDouble.at (clusterIndex).at ("z") - trackIDToDouble.at (*track).at ("z")) / hypot (clusterIDToDouble.at (clusterIndex).at ("ez"), trackIDToDouble.at (*track).at ("ez"));
299 if (nearestDistance < 0.0 || distance < nearestDistance)
300 {
301 nearestID = *track;
302 nearestDistance = distance;
303 }
304 }
305 }
306
307 // If no tracks within Sigma of the cluster were found, stop growing.
308 done = nearestDistance > fSigma || nearestDistance < 0.0;
309 if (done)
310 {
311 continue;
312 }
313
314 // Add the nearest track within Sigma to the cluster. If it already
315 // belonged to another cluster, remove it from that cluster first.
316 if (nearestDistance < fSigma)
317 {
318 oldClusterIndex = trackIDToInt.at (nearestID).at ("clusterIndex");
319 if (oldClusterIndex >= 0)
320 removeTrackFromCluster (nearestID, oldClusterIndex);
321
322 trackIDToBool[nearestID]["claimed"] = true;
323 addTrackToCluster (nearestID, clusterIndex);
324 }
325 }
326////////////////////////////////////////////////////////////////////////////////
327}
328
329double
330VertexFinder::weight (const unsigned trackID)
331{
332 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"))));
333}
334
335void
336VertexFinder::removeTrackFromCluster (const unsigned trackID, const unsigned clusterID)
337{
338 double wz = weight (trackID);
339
340 trackIDToInt[trackID]["clusterIndex"] = -1;
341 clusterIDToInt[clusterID]["ndf"]--;
342
343 clusterIDToDouble[clusterID]["sumZ"] -= wz * trackIDToDouble.at (trackID).at ("z");
344 clusterIDToDouble[clusterID]["errorSumZ"] -= wz * trackIDToDouble.at (trackID).at ("ez") * trackIDToDouble.at (trackID).at ("ez");
345 clusterIDToDouble[clusterID]["sumOfWeightsZ"] -= wz;
346 clusterIDToDouble[clusterID]["z"] = clusterIDToDouble.at (clusterID).at ("sumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ");
347 clusterIDToDouble[clusterID]["ez"] = sqrt ((1.0 / clusterIDToInt.at (clusterID).at ("ndf")) * (clusterIDToDouble.at (clusterID).at ("errorSumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ")));
348 clusterIDToDouble[clusterID]["sumPT2"] -= trackIDToDouble.at (trackID).at ("pt") * trackIDToDouble.at (trackID).at ("pt");
349}
350
351void
352VertexFinder::addTrackToCluster (const unsigned trackID, const unsigned clusterID)
353{
354 double wz = weight (trackID);
355
356 if (!clusterIDToInt.count (clusterID))
357 {
358 clusterIDToInt[clusterID]["ndf"] = 0;
359 clusterIDToInt[clusterID]["seed"] = trackID;
360 clusterIDToDouble[clusterID]["sumZ"] = 0.0;
361 clusterIDToDouble[clusterID]["errorSumZ"] = 0.0;
362 clusterIDToDouble[clusterID]["sumOfWeightsZ"] = 0.0;
363 clusterIDToDouble[clusterID]["sumPT2"] = 0.0;
364 }
365
366 trackIDToInt[trackID]["clusterIndex"] = clusterID;
367 clusterIDToInt[clusterID]["ndf"]++;
368
369 clusterIDToDouble[clusterID]["sumZ"] += wz * trackIDToDouble.at (trackID).at ("z");
370 clusterIDToDouble[clusterID]["errorSumZ"] += wz * trackIDToDouble.at (trackID).at ("ez") * trackIDToDouble.at (trackID).at ("ez");
371 clusterIDToDouble[clusterID]["sumOfWeightsZ"] += wz;
372 clusterIDToDouble[clusterID]["z"] = clusterIDToDouble.at (clusterID).at ("sumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ");
373 clusterIDToDouble[clusterID]["ez"] = sqrt ((1.0 / clusterIDToInt.at (clusterID).at ("ndf")) * (clusterIDToDouble.at (clusterID).at ("errorSumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ")));
374 clusterIDToDouble[clusterID]["sumPT2"] += trackIDToDouble.at (trackID).at ("pt") * trackIDToDouble.at (trackID).at ("pt");
375}
376
377//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.