Fork me on GitHub

source: git/modules/VertexFinder.cc@ 2600216

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

added variables for vertexing

  • Property mode set to 100644
File size: 15.3 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 // Grow the cluster if GrowSeeds is true
124 if (fGrowSeeds)
125 growCluster (cluster->first);
126
127 // If the cluster still has fewer than MinNDF tracks, release the tracks;
128 // otherwise, mark the seed track as claimed
129 if ((int) clusterIDToInt.at (cluster->first).at ("ndf") < fMinNDF)
130 {
131 for (map<unsigned, map<string, int> >::iterator track = trackIDToInt.begin (); track != trackIDToInt.end (); track++)
132 {
133 if (track->second.at ("clusterIndex") != (int) cluster->first)
134 continue;
135 track->second["clusterIndex"] = -1;
136 trackIDToBool[track->first]["claimed"] = false;
137 }
138 }
139 else
140 trackIDToBool[clusterIDToInt.at (cluster->first).at ("seed")]["claimed"] = true;
141 }
142////////////////////////////////////////////////////////////////////////////////
143
144////////////////////////////////////////////////////////////////////////////////
145// Add tracks to the output array after updating their ClusterIndex.
146////////////////////////////////////////////////////////////////////////////////
147 fItInputArray->Reset ();
148 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
149 {
150 if (candidate->Momentum.Pt () < fMinPT || fabs (candidate->Momentum.Eta ()) > fMaxEta)
151 continue;
152 candidate->ClusterIndex = trackIDToInt.at (candidate->GetUniqueID ()).at ("clusterIndex");
153 fOutputArray->Add(candidate);
154 }
155////////////////////////////////////////////////////////////////////////////////
156
157////////////////////////////////////////////////////////////////////////////////
158// Add clusters with at least MinNDF tracks to the output array in order of
159// descending sum(pt**2).
160////////////////////////////////////////////////////////////////////////////////
161 clusterSumPT2.clear ();
162 for (map<unsigned, map<string, int> >::const_iterator cluster = clusterIDToInt.begin (); cluster != clusterIDToInt.end (); cluster++)
163 {
164 if (cluster->second.at ("ndf") < fMinNDF)
165 continue;
166 clusterSumPT2.push_back (make_pair (cluster->first, clusterIDToDouble.at (cluster->first).at ("sumPT2")));
167 }
168 sort (clusterSumPT2.begin (), clusterSumPT2.end (), secondDescending);
169 for (vector<pair<unsigned, double> >::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
186void
187VertexFinder::createSeeds ()
188{
189 Candidate *candidate;
190 unsigned 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<unsigned, double> >::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<unsigned, double> >::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
248void
249VertexFinder::growCluster (const unsigned clusterIndex)
250{
251 bool done = false;
252 unsigned nearestID;
253 int oldClusterIndex;
254 double nearestDistance;
255 vector<unsigned> 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 for (map<unsigned, map<string, double> >::const_iterator track = trackIDToDouble.begin (); track != trackIDToDouble.end (); track++)
275 {
276 if (trackIDToBool.at (track->first).at ("claimed") || trackIDToInt.at (track->first).at ("clusterIndex") == (int) clusterIndex)
277 continue;
278 double distance = fabs (clusterIDToDouble.at (clusterIndex).at ("z") - track->second.at ("z")) / hypot (clusterIDToDouble.at (clusterIndex).at ("ez"), track->second.at ("ez"));
279 if (nearestDistance < 0.0 || distance < nearestDistance)
280 {
281 nearestID = track->first;
282 nearestDistance = distance;
283 }
284 if (distance < 10.0 * fSigma)
285 nearTracks.push_back (track->first);
286 }
287 }
288 else
289 {
290 for (vector<unsigned>::const_iterator track = nearTracks.begin (); track != nearTracks.end (); track++)
291 {
292 if (trackIDToBool.at (*track).at ("claimed") || trackIDToInt.at (*track).at ("clusterIndex") == (int) clusterIndex)
293 continue;
294 double distance = fabs (clusterIDToDouble.at (clusterIndex).at ("z") - trackIDToDouble.at (*track).at ("z")) / hypot (clusterIDToDouble.at (clusterIndex).at ("ez"), trackIDToDouble.at (*track).at ("ez"));
295 if (nearestDistance < 0.0 || distance < nearestDistance)
296 {
297 nearestID = *track;
298 nearestDistance = distance;
299 }
300 }
301 }
302
303 // If no tracks within Sigma of the cluster were found, stop growing.
304 done = nearestDistance > fSigma || nearestDistance < 0.0;
305 if (done)
306 {
307 continue;
308 }
309
310 // Add the nearest track within Sigma to the cluster. If it already
311 // belonged to another cluster, remove it from that cluster first.
312 if (nearestDistance < fSigma)
313 {
314 oldClusterIndex = trackIDToInt.at (nearestID).at ("clusterIndex");
315 if (oldClusterIndex >= 0)
316 removeTrackFromCluster (nearestID, oldClusterIndex);
317
318 trackIDToBool[nearestID]["claimed"] = true;
319 addTrackToCluster (nearestID, clusterIndex);
320 }
321 }
322////////////////////////////////////////////////////////////////////////////////
323}
324
325double
326VertexFinder::weight (const unsigned trackID)
327{
328 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"))));
329}
330
331void
332VertexFinder::removeTrackFromCluster (const unsigned trackID, const unsigned clusterID)
333{
334 double wz = weight (trackID);
335
336 trackIDToInt[trackID]["clusterIndex"] = -1;
337 clusterIDToInt[clusterID]["ndf"]--;
338
339 clusterIDToDouble[clusterID]["sumZ"] -= wz * trackIDToDouble.at (trackID).at ("z");
340 clusterIDToDouble[clusterID]["errorSumZ"] -= wz * trackIDToDouble.at (trackID).at ("ez") * trackIDToDouble.at (trackID).at ("ez");
341 clusterIDToDouble[clusterID]["sumOfWeightsZ"] -= wz;
342 clusterIDToDouble[clusterID]["z"] = clusterIDToDouble.at (clusterID).at ("sumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ");
343 clusterIDToDouble[clusterID]["ez"] = sqrt ((1.0 / clusterIDToInt.at (clusterID).at ("ndf")) * (clusterIDToDouble.at (clusterID).at ("errorSumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ")));
344 clusterIDToDouble[clusterID]["sumPT2"] -= trackIDToDouble.at (trackID).at ("pt") * trackIDToDouble.at (trackID).at ("pt");
345}
346
347void
348VertexFinder::addTrackToCluster (const unsigned trackID, const unsigned clusterID)
349{
350 double wz = weight (trackID);
351
352 if (!clusterIDToInt.count (clusterID))
353 {
354 clusterIDToInt[clusterID]["ndf"] = 0;
355 clusterIDToInt[clusterID]["seed"] = trackID;
356 clusterIDToDouble[clusterID]["sumZ"] = 0.0;
357 clusterIDToDouble[clusterID]["errorSumZ"] = 0.0;
358 clusterIDToDouble[clusterID]["sumOfWeightsZ"] = 0.0;
359 clusterIDToDouble[clusterID]["sumPT2"] = 0.0;
360 }
361
362 trackIDToInt[trackID]["clusterIndex"] = clusterID;
363 clusterIDToInt[clusterID]["ndf"]++;
364
365 clusterIDToDouble[clusterID]["sumZ"] += wz * trackIDToDouble.at (trackID).at ("z");
366 clusterIDToDouble[clusterID]["errorSumZ"] += wz * trackIDToDouble.at (trackID).at ("ez") * trackIDToDouble.at (trackID).at ("ez");
367 clusterIDToDouble[clusterID]["sumOfWeightsZ"] += wz;
368 clusterIDToDouble[clusterID]["z"] = clusterIDToDouble.at (clusterID).at ("sumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ");
369 clusterIDToDouble[clusterID]["ez"] = sqrt ((1.0 / clusterIDToInt.at (clusterID).at ("ndf")) * (clusterIDToDouble.at (clusterID).at ("errorSumZ") / clusterIDToDouble.at (clusterID).at ("sumOfWeightsZ")));
370 clusterIDToDouble[clusterID]["sumPT2"] += trackIDToDouble.at (trackID).at ("pt") * trackIDToDouble.at (trackID).at ("pt");
371}
372
373//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.