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