Fork me on GitHub

source: git/modules/VertexFinderDA4D.cc@ 6fc566b

Timing
Last change on this file since 6fc566b was 84a097e, checked in by Michele Selvaggi <michele.selvaggi@…>, 5 years ago

first iteration of adding timing in Delphes

  • Property mode set to 100644
File size: 39.3 KB
RevLine 
[29d662e]1/** \class VertexFinderDA4D
2 *
3 * Cluster vertices from tracks using deterministic annealing and timing information
4 *
[84a097e]5 * \authors O. Cerri
[29d662e]6 *
7 */
8
[84a097e]9
[29d662e]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
[341014c]16#include "ExRootAnalysis/ExRootResult.h"
[84a097e]17#include "ExRootAnalysis/ExRootFilter.h"
18#include "ExRootAnalysis/ExRootClassifier.h"
[29d662e]19
[84a097e]20#include "TMath.h"
21#include "TString.h"
[341014c]22#include "TFormula.h"
[84a097e]23#include "TRandom3.h"
24#include "TObjArray.h"
25#include "TDatabasePDG.h"
[29d662e]26#include "TLorentzVector.h"
27#include "TMatrixT.h"
[84a097e]28#include "TLatex.h"
[29d662e]29#include "TVector3.h"
30
[84a097e]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>
[95b4e9f]41#include <algorithm>
[341014c]42#include <stdexcept>
[84a097e]43#include <iostream>
[95b4e9f]44#include <vector>
45
46using namespace std;
47
[84a097e]48namespace vtx_DAZT
[4154bbd]49{
[84a097e]50 static const Double_t c_light = 2.99792458e+8; // [m/s]
[29d662e]51}
[84a097e]52using namespace vtx_DAZT;
[29d662e]53
54//------------------------------------------------------------------------------
55
[84a097e]56VertexFinderDA4D::VertexFinderDA4D()
[29d662e]57{
[84a097e]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;
[29d662e]74}
75
76//------------------------------------------------------------------------------
77
78VertexFinderDA4D::~VertexFinderDA4D()
79{
80}
81
82//------------------------------------------------------------------------------
83
84void VertexFinderDA4D::Init()
85{
[84a097e]86 fVerbose = GetInt("Verbose", 0);
87
88 fMaxIterations = GetInt("MaxIterations", 100);
89 fMaxVertexNumber = GetInt("MaxVertexNumber", 500);
[29d662e]90
[84a097e]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", ".");
[29d662e]114
115 fInputArray = ImportArray(GetString("InputArray", "TrackSmearing/tracks"));
116 fItInputArray = fInputArray->MakeIterator();
117
[84a097e]118 fTrackOutputArray = ExportArray(GetString("TrackOutputArray", "tracks"));
[29d662e]119 fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
[84a097e]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 }
[29d662e]141}
142
143//------------------------------------------------------------------------------
144
145void VertexFinderDA4D::Finish()
146{
147 if(fItInputArray) delete fItInputArray;
148}
149
150//------------------------------------------------------------------------------
151
152void VertexFinderDA4D::Process()
153{
154 fInputArray->Sort();
155
[84a097e]156 if (fVerbose)
[29d662e]157 {
[84a097e]158 cout<< endl << " Start processing vertices with VertexFinderDA4D" << endl;
159 cout<<" Found "<<fInputArray->GetEntriesFast()<<" input tracks"<<endl;
[29d662e]160 }
161
[84a097e]162 // clusterize tracks
163 TObjArray *ClusterArray = new TObjArray;
164 clusterize(*ClusterArray);
[29d662e]165
[84a097e]166 if(fVerbose>10)
[341014c]167 {
[84a097e]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();
[341014c]188 }
[29d662e]189
[84a097e]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();
[29d662e]194 ItClusterArray->Reset();
[84a097e]195 Candidate *candidate;
196 unsigned int k = 0;
197 while((candidate = static_cast<Candidate*>(ItClusterArray->Next())))
[29d662e]198 {
[84a097e]199 if(fVerbose)
[341014c]200 {
[84a097e]201 cout << Form("Cluster %d has %d tracks ", k, candidate->GetCandidates()->GetEntriesFast()) << endl;
[341014c]202 }
[84a097e]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;
[341014c]219
[84a097e]220 double wz = track->VertexingWeight/(track->ErrorDZ*track->ErrorDZ);
221 double wt = track->VertexingWeight/(track->ErrorT*track->ErrorT);
[341014c]222
[84a097e]223 sum_Dt_2 += wt*dt*dt;
224 sum_Dz_2 += wz*dz*dz;
225 sum_wt += wt;
226 sum_wz += wz;
227 }
[29d662e]228
[84a097e]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 }
[341014c]238
[84a097e]239 fVertexOutputArray->Add(candidate);
240 k++;
[29d662e]241 }
[84a097e]242 }// end of cluster loop
[29d662e]243
244 delete ClusterArray;
245}
246
247//------------------------------------------------------------------------------
248
[84a097e]249void VertexFinderDA4D::clusterize(TObjArray &clusters)
[29d662e]250{
[84a097e]251 tracks_t tks;
252 fill(tks);
253 unsigned int nt=tks.getSize();
[341014c]254 if(fVerbose)
255 {
[84a097e]256 cout << "Tracks added: " << nt << endl;
[29d662e]257 }
[84a097e]258 if (nt == 0) return;
[29d662e]259
260
261
[84a097e]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);
[29d662e]267
[84a097e]268 // Fit the vertex at T=inf and return the starting temperature
269 double beta=beta0(tks, vtx);
[29d662e]270
[84a097e]271 if( fVerbose > 1 )
[341014c]272 {
[84a097e]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;
[341014c]275 }
[29d662e]276
[84a097e]277 if( fVerbose > 10 ) plot_status(beta, vtx, tks, 0, "Ast");
[29d662e]278
[84a097e]279 if( fVerbose > 2){cout << "Cool down untill reaching the temperature to finish increasing the number of vertexes" << endl;}
[29d662e]280
[84a097e]281 double rho0=0.0; // start with no outlier rejection
282
283 unsigned int last_round = 0;
284 while(last_round < 2)
[341014c]285 {
[29d662e]286
[84a097e]287 unsigned int niter=0;
288 double delta2 = 0;
289 do {
290 delta2 = update(beta, tks, vtx, rho0);
[29d662e]291
[84a097e]292 if( fVerbose > 10 ) plot_status(beta, vtx, tks, niter, "Bup");
293 if (fVerbose > 3)
294 {
295 cout << "Update " << niter << " : " << delta2 << endl;
296 }
297 niter++;
298 }
299 while (delta2 > fD2UpdateLim && niter < fMaxIterations);
[29d662e]300
[4154bbd]301
[84a097e]302 unsigned int n_it = 0;
303 while(merge(vtx, fD2Merge) && n_it < fMaxIterations)
[4154bbd]304 {
[84a097e]305 unsigned int niter=0;
306 double delta2 = 0;
307 do {
308 delta2 = update(beta, tks, vtx, rho0);
309 niter++;
310 }
311 while (delta2 > fD2UpdateLim && niter < fMaxIterations);
312 n_it++;
[4154bbd]313
[84a097e]314 if( fVerbose > 10 ) plot_status(beta, vtx, tks, n_it, "Cme");
315 }
[4154bbd]316
[84a097e]317 beta /= fCoolingFactor;
318
319 if( beta < fBetaStop )
320 {
321 split(beta, vtx, tks);
322 if( fVerbose > 10 ) plot_status(beta, vtx, tks, 0, "Asp");
[4154bbd]323 }
324 else
325 {
[84a097e]326 beta = fBetaStop;
327 last_round++;
[4154bbd]328 }
329
[84a097e]330 if(fVerbose > 3)
[4154bbd]331 {
[84a097e]332 cout << endl << endl << " ----- Beta = " << beta << " --------" << endl;
333 cout << "Nv: " << vtx.getSize() << endl;
[4154bbd]334 }
335 }
[29d662e]336
[84a097e]337 if( fVerbose > 4)
[29d662e]338 {
[84a097e]339 for(unsigned int k = 0; k < vtx.getSize(); k++)
[341014c]340 {
[84a097e]341 cout << Form("Vertex %d next beta_c = %.3f", k, vtx.beta_c[k]) << endl;
[341014c]342 }
[29d662e]343 }
344
[84a097e]345 if(fVerbose > 2) {cout << "Adiabatic switch on of outlayr rejection" << endl;}
346 rho0 = 1./nt;
347 const double N_cycles = 10;
348 for(unsigned int f = 1; f <= N_cycles; f++)
[341014c]349 {
[84a097e]350 unsigned int niter=0;
351 double delta2 = 0;
352 do {
353 delta2 = update(beta, tks, vtx, rho0 * f/N_cycles);
354 niter++;
355 }
356 while (delta2 > 0.3*fD2UpdateLim && niter < fMaxIterations);
357 if( fVerbose > 10 ) plot_status(beta, vtx, tks, f, "Dadout");
[341014c]358 }
[29d662e]359
[84a097e]360 do {
361 beta /= fCoolingFactor;
362 if(beta > fBetaPurge) beta = fBetaPurge;
363 unsigned int i_pu = 0;
364 for(int min_trk = 2; min_trk<=fMinNTrack; min_trk++)
[341014c]365 {
[84a097e]366 while( purge(vtx, tks, rho0, beta, fMinTrackProb, min_trk) )
[341014c]367 {
[84a097e]368 unsigned int niter=0;
369 double delta2 = 0;
370 do {
371 delta2 = update(beta, tks, vtx, rho0);
372 niter++;
373 }
374 while (delta2 > fD2UpdateLim && niter < fMaxIterations);
375 if( fVerbose > 10 ) plot_status(beta, vtx, tks, i_pu, Form("Eprg%d",min_trk));
376 i_pu++;
[341014c]377 }
[29d662e]378 }
379
[84a097e]380 unsigned int n_it = 0;
381 while(merge(vtx, fD2Merge) && n_it < fMaxIterations)
[341014c]382 {
[84a097e]383 unsigned int niter=0;
384 double delta2 = 0;
385 do {
386 delta2 = update(beta, tks, vtx, rho0);
387 niter++;
[341014c]388 }
[84a097e]389 while (delta2 > fD2UpdateLim && niter < fMaxIterations);
390 n_it++;
391
392 if( fVerbose > 10 ) plot_status(beta, vtx, tks, n_it, "Cme");
[341014c]393 }
[84a097e]394 } while( beta < fBetaPurge );
395
396
397 if(fVerbose > 2){cout << "Cooldown untill the limit before assigning track to vertices" << endl;}
398 last_round = 0;
399 while(last_round < 2)
[341014c]400 {
[84a097e]401 unsigned int niter=0;
402 double delta2 = 0;
403 do {
404 delta2 = update(beta, tks, vtx, rho0);
405 niter++;
406 if( fVerbose > 10 ) plot_status(beta, vtx, tks, 0, "Bup");
[341014c]407 }
[84a097e]408 while (delta2 > 0.3*fD2UpdateLim && niter < fMaxIterations);
409
410 beta /= fCoolingFactor;
411 if ( beta >= fBetaMax )
[341014c]412 {
[84a097e]413 beta = fBetaMax;
414 last_round++;
[341014c]415 }
[29d662e]416 }
417
418
[84a097e]419 // Build the cluster candidates
420 for(unsigned int k = 0; k < vtx.getSize(); k++)
[341014c]421 {
[84a097e]422 DelphesFactory *factory = GetFactory();
423 Candidate * candidate = factory->NewCandidate();
424
425 candidate->ClusterIndex = k;
426 candidate->Position.SetXYZT(0.0, 0.0, vtx.z[k] , vtx.t[k]*1E-9*c_light);
427 candidate->PositionError.SetXYZT(0.0, 0.0, fVertexZSize , fVertexTSize*1E-9*c_light);
428 candidate->SumPT2 = 0;
429 candidate->SumPt = 0;
430 candidate->ClusterNDF = 0;
431
432 clusters.Add(candidate);
[341014c]433 }
[29d662e]434
[84a097e]435
436 // Assign each track to the most probable vertex
437 double Z_init = rho0 * exp(-beta * fMuOutlayer * fMuOutlayer); // Add fDtCutOff here toghether with this
438 vector<double> pk_exp_mBetaE = Compute_pk_exp_mBetaE(beta, vtx, tks, Z_init);
439 for(unsigned int i = 0; i< tks.getSize(); i++)
[341014c]440 {
[84a097e]441 if(tks.w[i] <= 0) continue;
442
443 double p_max = 0;
444 unsigned int k_max = 0;
445
446 for(unsigned int k = 0; k < vtx.getSize(); k++)
[341014c]447 {
[84a097e]448 unsigned int idx = k*nt + i;
449 if(pk_exp_mBetaE[idx] == 0 || tks.Z[i] == 0 || vtx.pk[k] == 0)
450 {
451 continue;
452 }
453
454 double pv_max = vtx.pk[k] / (vtx.pk[k] + rho0 * exp(-beta * fMuOutlayer* fMuOutlayer));
455 double p = pk_exp_mBetaE[idx] / tks.Z[i];
456
457 p /= pv_max;
458
459 if(p > p_max)
[341014c]460 {
[84a097e]461 p_max = p;
462 k_max = k;
[341014c]463 }
464 }
[84a097e]465
466 if(p_max > fMinTrackProb)
[341014c]467 {
[84a097e]468 tks.tt[i]->ClusterIndex = k_max;
469 tks.tt[i]->InitialPosition.SetT(1E-9*vtx.t[k_max]*c_light);
470 tks.tt[i]->InitialPosition.SetZ(vtx.z[k_max]);
471
472 ((Candidate *) clusters.At(k_max))->AddCandidate(tks.tt[i]);
473 ((Candidate *) clusters.At(k_max))->SumPT2 += tks.tt[i]->Momentum.Pt()*tks.tt[i]->Momentum.Pt();
474 ((Candidate *) clusters.At(k_max))->SumPt += tks.tt[i]->Momentum.Pt();
475 ((Candidate *) clusters.At(k_max))->ClusterNDF += 1;
[29d662e]476 }
[84a097e]477 else
478 {
479 tks.tt[i]->ClusterIndex = -1;
480 tks.tt[i]->InitialPosition.SetT(1E3*1000000*c_light);
481 tks.tt[i]->InitialPosition.SetZ(1E8);
482 }
483 fTrackOutputArray->Add(tks.tt[i]);
[29d662e]484 }
485
[84a097e]486 if(fVerbose > 10) plot_status_end(vtx, tks);
[29d662e]487
[84a097e]488}
[29d662e]489
[84a097e]490//------------------------------------------------------------------------------
491// Definition of the distance metrci between track and vertex
492double VertexFinderDA4D::Energy(double t_z, double v_z, double dz2_o, double t_t, double v_t, double dt2_o)
493{
494 return (t_z - v_z)*(t_z - v_z)* dz2_o + (t_t - v_t)*(t_t - v_t)*dt2_o;
495}
[29d662e]496
[84a097e]497//------------------------------------------------------------------------------
498// Fill tks with the input candidates array
499void VertexFinderDA4D::fill(tracks_t &tks)
500{
501 tks.sum_w_o_dt2 = 0;
502 tks.sum_w_o_dz2 = 0;
503 tks.sum_w = 0;
504
505 Candidate *candidate;
506
507 fItInputArray->Reset();
508 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
[341014c]509 {
[84a097e]510 unsigned int discard = 0;
511
512 double pt = candidate->Momentum.Pt();
513 if(pt<fPtMin || pt>fPtMax) discard = 1;
514
515 // ------------- Compute cloasest approach Z ----------------
516 double z = candidate->DZ; // [mm]
517
518 candidate->Zd = candidate->DZ; //Set the cloasest approach z
519 if(fabs(z) > 3*fDzCutOff) discard = 1;
520
521 // ------------- Compute cloasest approach T ----------------
522 //Asumme pion mass which is the most common particle
523 double M = 0.139570;
524 candidate->Mass = M;
525 double p = pt * sqrt(1 + candidate->CtgTheta*candidate->CtgTheta);
526 double e = sqrt(p*p + M*M);
527
528 double t = candidate->Position.T()*1.E9/c_light; // from [mm] to [ps]
529 if(t <= -9999) discard = 1; // Means that the time information has not been added
530
531 // DEBUG Here backpropagete for the whole length and not noly for z. Could improve resolution
532 // double bz = pt * candidate->CtgTheta/e;
533 // t += (z - candidate->Position.Z())*1E9/(c_light*bz);
534
535 // Use full path Length
536 t -= candidate->L*1E9/(c_light*p/e);
537
538 candidate->Td = t*1E-9*c_light;
539 if(fabs(t) > 3*fDtCutOff) discard = 1;
540
541 // auto genp = (Candidate*) candidate->GetCandidates()->At(0);
542 // cout << "Eta: " << candidate->Position.Eta() << endl;
543 // cout << genp->Momentum.Pt() << " -- " << candidate->Momentum.Pt() << endl;
544 // cout << genp->Momentum.Pz() << " -- " << candidate->Momentum.Pz() << endl;
545 // cout << genp->Momentum.P() << " -- " << p << endl;
546 // cout << genp->Momentum.E() << " -- " << e << endl;
547 // cout << Form("bz_true: %.4f -- bz_gen: %.4f", genp->Momentum.Pz()/genp->Momentum.E(), bz) << endl;
548
549 double dz2_o = candidate->ErrorDZ*candidate->ErrorDZ;
550 dz2_o += fVertexZSize*fVertexZSize;
551 // when needed add beam spot width (x-y)?? mha?
552 dz2_o = 1/dz2_o; //Multipling is faster than dividing all the times
553
554 double dt2_o = candidate->ErrorT*1.E9/c_light; // [ps]
555 dt2_o *= dt2_o;
556 dt2_o += fVertexTSize*fVertexTSize; // [ps^2]
557 // 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.
558 dt2_o = 1/dt2_o;
559
560 double w;
561 if(fD0CutOff > 0 && candidate->ErrorD0 > 0)
[341014c]562 {
[84a097e]563 double d0_sig = fabs(candidate->D0/candidate->ErrorD0);
564 w = exp(d0_sig*d0_sig - fD0CutOff*fD0CutOff);
565 w = 1./(1. + w);
566 if (w < 1E-4) discard = 1;
[29d662e]567 }
[84a097e]568 else
569 {
570 w = 1;
571 }
572 candidate->VertexingWeight = w;
[29d662e]573
574
[84a097e]575 if(discard)
[341014c]576 {
[84a097e]577 candidate->ClusterIndex = -1;
578 candidate->InitialPosition.SetT(1E3*1000000*c_light);
579 candidate->InitialPosition.SetZ(1E8);
580 fTrackOutputArray->Add(candidate);
581 }
582 else
583 {
584 tks.sum_w_o_dt2 += w * dt2_o;
585 tks.sum_w_o_dz2 += w * dz2_o;
586 tks.sum_w += w;
587 tks.addItem(z, t, dz2_o, dt2_o, &(*candidate), w, candidate->PID); //PROVA: rimuovi &(*---)
588 }
[29d662e]589
[84a097e]590 }
[29d662e]591
[84a097e]592 if(fVerbose > 1)
593 {
594 cout << "----->Filled tracks" << endl;
595 cout << "M z dz t dt w" << endl;
596 for(unsigned int i = 0; i < tks.getSize(); i++)
597 {
598 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;
[29d662e]599 }
600 }
601
[84a097e]602 return;
[29d662e]603}
604
605//------------------------------------------------------------------------------
[84a097e]606// Compute higher phase transition temperature
607double VertexFinderDA4D::beta0(tracks_t & tks, vertex_t &vtx)
[4154bbd]608{
[84a097e]609 if(vtx.getSize() != 1)
610 {
611 throw std::invalid_argument( "Unexpected number of vertices" );
612 }
[29d662e]613
[84a097e]614 unsigned int nt = tks.getSize();
[29d662e]615
[84a097e]616 //Set vertex position at T=inf as the weighted average of the tracks
617 double sum_wz = 0, sum_wt = 0;
618 for(unsigned int i = 0; i < nt; i++)
[341014c]619 {
[84a097e]620 sum_wz += tks.w[i] * tks.z[i] * tks.dz2_o[i];
621 sum_wt += tks.w[i] * tks.t[i] * tks.dt2_o[i];
[341014c]622 }
[84a097e]623 vtx.t[0] = sum_wt / tks.sum_w_o_dt2;
624 vtx.z[0] = sum_wz / tks.sum_w_o_dz2;
[29d662e]625
[84a097e]626 // Compute the posterior distribution covariance matrix elements
627 double s_zz = 0, s_tt = 0, s_tz = 0;
628 for(unsigned int i = 0; i < nt; i++)
[341014c]629 {
[84a097e]630 double dz = (tks.z[i] - vtx.z[0]) * tks.dz_o[i];
631 double dt = (tks.t[i] - vtx.t[0]) * tks.dt_o[i];
632
633 s_zz += tks.w[i] * dz * dz;
634 s_tt += tks.w[i] * dt * dt;
635 s_tz += tks.w[i] * dt * dz;
[29d662e]636 }
[84a097e]637 s_tt /= tks.sum_w;
638 s_zz /= tks.sum_w;
639 s_tz /= tks.sum_w;
640
641 // Copute the max eighenvalue
642 double beta_c = (s_tt - s_zz)*(s_tt - s_zz) + 4*s_tz*s_tz;
643 beta_c = 1. / (s_tt + s_zz + sqrt(beta_c));
644
645 double out;
646 if (beta_c < fBetaMax)
[341014c]647 {
[84a097e]648 // Cool down up to a step before the phase transition
649 out = beta_c * sqrt(fCoolingFactor);
[29d662e]650 }
[84a097e]651 else
[341014c]652 {
[84a097e]653 out = fBetaMax * fCoolingFactor;
[29d662e]654 }
655
[84a097e]656 return out;
657}
[29d662e]658
[84a097e]659//------------------------------------------------------------------------------
660// Compute the new vertexes position and mass (probability) -- mass constrained annealing without noise
661// Compute and store the posterior covariance matrix elements
662// Returns the squared sum of changes of vertexex position normalized by the vertex size declared in the init
663double VertexFinderDA4D::update(double beta, tracks_t &tks, vertex_t &vtx, double rho0)
664{
665 unsigned int nt = tks.getSize();
666 unsigned int nv = vtx.getSize();
667
668 //initialize sums
669 double Z_init = rho0 * exp(-beta * fMuOutlayer * fMuOutlayer);
670
671 // Compute all the energies (aka distances) and normalization partition function
672 vector<double> pk_exp_mBetaE = Compute_pk_exp_mBetaE(beta, vtx, tks, Z_init);
673
674 double sum_pk = 0;
675 double delta2_max = 0;
676 for (unsigned int k = 0; k < nv; k++)
[341014c]677 {
[84a097e]678 // Compute the new vertex positions and masses
679 double pk_new = 0;
680 double sw_z = 0, sw_t = 0;
681 // Compute the posterior covariance matrix Elements
682 double szz = 0, stt = 0, stz = 0;
683 double sum_wt = 0, sum_wz = 0;
684 double sum_ptt = 0, sum_pzz = 0, sum_ptz = 0;
[4154bbd]685
[84a097e]686
687 for (unsigned int i = 0; i < nt; i++)
[341014c]688 {
[84a097e]689 unsigned int idx = k*nt + i;
690
691 if(pk_exp_mBetaE[idx] == 0 || tks.Z[i] == 0)
[341014c]692 {
[84a097e]693 continue;
[341014c]694 }
[84a097e]695
696 double p_ygx = pk_exp_mBetaE[idx] / tks.Z[i]; //p(y|x), Gibbs distribution
697 if(std::isnan(p_ygx) || std::isinf(p_ygx) || p_ygx > 1)
[341014c]698 {
[84a097e]699 cout << Form("%1.6e %1.6e", pk_exp_mBetaE[idx], tks.Z[i]);
700 throw std::invalid_argument(Form("p_ygx is %.8f", p_ygx));
[29d662e]701 }
[84a097e]702 pk_new += tks.w[i] * p_ygx;
[29d662e]703
[84a097e]704 double wt = tks.w[i] * p_ygx * tks.dt2_o[i];
705 sw_t += wt * tks.t[i];
706 sum_wt += wt;
[29d662e]707
[84a097e]708 double wz = tks.w[i] * p_ygx * tks.dz2_o[i];
709 sw_z += wz * tks.z[i];
710 sum_wz += wz;
[29d662e]711
[84a097e]712 // Add the track contribution to the covariance matrix
713 double p_xgy = p_ygx * tks.w[i] / vtx.pk[k];
714 double dt = (tks.t[i] - vtx.t[k]) * tks.dt_o[i];
715 double dz = (tks.z[i] - vtx.z[k]) * tks.dz_o[i];
[29d662e]716
[84a097e]717 double wtt = p_xgy * tks.dt2_o[i];
718 double wzz = p_xgy * tks.dz2_o[i];
719 double wtz = p_xgy * tks.dt_o[i] * tks.dz_o[i];
[29d662e]720
[84a097e]721 stt += wtt * dt * dt;
722 szz += wzz * dz * dz;
723 stz += wtz * dt * dz;
[29d662e]724
[84a097e]725 sum_ptt += wtt;
726 sum_pzz += wzz;
727 sum_ptz += wtz;
[29d662e]728 }
[84a097e]729 if(pk_new == 0)
[341014c]730 {
[84a097e]731 vtx.removeItem(k);
732 k--;
733 // throw std::invalid_argument(Form("pk_new is %.8f", pk_new));
[341014c]734 }
735 else
736 {
[84a097e]737 pk_new /= tks.sum_w;
738 sum_pk += pk_new;
[29d662e]739
[84a097e]740 stt /= sum_ptt;
741 szz /= sum_pzz;
742 stz /= sum_ptz;
[29d662e]743
[84a097e]744 double new_t = sw_t/sum_wt;
745 double new_z = sw_z/sum_wz;
746 if(std::isnan(new_z) || std::isnan(new_t))
747 {
748 cout << endl << endl;
749 cout << Form("t: %.3e / %.3e", sw_t, sum_wt) << endl;
750 cout << Form("z: %.3e / %.3e", sw_z, sum_wz) << endl;
751 cout << "pk " << k << " " << vtx.pk[k] << endl;
752 throw std::invalid_argument("new_z is nan");
753 }
754
755 double z_displ = (new_z - vtx.z[k])/fVertexZSize;
756 double t_displ = (new_t - vtx.t[k])/fVertexTSize;
757 double delta2 = z_displ*z_displ + t_displ*t_displ;
758
759 if (delta2 > delta2_max) delta2_max = delta2;
760
761 vtx.z[k] = new_z;
762 vtx.t[k] = new_t;
763 vtx.pk[k] = pk_new;
764 vtx.szz[k] = szz;
765 vtx.stt[k] = stt;
766 vtx.stz[k] = stz;
[341014c]767 }
[84a097e]768 }
769
770 if(fabs((sum_pk - 1.) > 1E-4))
771 {
772 cout << "sum_pk " << sum_pk << endl;
773 for (unsigned int k = 0; k < nv; k++)
[341014c]774 {
[84a097e]775 cout << Form("%d: %1.4e", k, vtx.pk[k]) << endl;
[29d662e]776 }
[84a097e]777 throw std::invalid_argument("Sum of masses not unitary");
[29d662e]778 }
[84a097e]779 // if(fVerbose > 3)
780 // {
781 // cout << "===Update over" << endl;
782 // for (unsigned int k = 0; k < nv; k++)
783 // {
784 // cout << k << endl;
785 // cout << "z: " << vtx.z[k] << " , t: " << vtx.t[k] << " , p: " << vtx.pk[k] << endl;
786 // cout << " | " << vtx.szz[k] << " " << vtx.stz[k] << "|" << endl;
787 // cout << " | " << vtx.stz[k] << " " << vtx.stt[k] << "|" << endl << endl;
788 // }
789 // cout << "=======" << endl;
790 // }
[29d662e]791
[84a097e]792 return delta2_max;
[29d662e]793}
794
795//------------------------------------------------------------------------------
[84a097e]796// Split critical vertices (beta_c < beta)
797// Returns true if at least one cluster was split
798bool VertexFinderDA4D::split(double &beta, vertex_t &vtx, tracks_t & tks)
[4154bbd]799{
[84a097e]800 bool split = false;
[29d662e]801
[84a097e]802 auto pair_bc_k = vtx.ComputeAllBeta_c(fVerbose);
[29d662e]803
[84a097e]804 // If minimum beta_c is higher than beta, no split is necessaire
805 if( pair_bc_k.first > beta )
[341014c]806 {
[84a097e]807 split = false;
[29d662e]808 }
[84a097e]809 else
[341014c]810 {
[84a097e]811 const unsigned int nv = vtx.getSize();
812 for(unsigned int k = 0; k < nv; k++)
[341014c]813 {
[84a097e]814 if( fVerbose > 3 )
[341014c]815 {
[84a097e]816 cout << "vtx " << k << " beta_c = " << vtx.beta_c[k] << endl;
[29d662e]817 }
[84a097e]818 if(vtx.beta_c[k] <= beta)
819 {
820 double z_old = vtx.z[k];
821 double t_old = vtx.t[k];
822 double pk_old = vtx.pk[k];
823
824 // Compute splitting direction: given by the max eighenvalue eighenvector
825 double zn = (vtx.szz[k] - vtx.stt[k])*(vtx.szz[k] - vtx.stt[k]) + 4*vtx.stz[k]*vtx.stz[k];
826 zn = vtx.szz[k] - vtx.stt[k] + sqrt(zn);
827 double tn = 2*vtx.stz[k];
828 double norm = hypot(zn, tn);
829 tn /= norm;
830 zn /= norm;
831
832 // Estimate subcluster positions and weight
833 double p1=0, z1=0, t1=0, wz1=0, wt1=0;
834 double p2=0, z2=0, t2=0, wz2=0, wt2=0;
835 const unsigned int nt = tks.getSize();
836 for(unsigned int i=0; i<nt; ++i)
837 {
838 if (tks.Z[i] > 0)
839 {
840 double lr = (tks.t[i] - vtx.t[k]) * tn + (tks.z[i]-vtx.z[k]) * zn;
841 // winner-takes-all, usually overestimates splitting
842 double tl = lr < 0 ? 1.: 0.;
843 double tr = 1. - tl;
844
845 // soften it, especially at low T
846 // double arg = lr * sqrt(beta * ( zn*zn*tks.dz2_o[i] + tn*tn*tks.dt2_o[i] ) );
847 // if(abs(arg) < 20)
848 // {
849 // double t = exp(-arg);
850 // tl = t/(t+1.);
851 // tr = 1/(t+1.);
852 // }
853
854 double p = vtx.pk[k] * tks.w[i];
855 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];
856 double wt = p*tks.dt2_o[i];
857 double wz = p*tks.dz2_o[i];
858 p1 += p*tl; z1 += wz*tl*tks.z[i]; t1 += wt*tl*tks.t[i]; wz1 += wz*tl; wt1 += wt*tl;
859 p2 += p*tr; z2 += wz*tr*tks.z[i]; t2 += wt*tr*tks.t[i]; wz2 += wz*tr; wt2 += wt*tr;
860 }
861 }
862
863 if(wz1 > 0 && wt1 > 0 && wz2 > 0 && wt2 > 0)
864 {
865 t1 /= wt1;
866 z1 /= wz1;
867 t2 /= wt2;
868 z2 /= wz2;
869
870 if( fVerbose > 3 )
871 {
872 double aux = (z1-z2)*(z1-z2)/(fVertexZSize*fVertexZSize) + (t1-t2)*(t1-t2)/(fVertexTSize*fVertexTSize);
873 cout << "weighted split: delta = " << sqrt(aux) << endl;
874 }
875 }
876 else
877 {
878 continue;
879 // plot_split_crush(zn, tn, vtx, tks, k);
880 // throw std::invalid_argument( "0 division" );
881 }
[29d662e]882
[84a097e]883 while(vtx.NearestCluster(t1, z1) != k || vtx.NearestCluster(t2, z2) != k)
884 {
885 t1 = 0.5 * (t1 + t_old);
886 z1 = 0.5 * (z1 + z_old);
887 t2 = 0.5 * (t2 + t_old);
888 z2 = 0.5 * (z2 + z_old);
889 }
[29d662e]890
[84a097e]891 // Compute final distance and split if the distance is enough
892 double delta2 = (z1-z2)*(z1-z2)/(fVertexZSize*fVertexZSize) + (t1-t2)*(t1-t2)/(fVertexTSize*fVertexTSize);
893 if(delta2 > fD2Merge)
894 {
895 split = true;
896 vtx.t[k] = t1;
897 vtx.z[k] = z1;
898 vtx.pk[k] = p1 * pk_old/(p1+p2);
899
900 double new_t = t2;
901 double new_z = z2;
902 double new_pk = p2 * pk_old/(p1+p2);
903
904 vtx.addItem(new_z, new_t, new_pk);
905
906 if( fVerbose > 3 )
907 {
908 cout << "===Split happened on vtx " << k << endl;
909 cout << "OLD z: " << z_old << " , t: " << t_old << " , pk: " << pk_old << endl;
910 cout << "NEW+ z: " << vtx.z[k] << " , t: " << vtx.t[k] << " , pk: " << vtx.pk[k] << endl;
911 cout << "NEW- z: " << new_z << " , t: " << new_t << " , pk: " << new_pk << endl;
912 }
913 }
914 }
[29d662e]915 }
916 }
[84a097e]917 return split;
[29d662e]918}
919
920
[84a097e]921//------------------------------------------------------------------------------
922// Merge vertexes closer than declared dimensions
923bool VertexFinderDA4D::merge(vertex_t & vtx, double d2_merge = 2)
[4154bbd]924{
[84a097e]925 bool merged = false;
[29d662e]926
[84a097e]927 if(vtx.getSize() < 2) return merged;
[341014c]928
[84a097e]929 bool last_merge = false;
930 do {
931 double min_d2 = d2_merge;
932 unsigned int k1_min, k2_min;
933 for(unsigned int k1 = 0; k1 < vtx.getSize(); k1++)
934 {
935 for(unsigned int k2 = k1+1; k2 < vtx.getSize();k2++)
[341014c]936 {
[84a097e]937 double d2_tmp = vtx.DistanceSquare(k1, k2);
938 if(d2_tmp < min_d2)
939 {
940 min_d2 = d2_tmp;
941 k1_min = k1;
942 k2_min = k2;
943 }
[29d662e]944 }
[84a097e]945 }
[29d662e]946
[84a097e]947 if(min_d2 < d2_merge)
948 {
949 vtx.mergeItems(k1_min, k2_min);
950 last_merge = true;
951 merged = true;
[29d662e]952 }
[84a097e]953 else last_merge = false;
954 } while(last_merge);
[29d662e]955
[84a097e]956 return merged;
[29d662e]957}
958
[84a097e]959// -----------------------------------------------------------------------------
960// Compute all the energies and set the partition function normalization for each track
961vector<double> VertexFinderDA4D::Compute_pk_exp_mBetaE(double beta, vertex_t &vtx, tracks_t &tks, double Z_init)
[4154bbd]962{
[84a097e]963 unsigned int nt = tks.getSize();
964 unsigned int nv = vtx.getSize();
[341014c]965
[84a097e]966 vector<double> pk_exp_mBetaE(nt * nv);
967 for (unsigned int k = 0; k < nv; k++)
[341014c]968 {
[84a097e]969 for (unsigned int i = 0; i < nt; i++)
[341014c]970 {
[84a097e]971 if(k == 0) tks.Z[i] = Z_init;
[341014c]972
[84a097e]973 double aux = Energy(tks.z[i], vtx.z[k], tks.dz2_o[i], tks.t[i], vtx.t[k], tks.dt2_o[i]);
974 aux = vtx.pk[k] * exp(-beta * aux);
975 // if(aux < 1E-10) continue;
976 tks.Z[i] += aux;
977
978 unsigned int idx = k*nt + i;
979 pk_exp_mBetaE[idx] = aux;
[29d662e]980 }
981 }
[84a097e]982 return pk_exp_mBetaE;
[29d662e]983}
984
985//------------------------------------------------------------------------------
[84a097e]986// Eliminate clusters with only one significant/unique track
987bool VertexFinderDA4D::purge(vertex_t & vtx, tracks_t & tks, double & rho0, const double beta, double min_prob, double min_trk)
[4154bbd]988{
[84a097e]989 const unsigned int nv = vtx.getSize();
990 const unsigned int nt = tks.getSize();
991
992 if (nv < 2)
993 return false;
[341014c]994
995 double sumpmin = nt;
[84a097e]996 unsigned int k0 = nv;
997
998 int nUnique = 0;
999 double sump = 0;
1000
1001 double Z_init = rho0 * exp(-beta * fMuOutlayer * fMuOutlayer); // Add fDtCutOff here toghether with this
1002 vector<double> pk_exp_mBetaE = Compute_pk_exp_mBetaE(beta, vtx, tks, Z_init);
1003
1004 for (unsigned int k = 0; k < nv; ++k) {
1005
1006 nUnique = 0;
1007 sump = 0;
1008
1009 double pmax = vtx.pk[k] / (vtx.pk[k] + rho0 * exp(-beta * fMuOutlayer* fMuOutlayer));
1010 double pcut = min_prob * pmax;
1011
1012 for (unsigned int i = 0; i < nt; ++i) {
1013 unsigned int idx = k*nt + i;
1014
1015 if(pk_exp_mBetaE[idx] == 0 || tks.Z[i] == 0)
[341014c]1016 {
[84a097e]1017 continue;
[29d662e]1018 }
[84a097e]1019
1020 double p = pk_exp_mBetaE[idx] / tks.Z[i];
1021 sump += p;
1022 if( ( p > pcut ) & ( tks.w[i] > 0 ) ) nUnique++;
[29d662e]1023 }
1024
[84a097e]1025 if ((nUnique < min_trk) && (sump < sumpmin)) {
[341014c]1026 sumpmin = sump;
1027 k0 = k;
[29d662e]1028 }
[84a097e]1029
[29d662e]1030 }
1031
[84a097e]1032 if (k0 != nv) {
1033 if (fVerbose > 5) {
1034 std::cout << Form("eliminating prototype at z = %.3f mm, t = %.0f ps", vtx.z[k0], vtx.t[k0]) << " with sump=" << sumpmin
1035 << " rho*nt =" << vtx.pk[k0]*nt
1036 << endl;
1037 }
1038 vtx.removeItem(k0);
[29d662e]1039 return true;
[84a097e]1040 } else {
[29d662e]1041 return false;
1042 }
1043}
1044
1045
[84a097e]1046// -----------------------------------------------------------------------------
1047// Plot status
1048void VertexFinderDA4D::plot_status(double beta, vertex_t &vtx, tracks_t &tks, int n_it, const char* flag)
[4154bbd]1049{
[84a097e]1050 vector<int> vtx_color = {2,4,8,1,5,6,9,14,46,3};
1051 while(vtx.getSize() > vtx_color.size()) vtx_color.push_back(40);
[29d662e]1052
[84a097e]1053 vector<double> t_PV, dt_PV, z_PV, dz_PV;
1054 vector<double> t_PU, dt_PU, z_PU, dz_PU;
[29d662e]1055
[84a097e]1056 double ETot = 0;
1057 vector<double> pk_exp_mBetaE = Compute_pk_exp_mBetaE(beta, vtx, tks, 0);
[29d662e]1058
[84a097e]1059 for(unsigned int i = 0; i < tks.getSize(); i++)
1060 {
1061 for(unsigned int k = 0; k < vtx.getSize(); k++)
[341014c]1062 {
[84a097e]1063 unsigned int idx = k*tks.getSize() + i;
1064 if(pk_exp_mBetaE[idx] == 0) continue;
1065
1066 double p_ygx = pk_exp_mBetaE[idx] / tks.Z[i];
1067
1068 ETot += tks.w[i] * p_ygx * Energy(tks.z[i], vtx.z[k], tks.dz2_o[i], tks.t[i], vtx.t[k], tks.dt2_o[i]);
[29d662e]1069 }
1070
[84a097e]1071 if(tks.tt[i]->IsPU)
[341014c]1072 {
[84a097e]1073 t_PU.push_back(tks.t[i]);
1074 dt_PU.push_back(1./tks.dt_o[i]);
1075 z_PU.push_back(tks.z[i]);
1076 dz_PU.push_back(1./tks.dz_o[i]);
1077 }
1078 else
1079 {
1080 t_PV.push_back(tks.t[i]);
1081 dt_PV.push_back(1./tks.dt_o[i]);
1082 z_PV.push_back(tks.z[i]);
1083 dz_PV.push_back(1./tks.dz_o[i]);
[29d662e]1084 }
[341014c]1085 }
[84a097e]1086
1087
1088 ETot /= tks.sum_w;
1089 fEnergy_rec.push_back(ETot);
1090 fBeta_rec.push_back(beta);
1091 fNvtx_rec.push_back(vtx.getSize());
1092
1093 double t_min = TMath::Min( TMath::MinElement(t_PV.size(), &t_PV[0]), TMath::MinElement(t_PU.size(), &t_PU[0]) );
1094 t_min = TMath::Min(t_min, TMath::MinElement(vtx.getSize(), &(vtx.t[0])) ) - fVertexTSize;
1095 double t_max = TMath::Max( TMath::MaxElement(t_PV.size(), &t_PV[0]), TMath::MaxElement(t_PU.size(), &t_PU[0]) );
1096 t_max = TMath::Max(t_max, TMath::MaxElement(vtx.getSize(), &(vtx.t[0])) ) + fVertexTSize;
1097
1098 double z_min = TMath::Min( TMath::MinElement(z_PV.size(), &z_PV[0]), TMath::MinElement(z_PU.size(), &z_PU[0]) );
1099 z_min = TMath::Min(z_min, TMath::MinElement(vtx.getSize(), &(vtx.z[0])) ) - 5;
1100 double z_max = TMath::Max( TMath::MaxElement(z_PV.size(), &z_PV[0]), TMath::MaxElement(z_PU.size(), &z_PU[0]) );
1101 z_max = TMath::Max(z_max, TMath::MaxElement(vtx.getSize(), &(vtx.z[0])) ) + 5;
1102
1103 auto c_2Dspace = new TCanvas("c_2Dspace", "c_2Dspace", 800, 600);
1104
1105 TGraphErrors* gr_PVtks = new TGraphErrors(t_PV.size(), &t_PV[0], &z_PV[0], &dt_PV[0], &dz_PV[0]);
1106 gr_PVtks->SetTitle(Form("Clustering space - #beta = %.6f", beta));
1107 gr_PVtks->GetXaxis()->SetTitle("t CA [ps]");
1108 gr_PVtks->GetXaxis()->SetLimits(t_min, t_max);
1109 gr_PVtks->GetYaxis()->SetTitle("z CA [mm]");
1110 gr_PVtks->GetYaxis()->SetRangeUser(z_min, z_max);
1111 gr_PVtks->SetMarkerStyle(4);
1112 gr_PVtks->SetMarkerColor(8);
1113 gr_PVtks->SetLineColor(8);
1114 gr_PVtks->Draw("APE1");
1115
1116 TGraphErrors* gr_PUtks = new TGraphErrors(t_PU.size(), &t_PU[0], &z_PU[0], &dt_PU[0], &dz_PU[0]);
1117 gr_PUtks->SetMarkerStyle(3);
1118 gr_PUtks->Draw("PE1");
1119
1120 TGraph* gr_vtx = new TGraph(vtx.getSize(), &(vtx.t[0]), &(vtx.z[0]));
1121 gr_vtx->SetMarkerStyle(28);
1122 gr_vtx->SetMarkerColor(2);
1123 gr_vtx->SetMarkerSize(2.);
1124 gr_vtx->Draw("PE1");
1125
1126 fItInputGenVtx->Reset();
1127 TGraph* gr_genvtx = new TGraph(fInputGenVtx->GetEntriesFast());
1128 Candidate *candidate;
1129 unsigned int k = 0;
1130 while((candidate = static_cast<Candidate*>(fItInputGenVtx->Next())))
[341014c]1131 {
[84a097e]1132 gr_genvtx->SetPoint(k, candidate->Position.T()*1E9/c_light, candidate->Position.Z());
1133 k++;
[29d662e]1134 }
[84a097e]1135 gr_genvtx->SetMarkerStyle(33);
1136 gr_genvtx->SetMarkerColor(6);
1137 gr_genvtx->SetMarkerSize(2.);
1138 gr_genvtx->Draw("PE1");
[29d662e]1139
[84a097e]1140 // auto leg = new TLegend(0.1, 0.1);
1141 // leg->AddEntry(gr_PVtks, "PV tks", "ep");
1142 // leg->AddEntry(gr_PUtks, "PU tks", "ep");
1143 // leg->AddEntry(gr_vtx, "Cluster center", "p");
1144 // leg->Draw();
1145
1146 c_2Dspace->SetGrid();
1147 c_2Dspace->SaveAs(fFigFolderPath + Form("/c_2Dspace_beta%010.0f-%s%d.png", 1E7*beta, flag, n_it));
1148
1149 delete c_2Dspace;
1150}
[29d662e]1151
[84a097e]1152// -----------------------------------------------------------------------------
1153// Plot status at the end
1154void VertexFinderDA4D::plot_status_end(vertex_t &vtx, tracks_t &tks)
[4154bbd]1155{
[84a097e]1156 unsigned int nv = vtx.getSize();
[29d662e]1157
[84a097e]1158 // Define colors in a meaningfull way
1159 vector<int> MyPalette(nv);
1160
1161 const int Number = 3;
1162 double Red[Number] = { 1.00, 0.00, 0.00};
1163 double Green[Number] = { 0.00, 1.00, 0.00};
1164 double Blue[Number] = { 1.00, 0.00, 1.00};
1165 double Length[Number] = { 0.00, 0.50, 1.00 };
1166 int FI = TColor::CreateGradientColorTable(Number,Length,Red,Green,Blue,nv);
1167 for (unsigned int i=0;i<nv;i++) MyPalette[i] = FI+i;
1168
1169 TCanvas * c_out = new TCanvas("c_out", "c_out", 800, 600);
1170 double t_min = TMath::Min( TMath::MinElement(tks.getSize(), &tks.t[0]), TMath::MinElement(vtx.getSize(), &(vtx.t[0])) ) - 2*fVertexTSize;
1171 double t_max = TMath::Max(TMath::MaxElement(tks.getSize(), &tks.t[0]), TMath::MaxElement(vtx.getSize(), &(vtx.t[0])) ) + 2*fVertexTSize;
[29d662e]1172
[84a097e]1173 double z_min = TMath::Min( TMath::MinElement(tks.getSize(), &tks.z[0]), TMath::MinElement(vtx.getSize(), &(vtx.z[0])) ) - 15;
1174 double z_max = TMath::Max( TMath::MaxElement(tks.getSize(), &tks.z[0]), TMath::MaxElement(vtx.getSize(), &(vtx.z[0])) ) + 15;
[29d662e]1175
[84a097e]1176 // Draw tracks
1177 for(unsigned int i = 0; i < tks.getSize(); i++)
[341014c]1178 {
[84a097e]1179 double dt[] = {1./tks.dt_o[i]};
1180 double dz[] = {1./tks.dz_o[i]};
1181 TGraphErrors* gr = new TGraphErrors(1, &(tks.t[i]), &(tks.z[i]), dt, dz);
1182
1183 gr->SetNameTitle(Form("gr%d",i), Form("gr%d",i));
1184
1185 int marker = tks.tt[i]->IsPU? 1 : 4;
1186 gr->SetMarkerStyle(marker);
1187
1188 int idx = tks.tt[i]->ClusterIndex;
1189 int color = idx>=0 ? MyPalette[idx] : 13;
1190 gr->SetMarkerColor(color);
1191 gr->SetLineColor(color);
1192
1193 int line_style = idx>=0 ? 1 : 3;
1194 gr->SetLineStyle(line_style);
1195
1196 if(i==0)
[341014c]1197 {
[84a097e]1198 gr->SetTitle(Form("Clustering space - Tot Vertexes = %d", nv));
1199 gr->GetXaxis()->SetTitle("t CA [ps]");
1200 gr->GetXaxis()->SetLimits(t_min, t_max);
1201 gr->GetYaxis()->SetTitle("z CA [mm]");
1202 gr->GetYaxis()->SetRangeUser(z_min, z_max);
1203 gr->Draw("APE1");
[29d662e]1204 }
[84a097e]1205 else gr->Draw("PE1");
[29d662e]1206 }
1207
[84a097e]1208 // Draw vertices
1209 for(unsigned int k = 0; k < vtx.getSize(); k++)
[341014c]1210 {
[84a097e]1211 TGraph* gr = new TGraph(1, &(vtx.t[k]), &(vtx.z[k]));
[29d662e]1212
[84a097e]1213 gr->SetNameTitle(Form("grv%d",k), Form("grv%d",k));
[29d662e]1214
[84a097e]1215 gr->SetMarkerStyle(41);
1216 gr->SetMarkerSize(2.);
1217 gr->SetMarkerColor(MyPalette[k]);
1218
1219 gr->Draw("P");
[29d662e]1220 }
1221
[84a097e]1222 fItInputGenVtx->Reset();
1223 TGraph* gr_genvtx = new TGraph(fInputGenVtx->GetEntriesFast());
1224 TGraph* gr_genPV = new TGraph(1);
1225 Candidate *candidate;
1226 unsigned int k = 0;
1227 while((candidate = static_cast<Candidate*>(fItInputGenVtx->Next())))
1228 {
1229 if(k == 0 ) {
1230 gr_genPV->SetPoint(k, candidate->Position.T()*1E9/c_light, candidate->Position.Z());
1231 }
1232 else gr_genvtx->SetPoint(k, candidate->Position.T()*1E9/c_light, candidate->Position.Z());
[29d662e]1233
[84a097e]1234 k++;
1235 }
1236 gr_genvtx->SetMarkerStyle(20);
1237 gr_genvtx->SetMarkerColorAlpha(kBlack, 0.8);
1238 gr_genvtx->SetMarkerSize(.8);
1239 gr_genvtx->Draw("PE1");
1240 gr_genPV->SetMarkerStyle(33);
1241 gr_genPV->SetMarkerColorAlpha(kBlack, 1);
1242 gr_genPV->SetMarkerSize(2.5);
1243 gr_genPV->Draw("PE1");
1244
1245 // auto note = new TLatex();
1246 // note->DrawLatexNDC(0.5, 0.8, Form("#splitline{Vertexes Reco = %d }{Vertexes gen = %d}", vtx.getSize(), k) );
1247
1248 c_out->SetGrid();
1249 c_out->SaveAs(fFigFolderPath + Form("/c_final.root"));
1250 delete c_out;
1251}
[29d662e]1252
[84a097e]1253// -----------------------------------------------------------------------------
1254// Plot splitting
1255void VertexFinderDA4D::plot_split_crush(double zn, double tn, vertex_t &vtx, tracks_t &tks, int i_vtx)
[4154bbd]1256{
[84a097e]1257 vector<double> t, dt, z, dz;
[29d662e]1258
[84a097e]1259 for(unsigned int i = 0; i < tks.getSize(); i++)
1260 {
1261 t.push_back(tks.t[i]);
1262 dt.push_back(1./tks.dt_o[i]);
1263 z.push_back(tks.z[i]);
1264 dz.push_back(1./tks.dz_o[i]);
1265 }
[29d662e]1266
1267
[84a097e]1268 double t_min = TMath::Min(TMath::MinElement(t.size(), &t[0]), TMath::MinElement(vtx.getSize(), &(vtx.t[0])) ) - 50;
1269 double t_max = TMath::Max(TMath::MaxElement(t.size(), &t[0]), TMath::MaxElement(vtx.getSize(), &(vtx.t[0])) ) + 50;
1270
1271 double z_min = TMath::Min(TMath::MinElement(z.size(), &z[0]), TMath::MinElement(vtx.getSize(), &(vtx.z[0])) ) - 5;
1272 double z_max = TMath::Max(TMath::MaxElement(z.size(), &z[0]), TMath::MaxElement(vtx.getSize(), &(vtx.z[0])) ) + 5;
1273
1274 auto c_2Dspace = new TCanvas("c_2Dspace", "c_2Dspace", 800, 600);
1275
1276 TGraphErrors* gr_PVtks = new TGraphErrors(t.size(), &t[0], &z[0], &dt[0], &dz[0]);
1277 gr_PVtks->SetTitle(Form("Clustering space"));
1278 gr_PVtks->GetXaxis()->SetTitle("t CA [ps]");
1279 gr_PVtks->GetXaxis()->SetLimits(t_min, t_max);
1280 gr_PVtks->GetYaxis()->SetTitle("z CA [mm]");
1281 gr_PVtks->GetYaxis()->SetRangeUser(z_min, z_max);
1282 gr_PVtks->SetMarkerStyle(4);
1283 gr_PVtks->SetMarkerColor(1);
1284 gr_PVtks->SetLineColor(1);
1285 gr_PVtks->Draw("APE1");
1286
1287 TGraph* gr_vtx = new TGraph(1, &(vtx.t[i_vtx]), &(vtx.z[i_vtx]));
1288 gr_vtx->SetMarkerStyle(28);
1289 gr_vtx->SetMarkerColor(2);
1290 gr_vtx->SetMarkerSize(2.);
1291 gr_vtx->Draw("PE1");
1292
1293 double t_pos[] = {vtx.t[i_vtx], vtx.t[i_vtx]+100};
1294 double t_neg[] = {vtx.t[i_vtx], vtx.t[i_vtx]-100};
1295 double z_pos[] = {vtx.z[i_vtx], vtx.z[i_vtx]+(zn/tn)*100};
1296 double z_neg[] = {vtx.z[i_vtx], vtx.z[i_vtx]-(zn/tn)*100};
1297
1298 TGraph* gr_pos = new TGraph(2, &t_pos[0], &z_pos[0]);
1299 gr_pos->SetLineColor(8);
1300 gr_pos->SetMarkerColor(8);
1301 gr_pos->Draw("PL");
1302 TGraph* gr_neg = new TGraph(2, &t_neg[0], &z_neg[0]);
1303 gr_neg->SetLineColor(4);
1304 gr_neg->SetMarkerColor(4);
1305 gr_neg->Draw("PL");
1306
1307
1308 c_2Dspace->SetGrid();
1309 c_2Dspace->SaveAs(fFigFolderPath + Form("/crush_splitting.png"));
[29d662e]1310
[84a097e]1311 delete c_2Dspace;
[29d662e]1312}
Note: See TracBrowser for help on using the repository browser.