1 | /** \class VertexFinderDA4D
|
---|
2 | *
|
---|
3 | * Cluster vertices from tracks using deterministic annealing and timing information
|
---|
4 | *
|
---|
5 | * \authors O. Cerri
|
---|
6 | *
|
---|
7 | */
|
---|
8 |
|
---|
9 |
|
---|
10 | #include "modules/VertexFinderDA4D.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 "TLatex.h"
|
---|
29 | #include "TVector3.h"
|
---|
30 |
|
---|
31 | #include "TAxis.h"
|
---|
32 | #include "TGraphErrors.h"
|
---|
33 | #include "TCanvas.h"
|
---|
34 | #include "TString.h"
|
---|
35 | #include "TLegend.h"
|
---|
36 | #include "TFile.h"
|
---|
37 | #include "TColor.h"
|
---|
38 | #include "TLegend.h"
|
---|
39 |
|
---|
40 | #include <utility>
|
---|
41 | #include <algorithm>
|
---|
42 | #include <stdexcept>
|
---|
43 | #include <iostream>
|
---|
44 | #include <vector>
|
---|
45 |
|
---|
46 | using namespace std;
|
---|
47 |
|
---|
48 | namespace vtx_DAZT
|
---|
49 | {
|
---|
50 | static const Double_t c_light = 2.99792458e+8; // [m/s]
|
---|
51 | }
|
---|
52 | using namespace vtx_DAZT;
|
---|
53 |
|
---|
54 | //------------------------------------------------------------------------------
|
---|
55 |
|
---|
56 | VertexFinderDA4D::VertexFinderDA4D()
|
---|
57 | {
|
---|
58 | fVerbose = 0;
|
---|
59 | fMaxIterations = 0;
|
---|
60 | fBetaMax = 0;
|
---|
61 | fBetaStop = 0;
|
---|
62 | fBetaPurge = 0;
|
---|
63 | fVertexZSize = 0;
|
---|
64 | fVertexTSize = 0;
|
---|
65 | fCoolingFactor = 0;
|
---|
66 | fDzCutOff = 0;
|
---|
67 | fD0CutOff = 0;
|
---|
68 | fDtCutOff = 0;
|
---|
69 | fPtMin = 0;
|
---|
70 | fPtMax = 0;
|
---|
71 | fD2Merge = 0;
|
---|
72 | fMuOutlayer = 0;
|
---|
73 | fMinTrackProb = 0;
|
---|
74 | }
|
---|
75 |
|
---|
76 | //------------------------------------------------------------------------------
|
---|
77 |
|
---|
78 | VertexFinderDA4D::~VertexFinderDA4D()
|
---|
79 | {
|
---|
80 | }
|
---|
81 |
|
---|
82 | //------------------------------------------------------------------------------
|
---|
83 |
|
---|
84 | void VertexFinderDA4D::Init()
|
---|
85 | {
|
---|
86 | fVerbose = GetInt("Verbose", 0);
|
---|
87 |
|
---|
88 | fMaxIterations = GetInt("MaxIterations", 100);
|
---|
89 | fMaxVertexNumber = GetInt("MaxVertexNumber", 500);
|
---|
90 |
|
---|
91 | fBetaMax = GetDouble("BetaMax", 1.5);
|
---|
92 | fBetaPurge = GetDouble("BetaPurge", 1.);
|
---|
93 | fBetaStop = GetDouble("BetaStop", 0.2);
|
---|
94 |
|
---|
95 | fVertexZSize = GetDouble("VertexZSize", 0.1); //in mm
|
---|
96 | fVertexTSize = 1E12*GetDouble("VertexTimeSize", 15E-12); //Convert from [s] to [ps]
|
---|
97 |
|
---|
98 | fCoolingFactor = GetDouble("CoolingFactor", 0.8); // Multiply T so to cooldown must be <1
|
---|
99 |
|
---|
100 | fDzCutOff = GetDouble("DzCutOff", 40); // For the moment 3*DzCutOff is hard cut off for the considered tracks
|
---|
101 | fD0CutOff = GetDouble("D0CutOff", .5); // d0/sigma_d0, used to compute the pi (weight) of the track
|
---|
102 | fDtCutOff = GetDouble("DtCutOff", 160); // [ps], 3*DtCutOff is hard cut off for tracks
|
---|
103 | fPtMin = GetDouble("PtMin", 0.5); // Minimum pt accepted for tracks
|
---|
104 | fPtMax = GetDouble("PtMax", 50); // Maximum pt accepted for tracks
|
---|
105 |
|
---|
106 |
|
---|
107 | fD2UpdateLim = GetDouble("D2UpdateLim", .5); // ((dz/ZSize)^2+(dt/TSize)^2)/nv limit for merging vertices
|
---|
108 | fD2Merge = GetDouble("D2Merge", 4.0); // (dz/ZSize)^2+(dt/TSize)^2 limit for merging vertices
|
---|
109 | fMuOutlayer = GetDouble("MuOutlayer", 4); // Outlayer rejection exponent
|
---|
110 | fMinTrackProb = GetDouble("MinTrackProb", 0.6); // Minimum probability to be assigned at a vertex
|
---|
111 | fMinNTrack = GetInt("MinNTrack", 10); // Minimum number of tracks per vertex
|
---|
112 |
|
---|
113 | fFigFolderPath = GetString("DebugFigPath", ".");
|
---|
114 |
|
---|
115 | fInputArray = ImportArray(GetString("InputArray", "TrackSmearing/tracks"));
|
---|
116 | fItInputArray = fInputArray->MakeIterator();
|
---|
117 |
|
---|
118 | fTrackOutputArray = ExportArray(GetString("TrackOutputArray", "tracks"));
|
---|
119 | fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
|
---|
120 |
|
---|
121 | fInputGenVtx = ImportArray(GetString("InputGenVtx", "PileUpMerger/vertices"));
|
---|
122 | fItInputGenVtx = fInputGenVtx->MakeIterator();
|
---|
123 |
|
---|
124 | if (fBetaMax < fBetaPurge)
|
---|
125 | {
|
---|
126 | fBetaPurge = fBetaMax;
|
---|
127 | if (fVerbose)
|
---|
128 | {
|
---|
129 | cout << "BetaPurge set to " << fBetaPurge << endl;
|
---|
130 | }
|
---|
131 | }
|
---|
132 |
|
---|
133 | if (fBetaPurge < fBetaStop)
|
---|
134 | {
|
---|
135 | fBetaStop = fBetaPurge;
|
---|
136 | if (fVerbose)
|
---|
137 | {
|
---|
138 | cout << "BetaPurge set to " << fBetaPurge << endl;
|
---|
139 | }
|
---|
140 | }
|
---|
141 | }
|
---|
142 |
|
---|
143 | //------------------------------------------------------------------------------
|
---|
144 |
|
---|
145 | void VertexFinderDA4D::Finish()
|
---|
146 | {
|
---|
147 | if(fItInputArray) delete fItInputArray;
|
---|
148 | }
|
---|
149 |
|
---|
150 | //------------------------------------------------------------------------------
|
---|
151 |
|
---|
152 | void VertexFinderDA4D::Process()
|
---|
153 | {
|
---|
154 | fInputArray->Sort();
|
---|
155 |
|
---|
156 | if (fVerbose)
|
---|
157 | {
|
---|
158 | cout<< endl << " Start processing vertices with VertexFinderDA4D" << endl;
|
---|
159 | cout<<" Found "<<fInputArray->GetEntriesFast()<<" input tracks"<<endl;
|
---|
160 | }
|
---|
161 |
|
---|
162 | // clusterize tracks
|
---|
163 | TObjArray *ClusterArray = new TObjArray;
|
---|
164 | clusterize(*ClusterArray);
|
---|
165 |
|
---|
166 | if(fVerbose>10)
|
---|
167 | {
|
---|
168 | unsigned int N = fEnergy_rec.size();
|
---|
169 | TGraph* gr1 = new TGraph(N, &fBeta_rec[0], &fNvtx_rec[0]);
|
---|
170 | gr1->SetName("gr1");
|
---|
171 | gr1->GetXaxis()->SetTitle("beta");
|
---|
172 | gr1->GetYaxis()->SetTitle("# Vtx");
|
---|
173 | TGraph* gr2 = new TGraph(N, &fBeta_rec[0], &fEnergy_rec[0]);
|
---|
174 | gr2->SetName("gr2");
|
---|
175 | gr2->GetXaxis()->SetTitle("beta");
|
---|
176 | gr2->GetYaxis()->SetTitle("Total Energy");
|
---|
177 | TGraph* gr3 = new TGraph(N, &fNvtx_rec[0], &fEnergy_rec[0]);
|
---|
178 | gr3->SetName("gr3");
|
---|
179 | gr3->GetXaxis()->SetTitle("# Vtx");
|
---|
180 | gr3->GetYaxis()->SetTitle("Total Energy");
|
---|
181 |
|
---|
182 | auto f = new TFile("~/Desktop/debug/EnergyStat.root", "recreate");
|
---|
183 | gr1->Write("gr1");
|
---|
184 | gr2->Write("gr2");
|
---|
185 | gr3->Write("gr3");
|
---|
186 |
|
---|
187 | f->Close();
|
---|
188 | }
|
---|
189 |
|
---|
190 | if (fVerbose){std::cout << " clustering returned "<< ClusterArray->GetEntriesFast() << " clusters from " << fInputArray->GetEntriesFast() << " input tracks" <<std::endl;}
|
---|
191 |
|
---|
192 | // //loop over vertex candidates
|
---|
193 | TIterator * ItClusterArray = ClusterArray->MakeIterator();
|
---|
194 | ItClusterArray->Reset();
|
---|
195 | Candidate *candidate;
|
---|
196 | unsigned int k = 0;
|
---|
197 | while((candidate = static_cast<Candidate*>(ItClusterArray->Next())))
|
---|
198 | {
|
---|
199 | if(fVerbose)
|
---|
200 | {
|
---|
201 | cout << Form("Cluster %d has %d tracks ", k, candidate->GetCandidates()->GetEntriesFast()) << endl;
|
---|
202 | }
|
---|
203 | if(candidate->ClusterNDF>0)
|
---|
204 | {
|
---|
205 | // Estimate the vertex resolution
|
---|
206 | // loop over tracks belonging to this vertex
|
---|
207 | TIter it1(candidate->GetCandidates());
|
---|
208 | it1.Reset();
|
---|
209 |
|
---|
210 | Candidate *track;
|
---|
211 | double sum_Dt_2 = 0;
|
---|
212 | double sum_Dz_2 = 0;
|
---|
213 | double sum_wt = 0;
|
---|
214 | double sum_wz = 0;
|
---|
215 | while((track = static_cast<Candidate*>(it1.Next())))
|
---|
216 | {
|
---|
217 | double dz = candidate->Position.Z() - track->Zd;
|
---|
218 | double dt = candidate->Position.T() - track->Td;
|
---|
219 |
|
---|
220 | double wz = track->VertexingWeight/(track->ErrorDZ*track->ErrorDZ);
|
---|
221 | double wt = track->VertexingWeight/(track->ErrorT*track->ErrorT);
|
---|
222 |
|
---|
223 | sum_Dt_2 += wt*dt*dt;
|
---|
224 | sum_Dz_2 += wz*dz*dz;
|
---|
225 | sum_wt += wt;
|
---|
226 | sum_wz += wz;
|
---|
227 | }
|
---|
228 |
|
---|
229 | double sigma_z = sqrt(sum_Dz_2/sum_wz);
|
---|
230 | double sigma_t = sqrt(sum_Dt_2/sum_wt);
|
---|
231 | candidate->PositionError.SetXYZT(0.0, 0.0, sigma_z , sigma_t);
|
---|
232 | if(fVerbose > 3)
|
---|
233 | {
|
---|
234 | cout << "k: " << k << endl;
|
---|
235 | cout << "Sigma z: " << sigma_z*1E3 << " um" << endl;
|
---|
236 | cout << "Sigma t: " << sigma_t*1E9/c_light << " ps" << endl;
|
---|
237 | }
|
---|
238 |
|
---|
239 | fVertexOutputArray->Add(candidate);
|
---|
240 | k++;
|
---|
241 | }
|
---|
242 | }// end of cluster loop
|
---|
243 |
|
---|
244 | delete ClusterArray;
|
---|
245 | }
|
---|
246 |
|
---|
247 | //------------------------------------------------------------------------------
|
---|
248 |
|
---|
249 | void VertexFinderDA4D::clusterize(TObjArray &clusters)
|
---|
250 | {
|
---|
251 | tracks_t tks;
|
---|
252 | fill(tks);
|
---|
253 | unsigned int nt=tks.getSize();
|
---|
254 | if(fVerbose)
|
---|
255 | {
|
---|
256 | cout << "Tracks added: " << nt << endl;
|
---|
257 | }
|
---|
258 | if (nt == 0) return;
|
---|
259 |
|
---|
260 |
|
---|
261 |
|
---|
262 | vertex_t vtx; // the vertex prototypes
|
---|
263 | vtx.ZSize = fVertexZSize;
|
---|
264 | vtx.TSize = fVertexTSize;
|
---|
265 | // initialize:single vertex at infinite temperature
|
---|
266 | vtx.addItem(0, 0, 1);
|
---|
267 |
|
---|
268 | // Fit the vertex at T=inf and return the starting temperature
|
---|
269 | double beta=beta0(tks, vtx);
|
---|
270 |
|
---|
271 | if( fVerbose > 1 )
|
---|
272 | {
|
---|
273 | cout << "Cluster position at T=inf: z = " << vtx.z[0] << " mm , t = " << vtx.t[0] << " ps" << " pk = " << vtx.pk[0] << endl;
|
---|
274 | cout << Form("Beta Start = %2.1e", beta) << endl;
|
---|
275 | }
|
---|
276 |
|
---|
277 | if( fVerbose > 2){cout << "Cool down untill reaching the temperature to finish increasing the number of vertexes" << endl;}
|
---|
278 |
|
---|
279 | double rho0=0.0; // start with no outlier rejection
|
---|
280 |
|
---|
281 | unsigned int last_round = 0;
|
---|
282 | while(last_round < 2)
|
---|
283 | {
|
---|
284 |
|
---|
285 | unsigned int niter=0;
|
---|
286 | double delta2 = 0;
|
---|
287 | do {
|
---|
288 | delta2 = update(beta, tks, vtx, rho0);
|
---|
289 |
|
---|
290 | if (fVerbose > 3)
|
---|
291 | {
|
---|
292 | cout << "Update " << niter << " : " << delta2 << endl;
|
---|
293 | }
|
---|
294 | niter++;
|
---|
295 | }
|
---|
296 | while (delta2 > fD2UpdateLim && niter < fMaxIterations);
|
---|
297 |
|
---|
298 |
|
---|
299 | unsigned int n_it = 0;
|
---|
300 | while(merge(vtx, fD2Merge) && n_it < fMaxIterations)
|
---|
301 | {
|
---|
302 | unsigned int niter=0;
|
---|
303 | double delta2 = 0;
|
---|
304 | do {
|
---|
305 | delta2 = update(beta, tks, vtx, rho0);
|
---|
306 | niter++;
|
---|
307 | }
|
---|
308 | while (delta2 > fD2UpdateLim && niter < fMaxIterations);
|
---|
309 | n_it++;
|
---|
310 |
|
---|
311 | }
|
---|
312 |
|
---|
313 | beta /= fCoolingFactor;
|
---|
314 |
|
---|
315 | if( beta < fBetaStop )
|
---|
316 | {
|
---|
317 | split(beta, vtx, tks);
|
---|
318 | }
|
---|
319 | else
|
---|
320 | {
|
---|
321 | beta = fBetaStop;
|
---|
322 | last_round++;
|
---|
323 | }
|
---|
324 |
|
---|
325 | if(fVerbose > 3)
|
---|
326 | {
|
---|
327 | cout << endl << endl << " ----- Beta = " << beta << " --------" << endl;
|
---|
328 | cout << "Nv: " << vtx.getSize() << endl;
|
---|
329 | }
|
---|
330 | }
|
---|
331 |
|
---|
332 | if( fVerbose > 4)
|
---|
333 | {
|
---|
334 | for(unsigned int k = 0; k < vtx.getSize(); k++)
|
---|
335 | {
|
---|
336 | cout << Form("Vertex %d next beta_c = %.3f", k, vtx.beta_c[k]) << endl;
|
---|
337 | }
|
---|
338 | }
|
---|
339 |
|
---|
340 | if(fVerbose > 2) {cout << "Adiabatic switch on of outlayr rejection" << endl;}
|
---|
341 | rho0 = 1./nt;
|
---|
342 | const double N_cycles = 10;
|
---|
343 | for(unsigned int f = 1; f <= N_cycles; f++)
|
---|
344 | {
|
---|
345 | unsigned int niter=0;
|
---|
346 | double delta2 = 0;
|
---|
347 | do {
|
---|
348 | delta2 = update(beta, tks, vtx, rho0 * f/N_cycles);
|
---|
349 | niter++;
|
---|
350 | }
|
---|
351 | while (delta2 > 0.3*fD2UpdateLim && niter < fMaxIterations);
|
---|
352 | }
|
---|
353 |
|
---|
354 | do {
|
---|
355 | beta /= fCoolingFactor;
|
---|
356 | if(beta > fBetaPurge) beta = fBetaPurge;
|
---|
357 | unsigned int i_pu = 0;
|
---|
358 | for(int min_trk = 2; min_trk<=fMinNTrack; min_trk++)
|
---|
359 | {
|
---|
360 | while( purge(vtx, tks, rho0, beta, fMinTrackProb, min_trk) )
|
---|
361 | {
|
---|
362 | unsigned int niter=0;
|
---|
363 | double delta2 = 0;
|
---|
364 | do {
|
---|
365 | delta2 = update(beta, tks, vtx, rho0);
|
---|
366 | niter++;
|
---|
367 | }
|
---|
368 | while (delta2 > fD2UpdateLim && niter < fMaxIterations);
|
---|
369 | i_pu++;
|
---|
370 | }
|
---|
371 | }
|
---|
372 |
|
---|
373 | unsigned int n_it = 0;
|
---|
374 | while(merge(vtx, fD2Merge) && n_it < fMaxIterations)
|
---|
375 | {
|
---|
376 | unsigned int niter=0;
|
---|
377 | double delta2 = 0;
|
---|
378 | do {
|
---|
379 | delta2 = update(beta, tks, vtx, rho0);
|
---|
380 | niter++;
|
---|
381 | }
|
---|
382 | while (delta2 > fD2UpdateLim && niter < fMaxIterations);
|
---|
383 | n_it++;
|
---|
384 |
|
---|
385 | }
|
---|
386 | } while( beta < fBetaPurge );
|
---|
387 |
|
---|
388 |
|
---|
389 | if(fVerbose > 2){cout << "Cooldown untill the limit before assigning track to vertices" << endl;}
|
---|
390 | last_round = 0;
|
---|
391 | while(last_round < 2)
|
---|
392 | {
|
---|
393 | unsigned int niter=0;
|
---|
394 | double delta2 = 0;
|
---|
395 | do {
|
---|
396 | delta2 = update(beta, tks, vtx, rho0);
|
---|
397 | niter++;
|
---|
398 | }
|
---|
399 | while (delta2 > 0.3*fD2UpdateLim && niter < fMaxIterations);
|
---|
400 |
|
---|
401 | beta /= fCoolingFactor;
|
---|
402 | if ( beta >= fBetaMax )
|
---|
403 | {
|
---|
404 | beta = fBetaMax;
|
---|
405 | last_round++;
|
---|
406 | }
|
---|
407 | }
|
---|
408 |
|
---|
409 |
|
---|
410 | // Build the cluster candidates
|
---|
411 | for(unsigned int k = 0; k < vtx.getSize(); k++)
|
---|
412 | {
|
---|
413 | DelphesFactory *factory = GetFactory();
|
---|
414 | Candidate * candidate = factory->NewCandidate();
|
---|
415 |
|
---|
416 | candidate->ClusterIndex = k;
|
---|
417 | candidate->Position.SetXYZT(0.0, 0.0, vtx.z[k] , vtx.t[k]*1E-9*c_light);
|
---|
418 | candidate->InitialPosition.SetXYZT(0.0, 0.0, vtx.z[k] , vtx.t[k]*1E-9*c_light);
|
---|
419 | candidate->PositionError.SetXYZT(0.0, 0.0, fVertexZSize , fVertexTSize*1E-9*c_light);
|
---|
420 | candidate->SumPT2 = 0;
|
---|
421 | candidate->SumPt = 0;
|
---|
422 | candidate->ClusterNDF = 0;
|
---|
423 |
|
---|
424 | clusters.Add(candidate);
|
---|
425 | }
|
---|
426 |
|
---|
427 |
|
---|
428 | // Assign each track to the most probable vertex
|
---|
429 | double Z_init = rho0 * exp(-beta * fMuOutlayer * fMuOutlayer); // Add fDtCutOff here toghether with this
|
---|
430 | vector<double> pk_exp_mBetaE = Compute_pk_exp_mBetaE(beta, vtx, tks, Z_init);
|
---|
431 | for(unsigned int i = 0; i< tks.getSize(); i++)
|
---|
432 | {
|
---|
433 | if(tks.w[i] <= 0) continue;
|
---|
434 |
|
---|
435 | double p_max = 0;
|
---|
436 | unsigned int k_max = 0;
|
---|
437 |
|
---|
438 | for(unsigned int k = 0; k < vtx.getSize(); k++)
|
---|
439 | {
|
---|
440 | unsigned int idx = k*nt + i;
|
---|
441 | if(pk_exp_mBetaE[idx] == 0 || tks.Z[i] == 0 || vtx.pk[k] == 0)
|
---|
442 | {
|
---|
443 | continue;
|
---|
444 | }
|
---|
445 |
|
---|
446 | double pv_max = vtx.pk[k] / (vtx.pk[k] + rho0 * exp(-beta * fMuOutlayer* fMuOutlayer));
|
---|
447 | double p = pk_exp_mBetaE[idx] / tks.Z[i];
|
---|
448 |
|
---|
449 | p /= pv_max;
|
---|
450 |
|
---|
451 | if(p > p_max)
|
---|
452 | {
|
---|
453 | p_max = p;
|
---|
454 | k_max = k;
|
---|
455 | }
|
---|
456 | }
|
---|
457 |
|
---|
458 | if(p_max > fMinTrackProb)
|
---|
459 | {
|
---|
460 | tks.tt[i]->ClusterIndex = k_max;
|
---|
461 | tks.tt[i]->InitialPosition.SetT(1E-9*vtx.t[k_max]*c_light);
|
---|
462 | tks.tt[i]->InitialPosition.SetZ(vtx.z[k_max]);
|
---|
463 |
|
---|
464 | ((Candidate *) clusters.At(k_max))->AddCandidate(tks.tt[i]);
|
---|
465 | ((Candidate *) clusters.At(k_max))->SumPT2 += tks.tt[i]->Momentum.Pt()*tks.tt[i]->Momentum.Pt();
|
---|
466 | ((Candidate *) clusters.At(k_max))->SumPt += tks.tt[i]->Momentum.Pt();
|
---|
467 | ((Candidate *) clusters.At(k_max))->ClusterNDF += 1;
|
---|
468 | }
|
---|
469 | else
|
---|
470 | {
|
---|
471 | tks.tt[i]->ClusterIndex = -1;
|
---|
472 | tks.tt[i]->InitialPosition.SetT(1E3*1000000*c_light);
|
---|
473 | tks.tt[i]->InitialPosition.SetZ(1E8);
|
---|
474 | }
|
---|
475 | fTrackOutputArray->Add(tks.tt[i]);
|
---|
476 | }
|
---|
477 |
|
---|
478 | }
|
---|
479 |
|
---|
480 | //------------------------------------------------------------------------------
|
---|
481 | // Definition of the distance metrci between track and vertex
|
---|
482 | double VertexFinderDA4D::Energy(double t_z, double v_z, double dz2_o, double t_t, double v_t, double dt2_o)
|
---|
483 | {
|
---|
484 | return (t_z - v_z)*(t_z - v_z)* dz2_o + (t_t - v_t)*(t_t - v_t)*dt2_o;
|
---|
485 | }
|
---|
486 |
|
---|
487 | //------------------------------------------------------------------------------
|
---|
488 | // Fill tks with the input candidates array
|
---|
489 | void VertexFinderDA4D::fill(tracks_t &tks)
|
---|
490 | {
|
---|
491 | tks.sum_w_o_dt2 = 0;
|
---|
492 | tks.sum_w_o_dz2 = 0;
|
---|
493 | tks.sum_w = 0;
|
---|
494 |
|
---|
495 | Candidate *candidate;
|
---|
496 |
|
---|
497 | fItInputArray->Reset();
|
---|
498 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
499 | {
|
---|
500 | unsigned int discard = 0;
|
---|
501 |
|
---|
502 | double pt = candidate->Momentum.Pt();
|
---|
503 | if(pt<fPtMin || pt>fPtMax) discard = 1;
|
---|
504 |
|
---|
505 | // ------------- Compute cloasest approach Z ----------------
|
---|
506 | double z = candidate->DZ; // [mm]
|
---|
507 |
|
---|
508 | candidate->Zd = candidate->DZ; //Set the cloasest approach z
|
---|
509 | if(fabs(z) > 3*fDzCutOff) discard = 1;
|
---|
510 |
|
---|
511 | // ------------- Compute cloasest approach T ----------------
|
---|
512 | //Asumme pion mass which is the most common particle
|
---|
513 | double M = 0.139570;
|
---|
514 | candidate->Mass = M;
|
---|
515 | double p = pt * sqrt(1 + candidate->CtgTheta*candidate->CtgTheta);
|
---|
516 | double e = sqrt(p*p + M*M);
|
---|
517 |
|
---|
518 | double t = candidate->Position.T()*1.E9/c_light; // from [mm] to [ps]
|
---|
519 | if(t <= -9999) discard = 1; // Means that the time information has not been added
|
---|
520 |
|
---|
521 | // DEBUG Here backpropagete for the whole length and not noly for z. Could improve resolution
|
---|
522 | // double bz = pt * candidate->CtgTheta/e;
|
---|
523 | // t += (z - candidate->Position.Z())*1E9/(c_light*bz);
|
---|
524 |
|
---|
525 | // Use full path Length
|
---|
526 | t -= candidate->L*1E9/(c_light*p/e);
|
---|
527 |
|
---|
528 | candidate->Td = t*1E-9*c_light;
|
---|
529 | if(fabs(t) > 3*fDtCutOff) discard = 1;
|
---|
530 |
|
---|
531 | // auto genp = (Candidate*) candidate->GetCandidates()->At(0);
|
---|
532 | // cout << "Eta: " << candidate->Position.Eta() << endl;
|
---|
533 | // cout << genp->Momentum.Pt() << " -- " << candidate->Momentum.Pt() << endl;
|
---|
534 | // cout << genp->Momentum.Pz() << " -- " << candidate->Momentum.Pz() << endl;
|
---|
535 | // cout << genp->Momentum.P() << " -- " << p << endl;
|
---|
536 | // cout << genp->Momentum.E() << " -- " << e << endl;
|
---|
537 | // cout << Form("bz_true: %.4f -- bz_gen: %.4f", genp->Momentum.Pz()/genp->Momentum.E(), bz) << endl;
|
---|
538 |
|
---|
539 | double dz2_o = candidate->ErrorDZ*candidate->ErrorDZ;
|
---|
540 | dz2_o += fVertexZSize*fVertexZSize;
|
---|
541 | // when needed add beam spot width (x-y)?? mha?
|
---|
542 | dz2_o = 1/dz2_o; //Multipling is faster than dividing all the times
|
---|
543 |
|
---|
544 | double dt2_o = candidate->ErrorT*1.E9/c_light; // [ps]
|
---|
545 | dt2_o *= dt2_o;
|
---|
546 | dt2_o += fVertexTSize*fVertexTSize; // [ps^2]
|
---|
547 | // Ideally we should also add the induced uncertantiy from dz, z_out, pt, ctgthetaand all the other thing used above (total around 5ps). For the moment we compensae using a high value for vertex time.
|
---|
548 | dt2_o = 1/dt2_o;
|
---|
549 |
|
---|
550 | double w;
|
---|
551 | if(fD0CutOff > 0 && candidate->ErrorD0 > 0)
|
---|
552 | {
|
---|
553 | double d0_sig = fabs(candidate->D0/candidate->ErrorD0);
|
---|
554 | w = exp(d0_sig*d0_sig - fD0CutOff*fD0CutOff);
|
---|
555 | w = 1./(1. + w);
|
---|
556 | if (w < 1E-4) discard = 1;
|
---|
557 | }
|
---|
558 | else
|
---|
559 | {
|
---|
560 | w = 1;
|
---|
561 | }
|
---|
562 | candidate->VertexingWeight = w;
|
---|
563 |
|
---|
564 |
|
---|
565 | if(discard)
|
---|
566 | {
|
---|
567 | candidate->ClusterIndex = -1;
|
---|
568 | candidate->InitialPosition.SetT(1E3*1000000*c_light);
|
---|
569 | candidate->InitialPosition.SetZ(1E8);
|
---|
570 | fTrackOutputArray->Add(candidate);
|
---|
571 | }
|
---|
572 | else
|
---|
573 | {
|
---|
574 | tks.sum_w_o_dt2 += w * dt2_o;
|
---|
575 | tks.sum_w_o_dz2 += w * dz2_o;
|
---|
576 | tks.sum_w += w;
|
---|
577 | tks.addItem(z, t, dz2_o, dt2_o, &(*candidate), w, candidate->PID); //PROVA: rimuovi &(*---)
|
---|
578 | }
|
---|
579 |
|
---|
580 | }
|
---|
581 |
|
---|
582 | if(fVerbose > 1)
|
---|
583 | {
|
---|
584 | cout << "----->Filled tracks" << endl;
|
---|
585 | cout << "M z dz t dt w" << endl;
|
---|
586 | for(unsigned int i = 0; i < tks.getSize(); i++)
|
---|
587 | {
|
---|
588 | cout << Form("%d\t%1.1e\t%1.1e\t%1.1e\t%1.1e\t%1.1e", tks.PID[i], tks.z[i], 1/sqrt(tks.dz2_o[i]), tks.t[i], 1/sqrt(tks.dt2_o[i]), tks.w[i]) << endl;
|
---|
589 | }
|
---|
590 | }
|
---|
591 |
|
---|
592 | return;
|
---|
593 | }
|
---|
594 |
|
---|
595 | //------------------------------------------------------------------------------
|
---|
596 | // Compute higher phase transition temperature
|
---|
597 | double VertexFinderDA4D::beta0(tracks_t & tks, vertex_t &vtx)
|
---|
598 | {
|
---|
599 | if(vtx.getSize() != 1)
|
---|
600 | {
|
---|
601 | throw std::invalid_argument( "Unexpected number of vertices" );
|
---|
602 | }
|
---|
603 |
|
---|
604 | unsigned int nt = tks.getSize();
|
---|
605 |
|
---|
606 | //Set vertex position at T=inf as the weighted average of the tracks
|
---|
607 | double sum_wz = 0, sum_wt = 0;
|
---|
608 | for(unsigned int i = 0; i < nt; i++)
|
---|
609 | {
|
---|
610 | sum_wz += tks.w[i] * tks.z[i] * tks.dz2_o[i];
|
---|
611 | sum_wt += tks.w[i] * tks.t[i] * tks.dt2_o[i];
|
---|
612 | }
|
---|
613 | vtx.t[0] = sum_wt / tks.sum_w_o_dt2;
|
---|
614 | vtx.z[0] = sum_wz / tks.sum_w_o_dz2;
|
---|
615 |
|
---|
616 | // Compute the posterior distribution covariance matrix elements
|
---|
617 | double s_zz = 0, s_tt = 0, s_tz = 0;
|
---|
618 | for(unsigned int i = 0; i < nt; i++)
|
---|
619 | {
|
---|
620 | double dz = (tks.z[i] - vtx.z[0]) * tks.dz_o[i];
|
---|
621 | double dt = (tks.t[i] - vtx.t[0]) * tks.dt_o[i];
|
---|
622 |
|
---|
623 | s_zz += tks.w[i] * dz * dz;
|
---|
624 | s_tt += tks.w[i] * dt * dt;
|
---|
625 | s_tz += tks.w[i] * dt * dz;
|
---|
626 | }
|
---|
627 | s_tt /= tks.sum_w;
|
---|
628 | s_zz /= tks.sum_w;
|
---|
629 | s_tz /= tks.sum_w;
|
---|
630 |
|
---|
631 | // Copute the max eighenvalue
|
---|
632 | double beta_c = (s_tt - s_zz)*(s_tt - s_zz) + 4*s_tz*s_tz;
|
---|
633 | beta_c = 1. / (s_tt + s_zz + sqrt(beta_c));
|
---|
634 |
|
---|
635 | double out;
|
---|
636 | if (beta_c < fBetaMax)
|
---|
637 | {
|
---|
638 | // Cool down up to a step before the phase transition
|
---|
639 | out = beta_c * sqrt(fCoolingFactor);
|
---|
640 | }
|
---|
641 | else
|
---|
642 | {
|
---|
643 | out = fBetaMax * fCoolingFactor;
|
---|
644 | }
|
---|
645 |
|
---|
646 | return out;
|
---|
647 | }
|
---|
648 |
|
---|
649 | //------------------------------------------------------------------------------
|
---|
650 | // Compute the new vertexes position and mass (probability) -- mass constrained annealing without noise
|
---|
651 | // Compute and store the posterior covariance matrix elements
|
---|
652 | // Returns the squared sum of changes of vertexex position normalized by the vertex size declared in the init
|
---|
653 | double VertexFinderDA4D::update(double beta, tracks_t &tks, vertex_t &vtx, double rho0)
|
---|
654 | {
|
---|
655 | unsigned int nt = tks.getSize();
|
---|
656 | unsigned int nv = vtx.getSize();
|
---|
657 |
|
---|
658 | //initialize sums
|
---|
659 | double Z_init = rho0 * exp(-beta * fMuOutlayer * fMuOutlayer);
|
---|
660 |
|
---|
661 | // Compute all the energies (aka distances) and normalization partition function
|
---|
662 | vector<double> pk_exp_mBetaE = Compute_pk_exp_mBetaE(beta, vtx, tks, Z_init);
|
---|
663 |
|
---|
664 | double sum_pk = 0;
|
---|
665 | double delta2_max = 0;
|
---|
666 | for (unsigned int k = 0; k < nv; k++)
|
---|
667 | {
|
---|
668 | // Compute the new vertex positions and masses
|
---|
669 | double pk_new = 0;
|
---|
670 | double sw_z = 0, sw_t = 0;
|
---|
671 | // Compute the posterior covariance matrix Elements
|
---|
672 | double szz = 0, stt = 0, stz = 0;
|
---|
673 | double sum_wt = 0, sum_wz = 0;
|
---|
674 | double sum_ptt = 0, sum_pzz = 0, sum_ptz = 0;
|
---|
675 |
|
---|
676 |
|
---|
677 | for (unsigned int i = 0; i < nt; i++)
|
---|
678 | {
|
---|
679 | unsigned int idx = k*nt + i;
|
---|
680 |
|
---|
681 | if(pk_exp_mBetaE[idx] == 0 || tks.Z[i] == 0)
|
---|
682 | {
|
---|
683 | continue;
|
---|
684 | }
|
---|
685 |
|
---|
686 | double p_ygx = pk_exp_mBetaE[idx] / tks.Z[i]; //p(y|x), Gibbs distribution
|
---|
687 | if(std::isnan(p_ygx) || std::isinf(p_ygx) || p_ygx > 1)
|
---|
688 | {
|
---|
689 | cout << Form("%1.6e %1.6e", pk_exp_mBetaE[idx], tks.Z[i]);
|
---|
690 | throw std::invalid_argument(Form("p_ygx is %.8f", p_ygx));
|
---|
691 | }
|
---|
692 | pk_new += tks.w[i] * p_ygx;
|
---|
693 |
|
---|
694 | double wt = tks.w[i] * p_ygx * tks.dt2_o[i];
|
---|
695 | sw_t += wt * tks.t[i];
|
---|
696 | sum_wt += wt;
|
---|
697 |
|
---|
698 | double wz = tks.w[i] * p_ygx * tks.dz2_o[i];
|
---|
699 | sw_z += wz * tks.z[i];
|
---|
700 | sum_wz += wz;
|
---|
701 |
|
---|
702 | // Add the track contribution to the covariance matrix
|
---|
703 | double p_xgy = p_ygx * tks.w[i] / vtx.pk[k];
|
---|
704 | double dt = (tks.t[i] - vtx.t[k]) * tks.dt_o[i];
|
---|
705 | double dz = (tks.z[i] - vtx.z[k]) * tks.dz_o[i];
|
---|
706 |
|
---|
707 | double wtt = p_xgy * tks.dt2_o[i];
|
---|
708 | double wzz = p_xgy * tks.dz2_o[i];
|
---|
709 | double wtz = p_xgy * tks.dt_o[i] * tks.dz_o[i];
|
---|
710 |
|
---|
711 | stt += wtt * dt * dt;
|
---|
712 | szz += wzz * dz * dz;
|
---|
713 | stz += wtz * dt * dz;
|
---|
714 |
|
---|
715 | sum_ptt += wtt;
|
---|
716 | sum_pzz += wzz;
|
---|
717 | sum_ptz += wtz;
|
---|
718 | }
|
---|
719 | if(pk_new == 0)
|
---|
720 | {
|
---|
721 | vtx.removeItem(k);
|
---|
722 | k--;
|
---|
723 | // throw std::invalid_argument(Form("pk_new is %.8f", pk_new));
|
---|
724 | }
|
---|
725 | else
|
---|
726 | {
|
---|
727 | pk_new /= tks.sum_w;
|
---|
728 | sum_pk += pk_new;
|
---|
729 |
|
---|
730 | stt /= sum_ptt;
|
---|
731 | szz /= sum_pzz;
|
---|
732 | stz /= sum_ptz;
|
---|
733 |
|
---|
734 | double new_t = sw_t/sum_wt;
|
---|
735 | double new_z = sw_z/sum_wz;
|
---|
736 | if(std::isnan(new_z) || std::isnan(new_t))
|
---|
737 | {
|
---|
738 | cout << endl << endl;
|
---|
739 | cout << Form("t: %.3e / %.3e", sw_t, sum_wt) << endl;
|
---|
740 | cout << Form("z: %.3e / %.3e", sw_z, sum_wz) << endl;
|
---|
741 | cout << "pk " << k << " " << vtx.pk[k] << endl;
|
---|
742 | throw std::invalid_argument("new_z is nan");
|
---|
743 | }
|
---|
744 |
|
---|
745 | double z_displ = (new_z - vtx.z[k])/fVertexZSize;
|
---|
746 | double t_displ = (new_t - vtx.t[k])/fVertexTSize;
|
---|
747 | double delta2 = z_displ*z_displ + t_displ*t_displ;
|
---|
748 |
|
---|
749 | if (delta2 > delta2_max) delta2_max = delta2;
|
---|
750 |
|
---|
751 | vtx.z[k] = new_z;
|
---|
752 | vtx.t[k] = new_t;
|
---|
753 | vtx.pk[k] = pk_new;
|
---|
754 | vtx.szz[k] = szz;
|
---|
755 | vtx.stt[k] = stt;
|
---|
756 | vtx.stz[k] = stz;
|
---|
757 | }
|
---|
758 | }
|
---|
759 |
|
---|
760 | if(fabs((sum_pk - 1.) > 1E-4))
|
---|
761 | {
|
---|
762 | cout << "sum_pk " << sum_pk << endl;
|
---|
763 | for (unsigned int k = 0; k < nv; k++)
|
---|
764 | {
|
---|
765 | cout << Form("%d: %1.4e", k, vtx.pk[k]) << endl;
|
---|
766 | }
|
---|
767 | throw std::invalid_argument("Sum of masses not unitary");
|
---|
768 | }
|
---|
769 | // if(fVerbose > 3)
|
---|
770 | // {
|
---|
771 | // cout << "===Update over" << endl;
|
---|
772 | // for (unsigned int k = 0; k < nv; k++)
|
---|
773 | // {
|
---|
774 | // cout << k << endl;
|
---|
775 | // cout << "z: " << vtx.z[k] << " , t: " << vtx.t[k] << " , p: " << vtx.pk[k] << endl;
|
---|
776 | // cout << " | " << vtx.szz[k] << " " << vtx.stz[k] << "|" << endl;
|
---|
777 | // cout << " | " << vtx.stz[k] << " " << vtx.stt[k] << "|" << endl << endl;
|
---|
778 | // }
|
---|
779 | // cout << "=======" << endl;
|
---|
780 | // }
|
---|
781 |
|
---|
782 | return delta2_max;
|
---|
783 | }
|
---|
784 |
|
---|
785 | //------------------------------------------------------------------------------
|
---|
786 | // Split critical vertices (beta_c < beta)
|
---|
787 | // Returns true if at least one cluster was split
|
---|
788 | bool VertexFinderDA4D::split(double &beta, vertex_t &vtx, tracks_t & tks)
|
---|
789 | {
|
---|
790 | bool split = false;
|
---|
791 |
|
---|
792 | auto pair_bc_k = vtx.ComputeAllBeta_c(fVerbose);
|
---|
793 |
|
---|
794 | // If minimum beta_c is higher than beta, no split is necessaire
|
---|
795 | if( pair_bc_k.first > beta )
|
---|
796 | {
|
---|
797 | split = false;
|
---|
798 | }
|
---|
799 | else
|
---|
800 | {
|
---|
801 | const unsigned int nv = vtx.getSize();
|
---|
802 | for(unsigned int k = 0; k < nv; k++)
|
---|
803 | {
|
---|
804 | if( fVerbose > 3 )
|
---|
805 | {
|
---|
806 | cout << "vtx " << k << " beta_c = " << vtx.beta_c[k] << endl;
|
---|
807 | }
|
---|
808 | if(vtx.beta_c[k] <= beta)
|
---|
809 | {
|
---|
810 | double z_old = vtx.z[k];
|
---|
811 | double t_old = vtx.t[k];
|
---|
812 | double pk_old = vtx.pk[k];
|
---|
813 |
|
---|
814 | // Compute splitting direction: given by the max eighenvalue eighenvector
|
---|
815 | double zn = (vtx.szz[k] - vtx.stt[k])*(vtx.szz[k] - vtx.stt[k]) + 4*vtx.stz[k]*vtx.stz[k];
|
---|
816 | zn = vtx.szz[k] - vtx.stt[k] + sqrt(zn);
|
---|
817 | double tn = 2*vtx.stz[k];
|
---|
818 | double norm = hypot(zn, tn);
|
---|
819 | tn /= norm;
|
---|
820 | zn /= norm;
|
---|
821 |
|
---|
822 | // Estimate subcluster positions and weight
|
---|
823 | double p1=0, z1=0, t1=0, wz1=0, wt1=0;
|
---|
824 | double p2=0, z2=0, t2=0, wz2=0, wt2=0;
|
---|
825 | const unsigned int nt = tks.getSize();
|
---|
826 | for(unsigned int i=0; i<nt; ++i)
|
---|
827 | {
|
---|
828 | if (tks.Z[i] > 0)
|
---|
829 | {
|
---|
830 | double lr = (tks.t[i] - vtx.t[k]) * tn + (tks.z[i]-vtx.z[k]) * zn;
|
---|
831 | // winner-takes-all, usually overestimates splitting
|
---|
832 | double tl = lr < 0 ? 1.: 0.;
|
---|
833 | double tr = 1. - tl;
|
---|
834 |
|
---|
835 | // soften it, especially at low T
|
---|
836 | // double arg = lr * sqrt(beta * ( zn*zn*tks.dz2_o[i] + tn*tn*tks.dt2_o[i] ) );
|
---|
837 | // if(abs(arg) < 20)
|
---|
838 | // {
|
---|
839 | // double t = exp(-arg);
|
---|
840 | // tl = t/(t+1.);
|
---|
841 | // tr = 1/(t+1.);
|
---|
842 | // }
|
---|
843 |
|
---|
844 | double p = vtx.pk[k] * tks.w[i];
|
---|
845 | p *= exp(-beta * Energy(tks.z[i], vtx.z[k], tks.dz2_o[i], tks.t[i], vtx.t[k], tks.dt2_o[i])) / tks.Z[i];
|
---|
846 | double wt = p*tks.dt2_o[i];
|
---|
847 | double wz = p*tks.dz2_o[i];
|
---|
848 | p1 += p*tl; z1 += wz*tl*tks.z[i]; t1 += wt*tl*tks.t[i]; wz1 += wz*tl; wt1 += wt*tl;
|
---|
849 | p2 += p*tr; z2 += wz*tr*tks.z[i]; t2 += wt*tr*tks.t[i]; wz2 += wz*tr; wt2 += wt*tr;
|
---|
850 | }
|
---|
851 | }
|
---|
852 |
|
---|
853 | if(wz1 > 0 && wt1 > 0 && wz2 > 0 && wt2 > 0)
|
---|
854 | {
|
---|
855 | t1 /= wt1;
|
---|
856 | z1 /= wz1;
|
---|
857 | t2 /= wt2;
|
---|
858 | z2 /= wz2;
|
---|
859 |
|
---|
860 | if( fVerbose > 3 )
|
---|
861 | {
|
---|
862 | double aux = (z1-z2)*(z1-z2)/(fVertexZSize*fVertexZSize) + (t1-t2)*(t1-t2)/(fVertexTSize*fVertexTSize);
|
---|
863 | cout << "weighted split: delta = " << sqrt(aux) << endl;
|
---|
864 | }
|
---|
865 | }
|
---|
866 | else
|
---|
867 | {
|
---|
868 | continue;
|
---|
869 | // throw std::invalid_argument( "0 division" );
|
---|
870 | }
|
---|
871 |
|
---|
872 | while(vtx.NearestCluster(t1, z1) != k || vtx.NearestCluster(t2, z2) != k)
|
---|
873 | {
|
---|
874 | t1 = 0.5 * (t1 + t_old);
|
---|
875 | z1 = 0.5 * (z1 + z_old);
|
---|
876 | t2 = 0.5 * (t2 + t_old);
|
---|
877 | z2 = 0.5 * (z2 + z_old);
|
---|
878 | }
|
---|
879 |
|
---|
880 | // Compute final distance and split if the distance is enough
|
---|
881 | double delta2 = (z1-z2)*(z1-z2)/(fVertexZSize*fVertexZSize) + (t1-t2)*(t1-t2)/(fVertexTSize*fVertexTSize);
|
---|
882 | if(delta2 > fD2Merge)
|
---|
883 | {
|
---|
884 | split = true;
|
---|
885 | vtx.t[k] = t1;
|
---|
886 | vtx.z[k] = z1;
|
---|
887 | vtx.pk[k] = p1 * pk_old/(p1+p2);
|
---|
888 |
|
---|
889 | double new_t = t2;
|
---|
890 | double new_z = z2;
|
---|
891 | double new_pk = p2 * pk_old/(p1+p2);
|
---|
892 |
|
---|
893 | vtx.addItem(new_z, new_t, new_pk);
|
---|
894 |
|
---|
895 | if( fVerbose > 3 )
|
---|
896 | {
|
---|
897 | cout << "===Split happened on vtx " << k << endl;
|
---|
898 | cout << "OLD z: " << z_old << " , t: " << t_old << " , pk: " << pk_old << endl;
|
---|
899 | cout << "NEW+ z: " << vtx.z[k] << " , t: " << vtx.t[k] << " , pk: " << vtx.pk[k] << endl;
|
---|
900 | cout << "NEW- z: " << new_z << " , t: " << new_t << " , pk: " << new_pk << endl;
|
---|
901 | }
|
---|
902 | }
|
---|
903 | }
|
---|
904 | }
|
---|
905 | }
|
---|
906 | return split;
|
---|
907 | }
|
---|
908 |
|
---|
909 |
|
---|
910 | //------------------------------------------------------------------------------
|
---|
911 | // Merge vertexes closer than declared dimensions
|
---|
912 | bool VertexFinderDA4D::merge(vertex_t & vtx, double d2_merge = 2)
|
---|
913 | {
|
---|
914 | bool merged = false;
|
---|
915 |
|
---|
916 | if(vtx.getSize() < 2) return merged;
|
---|
917 |
|
---|
918 | bool last_merge = false;
|
---|
919 | do {
|
---|
920 | double min_d2 = d2_merge;
|
---|
921 | unsigned int k1_min, k2_min;
|
---|
922 | for(unsigned int k1 = 0; k1 < vtx.getSize(); k1++)
|
---|
923 | {
|
---|
924 | for(unsigned int k2 = k1+1; k2 < vtx.getSize();k2++)
|
---|
925 | {
|
---|
926 | double d2_tmp = vtx.DistanceSquare(k1, k2);
|
---|
927 | if(d2_tmp < min_d2)
|
---|
928 | {
|
---|
929 | min_d2 = d2_tmp;
|
---|
930 | k1_min = k1;
|
---|
931 | k2_min = k2;
|
---|
932 | }
|
---|
933 | }
|
---|
934 | }
|
---|
935 |
|
---|
936 | if(min_d2 < d2_merge)
|
---|
937 | {
|
---|
938 | vtx.mergeItems(k1_min, k2_min);
|
---|
939 | last_merge = true;
|
---|
940 | merged = true;
|
---|
941 | }
|
---|
942 | else last_merge = false;
|
---|
943 | } while(last_merge);
|
---|
944 |
|
---|
945 | return merged;
|
---|
946 | }
|
---|
947 |
|
---|
948 | // -----------------------------------------------------------------------------
|
---|
949 | // Compute all the energies and set the partition function normalization for each track
|
---|
950 | vector<double> VertexFinderDA4D::Compute_pk_exp_mBetaE(double beta, vertex_t &vtx, tracks_t &tks, double Z_init)
|
---|
951 | {
|
---|
952 | unsigned int nt = tks.getSize();
|
---|
953 | unsigned int nv = vtx.getSize();
|
---|
954 |
|
---|
955 | vector<double> pk_exp_mBetaE(nt * nv);
|
---|
956 | for (unsigned int k = 0; k < nv; k++)
|
---|
957 | {
|
---|
958 | for (unsigned int i = 0; i < nt; i++)
|
---|
959 | {
|
---|
960 | if(k == 0) tks.Z[i] = Z_init;
|
---|
961 |
|
---|
962 | double aux = Energy(tks.z[i], vtx.z[k], tks.dz2_o[i], tks.t[i], vtx.t[k], tks.dt2_o[i]);
|
---|
963 | aux = vtx.pk[k] * exp(-beta * aux);
|
---|
964 | // if(aux < 1E-10) continue;
|
---|
965 | tks.Z[i] += aux;
|
---|
966 |
|
---|
967 | unsigned int idx = k*nt + i;
|
---|
968 | pk_exp_mBetaE[idx] = aux;
|
---|
969 | }
|
---|
970 | }
|
---|
971 | return pk_exp_mBetaE;
|
---|
972 | }
|
---|
973 |
|
---|
974 | //------------------------------------------------------------------------------
|
---|
975 | // Eliminate clusters with only one significant/unique track
|
---|
976 | bool VertexFinderDA4D::purge(vertex_t & vtx, tracks_t & tks, double & rho0, const double beta, double min_prob, double min_trk)
|
---|
977 | {
|
---|
978 | const unsigned int nv = vtx.getSize();
|
---|
979 | const unsigned int nt = tks.getSize();
|
---|
980 |
|
---|
981 | if (nv < 2)
|
---|
982 | return false;
|
---|
983 |
|
---|
984 | double sumpmin = nt;
|
---|
985 | unsigned int k0 = nv;
|
---|
986 |
|
---|
987 | int nUnique = 0;
|
---|
988 | double sump = 0;
|
---|
989 |
|
---|
990 | double Z_init = rho0 * exp(-beta * fMuOutlayer * fMuOutlayer); // Add fDtCutOff here toghether with this
|
---|
991 | vector<double> pk_exp_mBetaE = Compute_pk_exp_mBetaE(beta, vtx, tks, Z_init);
|
---|
992 |
|
---|
993 | for (unsigned int k = 0; k < nv; ++k) {
|
---|
994 |
|
---|
995 | nUnique = 0;
|
---|
996 | sump = 0;
|
---|
997 |
|
---|
998 | double pmax = vtx.pk[k] / (vtx.pk[k] + rho0 * exp(-beta * fMuOutlayer* fMuOutlayer));
|
---|
999 | double pcut = min_prob * pmax;
|
---|
1000 |
|
---|
1001 | for (unsigned int i = 0; i < nt; ++i) {
|
---|
1002 | unsigned int idx = k*nt + i;
|
---|
1003 |
|
---|
1004 | if(pk_exp_mBetaE[idx] == 0 || tks.Z[i] == 0)
|
---|
1005 | {
|
---|
1006 | continue;
|
---|
1007 | }
|
---|
1008 |
|
---|
1009 | double p = pk_exp_mBetaE[idx] / tks.Z[i];
|
---|
1010 | sump += p;
|
---|
1011 | if( ( p > pcut ) & ( tks.w[i] > 0 ) ) nUnique++;
|
---|
1012 | }
|
---|
1013 |
|
---|
1014 | if ((nUnique < min_trk) && (sump < sumpmin)) {
|
---|
1015 | sumpmin = sump;
|
---|
1016 | k0 = k;
|
---|
1017 | }
|
---|
1018 |
|
---|
1019 | }
|
---|
1020 |
|
---|
1021 | if (k0 != nv) {
|
---|
1022 | if (fVerbose > 5) {
|
---|
1023 | std::cout << Form("eliminating prototype at z = %.3f mm, t = %.0f ps", vtx.z[k0], vtx.t[k0]) << " with sump=" << sumpmin
|
---|
1024 | << " rho*nt =" << vtx.pk[k0]*nt
|
---|
1025 | << endl;
|
---|
1026 | }
|
---|
1027 | vtx.removeItem(k0);
|
---|
1028 | return true;
|
---|
1029 | } else {
|
---|
1030 | return false;
|
---|
1031 | }
|
---|
1032 | }
|
---|