Fork me on GitHub

source: git/modules/VertexFinderDA4D.cc@ e57c062

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

add print out for the vertex position error

  • Property mode set to 100644
File size: 31.9 KB
Line 
1/** \class VertexFinderDA4D
2 *
3 * Cluster vertices from tracks using deterministic annealing and timing information
4 *
5 * \authors M. Selvaggi, L. Gray
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 "TVector3.h"
29
30#include <utility>
31#include <algorithm>
32#include <stdexcept>
33#include <iostream>
34#include <vector>
35
36using namespace std;
37
38static const Double_t mm = 1.;
39static const Double_t m = 1000.*mm;
40static const Double_t ns = 1.;
41static const Double_t s = 1.e+9 *ns;
42static const Double_t c_light = 2.99792458e+8 * m/s;
43
44struct track_t
45{
46 double z; // z-coordinate at point of closest approach to the beamline
47 double t; // t-coordinate at point of closest approach to the beamline
48 double dz2; // square of the error of z(pca)
49 double dtz; // covariance of z-t
50 double dt2; // square of the error of t(pca)
51 Candidate *tt; // a pointer to the Candidate Track
52 double Z; // Z[i] for DA clustering
53 double pi; // track weight
54 double pt;
55 double eta;
56 double phi;
57};
58
59struct vertex_t
60{
61 double z;
62 double t;
63 double pk; // vertex weight for "constrained" clustering
64 // --- temporary numbers, used during update
65 double ei;
66 double sw;
67 double swz;
68 double swt;
69 double se;
70 // ---for Tc
71 double swE;
72 double Tc;
73};
74
75static bool split(double beta, std::vector<track_t> &tks, std::vector<vertex_t> &y);
76static double update1(double beta, std::vector<track_t> &tks, std::vector<vertex_t> &y);
77static double update2(double beta, std::vector<track_t> &tks, std::vector<vertex_t> &y, double &rho0, const double dzCutOff);
78static void dump(const double beta, const std::vector<vertex_t> & y, const std::vector<track_t> & tks);
79static bool merge(std::vector<vertex_t> &);
80static bool merge(std::vector<vertex_t> &, double &);
81static bool purge(std::vector<vertex_t> &, std::vector<track_t> & , double &, const double, const double);
82static void splitAll(std::vector<vertex_t> &y);
83static double beta0(const double betamax, std::vector<track_t> &tks, std::vector<vertex_t> &y, const double coolingFactor);
84static double Eik(const track_t &t, const vertex_t &k);
85
86static bool recTrackLessZ1(const track_t & tk1, const track_t & tk2)
87{
88 return tk1.z < tk2.z;
89}
90
91using namespace std;
92
93//------------------------------------------------------------------------------
94
95VertexFinderDA4D::VertexFinderDA4D() :
96 fVerbose(0), fMinPT(0), fVertexSpaceSize(0), fVertexTimeSize(0),
97 fUseTc(0), fBetaMax(0), fBetaStop(0), fCoolingFactor(0),
98 fMaxIterations(0), fDzCutOff(0), fD0CutOff(0), fDtCutOff(0)
99{
100}
101
102//------------------------------------------------------------------------------
103
104VertexFinderDA4D::~VertexFinderDA4D()
105{
106}
107
108//------------------------------------------------------------------------------
109
110void VertexFinderDA4D::Init()
111{
112
113 fVerbose = GetBool("Verbose", 1);
114 fMinPT = GetDouble("MinPT", 0.1);
115 fVertexSpaceSize = GetDouble("VertexSpaceSize", 0.5); //in mm
116 fVertexTimeSize = GetDouble("VertexTimeSize", 10E-12); //in s
117 fUseTc = GetBool("UseTc", 1);
118 fBetaMax = GetDouble("BetaMax ", 0.1);
119 fBetaStop = GetDouble("BetaStop", 1.0);
120 fCoolingFactor = GetDouble("CoolingFactor", 0.8);
121 fMaxIterations = GetInt("MaxIterations", 100);
122 fDzCutOff = GetDouble("DzCutOff", 40); // Adaptive Fitter uses 30 mm but that appears to be a bit tight here sometimes
123 fD0CutOff = GetDouble("D0CutOff", 30);
124 fDtCutOff = GetDouble("DtCutOff", 100E-12); // dummy
125
126 // convert stuff in cm, ns
127 fVertexSpaceSize /= 10.0;
128 fVertexTimeSize *= 1E9;
129 fDzCutOff /= 10.0; // Adaptive Fitter uses 3.0 but that appears to be a bit tight here sometimes
130 fD0CutOff /= 10.0;
131
132 fInputArray = ImportArray(GetString("InputArray", "TrackSmearing/tracks"));
133 fItInputArray = fInputArray->MakeIterator();
134
135 fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
136 fVertexOutputArray = ExportArray(GetString("VertexOutputArray", "vertices"));
137}
138
139//------------------------------------------------------------------------------
140
141void VertexFinderDA4D::Finish()
142{
143 if(fItInputArray) delete fItInputArray;
144}
145
146//------------------------------------------------------------------------------
147
148void VertexFinderDA4D::Process()
149{
150 Candidate *candidate, *track;
151 TObjArray *ClusterArray;
152 ClusterArray = new TObjArray;
153 TIterator *ItClusterArray;
154 Int_t ivtx = 0;
155
156 fInputArray->Sort();
157
158 TLorentzVector pos, mom;
159 if (fVerbose)
160 {
161 cout<<" start processing vertices ..."<<endl;
162 cout<<" Found "<<fInputArray->GetEntriesFast()<<" input tracks"<<endl;
163 //loop over input tracks
164 fItInputArray->Reset();
165 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
166 {
167 pos = candidate->InitialPosition;
168 mom = candidate->Momentum;
169
170 cout<<"pt: "<<mom.Pt()<<", eta: "<<mom.Eta()<<", phi: "<<mom.Phi()<<", z: "<<candidate->DZ/10<<endl;
171 }
172 }
173
174 // clusterize tracks in Z
175 clusterize(*fInputArray, *ClusterArray);
176
177 if (fVerbose){std::cout << " clustering returned "<< ClusterArray->GetEntriesFast() << " clusters from " << fInputArray->GetEntriesFast() << " selected tracks" <<std::endl;}
178
179 //loop over vertex candidates
180 ItClusterArray = ClusterArray->MakeIterator();
181 ItClusterArray->Reset();
182 while((candidate = static_cast<Candidate*>(ItClusterArray->Next())))
183 {
184
185 double meantime = 0.;
186 double expv_x2 = 0.;
187 double normw = 0.;
188 double errtime = 0;
189
190 double meanpos = 0.;
191 double meanerr2 = 0.;
192 double normpos = 0.;
193 double errpos = 0.;
194
195 double sumpt2 = 0.;
196
197 int itr = 0;
198
199 if(fVerbose)cout<<"this vertex has: "<<candidate->GetCandidates()->GetEntriesFast()<<" tracks"<<endl;
200
201 // loop over tracks belonging to this vertex
202 TIter it1(candidate->GetCandidates());
203 it1.Reset();
204
205 while((track = static_cast<Candidate*>(it1.Next())))
206 {
207
208 itr++;
209 // TBC: the time is in ns for now TBC
210 double t = track->InitialPosition.T()/c_light;
211 double dt = track->ErrorT/c_light;
212 const double time = t;
213 const double inverr = 1.0/dt;
214 meantime += time*inverr;
215 expv_x2 += time*time*inverr;
216 normw += inverr;
217
218 // compute error position TBC
219 const double pt = track->Momentum.Pt();
220 const double z = track->DZ/10.0;
221 const double err_pt = track->ErrorPT;
222 const double err_z = track->ErrorDZ;
223
224 const double wi = (pt/(err_pt*err_z))*(pt/(err_pt*err_z));
225 meanpos += z*wi;
226
227 meanerr2 += err_z*err_z*wi;
228 normpos += wi;
229 sumpt2 += pt*pt;
230
231 // while we are here store cluster index in tracks
232 track->ClusterIndex = ivtx;
233 }
234
235 meantime = meantime/normw;
236 expv_x2 = expv_x2/normw;
237 errtime = TMath::Sqrt((expv_x2 - meantime*meantime)/itr);
238 meanpos = meanpos/normpos;
239 meanerr2 = meanerr2/normpos;
240 errpos = TMath::Sqrt(meanerr2/itr);
241
242 candidate->Position.SetXYZT(0.0, 0.0, meanpos*10.0 , meantime*c_light);
243 candidate->PositionError.SetXYZT(0.0, 0.0, errpos*10.0 , errtime*c_light);
244 candidate->SumPT2 = sumpt2;
245 candidate->ClusterNDF = itr;
246 candidate->ClusterIndex = ivtx;
247
248 fVertexOutputArray->Add(candidate);
249
250 ivtx++;
251
252 if (fVerbose){
253 std::cout << "x,y,z";
254 std::cout << ",t";
255 std::cout << "=" << candidate->Position.X()/10.0 <<" " << candidate->Position.Y()/10.0 << " " << candidate->Position.Z()/10.0;
256 std::cout << " " << candidate->Position.T()/c_light;
257
258 std::cout << std::endl;
259 std::cout << "sumpt2 " << candidate->SumPT2<<endl;
260
261 std::cout << "ex,ey,ez";
262 std::cout << ",et";
263 std::cout << "=" << candidate->PositionError.X()/10.0 <<" " << candidate->PositionError.Y()/10.0 << " " << candidate->PositionError.Z()/10.0;
264 std::cout << " " << candidate->PositionError.T()/c_light;
265 std::cout << std::endl;
266
267 }
268 }// end of cluster loop
269
270
271 if(fVerbose){
272 std::cout << "PrimaryVertexProducerAlgorithm::vertices candidates =" << ClusterArray->GetEntriesFast() << std::endl;
273 }
274
275 //TBC maybe this can be done later
276 // sort vertices by pt**2 vertex (aka signal vertex tagging)
277 /*if(pvs.size()>1){
278 sort(pvs.begin(), pvs.end(), VertexHigherPtSquared());
279 }
280 */
281
282 delete ClusterArray;
283
284}
285
286//------------------------------------------------------------------------------
287
288void VertexFinderDA4D::clusterize(const TObjArray &tracks, TObjArray &clusters)
289{
290 if(fVerbose) {
291 cout << "###################################################" << endl;
292 cout << "# VertexFinderDA4D::clusterize nt="<<tracks.GetEntriesFast() << endl;
293 cout << "###################################################" << endl;
294 }
295
296 vector< Candidate* > pv = vertices();
297
298 if(fVerbose){ cout << "# VertexFinderDA4D::clusterize pv.size="<<pv.size() << endl; }
299 if (pv.size()==0){ return; }
300
301 // convert into vector of candidates
302 //TObjArray *ClusterArray = pv.begin()->GetCandidates();
303 //Candidate *aCluster = static_cast<Candidate*>(&(pv.at(0)));
304 Candidate *aCluster = pv.at(0);
305
306 // fill into clusters and merge
307
308
309 if( fVerbose ) {
310 std::cout << '\t' << 0;
311 std::cout << ' ' << (*pv.begin())->Position.Z()/10.0 << ' ' << (*pv.begin())->Position.T()/c_light << std::endl;
312 }
313
314 for(vector<Candidate*>::iterator k=pv.begin()+1; k!=pv.end(); k++){
315 if( fVerbose ) {
316 std::cout << '\t' << std::distance(pv.begin(),k);
317 std::cout << ' ' << (*k)->Position.Z() << ' ' << (*k)->Position.T() << std::endl;
318 }
319
320
321 // TBC - check units here
322 if ( std::abs((*k)->Position.Z() - (*(k-1))->Position.Z())/10.0 > (2*fVertexSpaceSize) ||
323 std::abs((*k)->Position.T() - (*(k-1))->Position.Z())/c_light > 2*0.010 ) {
324 // close a cluster
325 clusters.Add(aCluster);
326 //aCluster.clear();
327 }
328 //for(unsigned int i=0; i<k->GetCandidates().GetEntriesFast(); i++){
329 aCluster = *k;
330 //}
331
332 }
333 clusters.Add(aCluster);
334
335 if(fVerbose) { std::cout << "# VertexFinderDA4D::clusterize clusters.size="<<clusters.GetEntriesFast() << std::endl; }
336
337}
338
339//------------------------------------------------------------------------------
340
341vector< Candidate* > VertexFinderDA4D::vertices()
342{
343 Candidate *candidate;
344 UInt_t clusterIndex = 0;
345 vector< Candidate* > clusters;
346
347 vector<track_t> tks;
348 track_t tr;
349 Double_t z, dz, t, l, dt, d0, d0error;
350
351 // loop over input tracks
352 fItInputArray->Reset();
353 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
354 {
355 //TBC everything in cm
356 z = candidate->DZ/10;
357 tr.z = z;
358 dz = candidate->ErrorDZ/10;
359 tr.dz2 = dz*dz // track error
360 //TBC: beamspot size induced error, take 0 for now.
361 // + (std::pow(beamspot.BeamWidthX()*cos(phi),2.)+std::pow(beamspot.BeamWidthY()*sin(phi),2.))/std::pow(tantheta,2.) // beam-width induced
362 + fVertexSpaceSize*fVertexSpaceSize; // intrinsic vertex size, safer for outliers and short lived decays
363
364 // TBC: the time is in ns for now TBC
365 //t = candidate->Position.T()/c_light;
366 t = candidate->InitialPosition.T()/c_light;
367 l = candidate->L/c_light;
368 double pt = candidate->Momentum.Pt();
369 double eta = candidate->Momentum.Eta();
370 double phi = candidate->Momentum.Phi();
371
372 tr.pt = pt;
373 tr.eta = eta;
374 tr.phi = phi;
375 tr.t = t; //
376 tr.dtz = 0.;
377 dt = candidate->ErrorT/c_light;
378 tr.dt2 = dt*dt + fVertexTimeSize*fVertexTimeSize; // the ~injected~ timing error plus a small minimum vertex size in time
379 if(fD0CutOff>0)
380 {
381
382 d0 = TMath::Abs(candidate->D0)/10.0;
383 d0error = candidate->ErrorD0/10.0;
384
385 tr.pi=1./(1.+exp((d0*d0)/(d0error*d0error) - fD0CutOff*fD0CutOff)); // reduce weight for high ip tracks
386
387 }
388 else
389 {
390 tr.pi=1.;
391 }
392 tr.tt=&(*candidate);
393 tr.Z=1.;
394
395 // TBC now putting track selection here (> fPTMin)
396 if(tr.pi > 1e-3 && tr.pt > fMinPT)
397 {
398 tks.push_back(tr);
399 }
400 }
401
402 //print out input tracks
403
404 if(fVerbose)
405 {
406 std::cout<<" start processing vertices ..."<<std::endl;
407 std::cout<<" Found "<<tks.size()<<" input tracks"<<std::endl;
408 //loop over input tracks
409
410
411 for(std::vector<track_t>::const_iterator it=tks.begin(); it!=tks.end(); it++){
412 double z = it->z;
413 double pt=it->pt;
414 double eta=it->eta;
415 double phi=it->phi;
416 double t = it->t;
417
418 std::cout<<"pt: "<<pt<<", eta: "<<eta<<", phi: "<<phi<<", z: "<<z<<", t: "<<t<<std::endl;
419 }
420 }
421
422 unsigned int nt=tks.size();
423 double rho0=0.0; // start with no outlier rejection
424
425 if (tks.empty()) return clusters;
426
427 vector<vertex_t> y; // the vertex prototypes
428
429 // initialize:single vertex at infinite temperature
430 vertex_t vstart;
431 vstart.z=0.;
432 vstart.t=0.;
433 vstart.pk=1.;
434 y.push_back(vstart);
435 int niter=0; // number of iterations
436
437 // estimate first critical temperature
438 double beta=beta0(fBetaMax, tks, y, fCoolingFactor);
439 niter=0; while((update1(beta, tks,y)>1.e-6) && (niter++ < fMaxIterations)){ }
440
441 // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
442 while(beta<fBetaMax){
443
444 if(fUseTc){
445 update1(beta, tks,y);
446 while(merge(y,beta)){update1(beta, tks,y);}
447 split(beta, tks,y);
448 beta=beta/fCoolingFactor;
449 }else{
450 beta=beta/fCoolingFactor;
451 splitAll(y);
452 }
453
454 // make sure we are not too far from equilibrium before cooling further
455 niter=0; while((update1(beta, tks,y)>1.e-6) && (niter++ < fMaxIterations)){ }
456
457 }
458
459 if(fUseTc){
460 // last round of splitting, make sure no critical clusters are left
461 update1(beta, tks,y);
462 while(merge(y,beta)){update1(beta, tks,y);}
463 unsigned int ntry=0;
464 while( split(beta, tks,y) && (ntry++<10) ){
465 niter=0;
466 while((update1(beta, tks,y)>1.e-6) && (niter++ < fMaxIterations)){}
467 merge(y,beta);
468 update1(beta, tks,y);
469 }
470 }else{
471 // merge collapsed clusters
472 while(merge(y,beta)){update1(beta, tks,y);}
473 if(fVerbose ){ cout << "dump after 1st merging " << endl; dump(beta,y,tks);}
474 }
475
476 // switch on outlier rejection
477 rho0=1./nt; for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){ k->pk =1.; } // democratic
478 niter=0; while((update2(beta, tks,y,rho0, fDzCutOff) > 1.e-8) && (niter++ < fMaxIterations)){ }
479 if(fVerbose ){ cout << "rho0=" << rho0 << " niter=" << niter << endl; dump(beta,y,tks);}
480
481
482 // merge again (some cluster split by outliers collapse here)
483 while(merge(y)){}
484 if(fVerbose ){ cout << "dump after 2nd merging " << endl; dump(beta,y,tks);}
485
486
487 // continue from freeze-out to Tstop (=1) without splitting, eliminate insignificant vertices
488 while(beta<=fBetaStop){
489 while(purge(y,tks,rho0, beta, fDzCutOff)){
490 niter=0; while((update2(beta, tks, y, rho0, fDzCutOff) > 1.e-6) && (niter++ < fMaxIterations)){ }
491 }
492 beta/=fCoolingFactor;
493 niter=0; while((update2(beta, tks, y, rho0, fDzCutOff) > 1.e-6) && (niter++ < fMaxIterations)){ }
494 }
495
496
497// // new, one last round of cleaning at T=Tstop
498// while(purge(y,tks,rho0, beta)){
499// niter=0; while((update2(beta, tks,y,rho0, fDzCutOff) > 1.e-6) && (niter++ < fMaxIterations)){ }
500// }
501
502
503 if(fVerbose){
504 cout << "Final result, rho0=" << rho0 << endl;
505 dump(beta,y,tks);
506 }
507
508
509 // select significant tracks and use a TransientVertex as a container
510 //GlobalError dummyError;
511
512 // ensure correct normalization of probabilities, should make double assginment reasonably impossible
513 for(unsigned int i=0; i<nt; i++){
514 tks[i].Z=rho0*exp(-beta*( fDzCutOff*fDzCutOff));
515 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
516 tks[i].Z += k->pk * exp(-beta*Eik(tks[i],*k));
517 }
518 }
519
520 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
521
522 DelphesFactory *factory = GetFactory();
523 candidate = factory->NewCandidate();
524
525 //cout<<"new vertex"<<endl;
526 //GlobalPoint pos(0, 0, k->z);
527 double time = k->t;
528 double z = k->z;
529 //vector< reco::TransientTrack > vertexTracks;
530 //double max_track_time_err2 = 0;
531 double mean = 0.;
532 double expv_x2 = 0.;
533 double normw = 0.;
534 for(unsigned int i=0; i<nt; i++){
535 const double invdt = 1.0/std::sqrt(tks[i].dt2);
536 if(tks[i].Z>0){
537 double p = k->pk * exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
538 if( (tks[i].pi>0) && ( p > 0.5 ) ){
539 //std::cout << "pushing back " << i << ' ' << tks[i].tt << std::endl;
540 //vertexTracks.push_back(*(tks[i].tt)); tks[i].Z=0;
541
542 candidate->AddCandidate(tks[i].tt); tks[i].Z=0;
543
544 mean += tks[i].t*invdt*p;
545 expv_x2 += tks[i].t*tks[i].t*invdt*p;
546 normw += invdt*p;
547 } // setting Z=0 excludes double assignment
548 }
549 }
550
551 mean = mean/normw;
552 expv_x2 = expv_x2/normw;
553 const double time_var = expv_x2 - mean*mean;
554 const double crappy_error_guess = std::sqrt(time_var);
555 /*GlobalError dummyErrorWithTime(0,
556 0,0,
557 0,0,0,
558 0,0,0,crappy_error_guess);*/
559 //TransientVertex v(pos, time, dummyErrorWithTime, vertexTracks, 5);
560
561
562 candidate->ClusterIndex = clusterIndex++;;
563 candidate->Position.SetXYZT(0.0, 0.0, z*10.0 , time*c_light);
564
565 // TBC - fill error later ...
566 candidate->PositionError.SetXYZT(0.0, 0.0, 0.0 , crappy_error_guess*c_light);
567
568 clusterIndex++;
569 clusters.push_back(candidate);
570 }
571
572
573 return clusters;
574
575}
576
577//------------------------------------------------------------------------------
578
579static double Eik(const track_t & t, const vertex_t &k)
580{
581 return std::pow(t.z-k.z,2.)/t.dz2 + std::pow(t.t - k.t,2.)/t.dt2;
582}
583
584//------------------------------------------------------------------------------
585
586static void dump(const double beta, const vector<vertex_t> &y, const vector<track_t> &tks0)
587{
588 // copy and sort for nicer printout
589 vector<track_t> tks;
590 for(vector<track_t>::const_iterator t=tks0.begin(); t!=tks0.end(); t++){tks.push_back(*t); }
591 std::stable_sort(tks.begin(), tks.end(), recTrackLessZ1);
592
593 cout << "-----DAClusterizerInZT::dump ----" << endl;
594 cout << " beta=" << beta << endl;
595 cout << " z= ";
596 cout.precision(4);
597 for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
598 //cout << setw(8) << fixed << k->z;
599 }
600 cout << endl << " t= ";
601 for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
602 //cout << setw(8) << fixed << k->t;
603 }
604 //cout << endl << "T=" << setw(15) << 1./beta <<" Tc= ";
605 for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
606 //cout << setw(8) << fixed << k->Tc ;
607 }
608
609 cout << endl << " pk=";
610 double sumpk=0;
611 for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
612 //cout << setw(8) << setprecision(3) << fixed << k->pk;
613 sumpk+=k->pk;
614 }
615 cout << endl;
616
617 double E=0, F=0;
618 cout << endl;
619 cout << "---- z +/- dz t +/- dt ip +/-dip pt phi eta weights ----" << endl;
620 cout.precision(4);
621 for(unsigned int i=0; i<tks.size(); i++){
622 if (tks[i].Z>0){ F-=log(tks[i].Z)/beta;}
623 double tz= tks[i].z;
624 double tt= tks[i].t;
625 //cout << setw (3)<< i << ")" << setw (8) << fixed << setprecision(4)<< tz << " +/-" << setw (6)<< sqrt(tks[i].dz2)
626 // << setw(8) << fixed << setprecision(4) << tt << " +/-" << setw(6) << std::sqrt(tks[i].dt2) ;
627
628 double sump=0.;
629 for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
630 if((tks[i].pi>0)&&(tks[i].Z>0)){
631 //double p=pik(beta,tks[i],*k);
632 double p=k->pk * std::exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
633 if( p > 0.0001){
634 //cout << setw (8) << setprecision(3) << p;
635 }else{
636 cout << " . ";
637 }
638 E+=p*Eik(tks[i],*k);
639 sump+=p;
640 }else{
641 cout << " ";
642 }
643 }
644 cout << endl;
645 }
646 cout << endl << "T=" << 1/beta << " E=" << E << " n="<< y.size() << " F= " << F << endl << "----------" << endl;
647}
648
649//------------------------------------------------------------------------------
650
651static double update1(double beta, vector<track_t> &tks, vector<vertex_t> &y)
652{
653 //update weights and vertex positions
654 // mass constrained annealing without noise
655 // returns the squared sum of changes of vertex positions
656
657 unsigned int nt=tks.size();
658
659 //initialize sums
660 double sumpi=0;
661 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); ++k){
662 k->sw=0.; k->swz=0.; k->swt = 0.; k->se=0.;
663 k->swE=0.; k->Tc=0.;
664 }
665
666
667 // loop over tracks
668 for(unsigned int i=0; i<nt; i++){
669
670 // update pik and Zi
671 double Zi = 0.;
672 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); ++k){
673 k->ei = std::exp(-beta*Eik(tks[i],*k));// cache exponential for one track at a time
674 Zi += k->pk * k->ei;
675 }
676 tks[i].Z=Zi;
677
678 // normalization for pk
679 if (tks[i].Z>0){
680 sumpi += tks[i].pi;
681 // accumulate weighted z and weights for vertex update
682 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); ++k){
683 k->se += tks[i].pi* k->ei / Zi;
684 const double w = k->pk * tks[i].pi* k->ei / ( Zi * ( tks[i].dz2 * tks[i].dt2 ) );
685 k->sw += w;
686 k->swz += w * tks[i].z;
687 k->swt += w * tks[i].t;
688 k->swE += w * Eik(tks[i],*k);
689 }
690 }else{
691 sumpi += tks[i].pi;
692 }
693
694
695 } // end of track loop
696
697
698 // now update z and pk
699 double delta=0;
700 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
701 if ( k->sw > 0){
702 const double znew = k->swz/k->sw;
703 const double tnew = k->swt/k->sw;
704 delta += std::pow(k->z-znew,2.) + std::pow(k->t-tnew,2.);
705 k->z = znew;
706 k->t = tnew;
707 k->Tc = 2.*k->swE/k->sw;
708 }else{
709 // cout << " a cluster melted away ? pk=" << k->pk << " sumw=" << k->sw << endl
710 k->Tc=-1;
711 }
712
713 k->pk = k->pk * k->se / sumpi;
714 }
715
716 // return how much the prototypes moved
717 return delta;
718}
719
720//------------------------------------------------------------------------------
721
722static double update2(double beta, vector<track_t> &tks, vector<vertex_t> &y, double &rho0, double dzCutOff)
723{
724 // MVF style, no more vertex weights, update tracks weights and vertex positions, with noise
725 // returns the squared sum of changes of vertex positions
726
727 unsigned int nt=tks.size();
728
729 //initialize sums
730 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
731 k->sw = 0.; k->swz = 0.; k->swt = 0.; k->se = 0.;
732 k->swE = 0.; k->Tc=0.;
733 }
734
735
736 // loop over tracks
737 for(unsigned int i=0; i<nt; i++){
738
739 // update pik and Zi and Ti
740 double Zi = rho0*std::exp(-beta*(dzCutOff*dzCutOff));// cut-off (eventually add finite size in time)
741 //double Ti = 0.; // dt0*std::exp(-beta*fDtCutOff);
742 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
743 k->ei = std::exp(-beta*Eik(tks[i],*k));// cache exponential for one track at a time
744 Zi += k->pk * k->ei;
745 }
746 tks[i].Z=Zi;
747
748 // normalization
749 if (tks[i].Z>0){
750 // accumulate weighted z and weights for vertex update
751 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
752 k->se += tks[i].pi* k->ei / Zi;
753 double w = k->pk * tks[i].pi * k->ei /( Zi * ( tks[i].dz2 * tks[i].dt2 ) );
754 k->sw += w;
755 k->swz += w * tks[i].z;
756 k->swt += w * tks[i].t;
757 k->swE += w * Eik(tks[i],*k);
758 }
759 }
760
761 } // end of track loop
762
763 // now update z
764 double delta=0;
765 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
766 if ( k->sw > 0){
767 const double znew=k->swz/k->sw;
768 const double tnew=k->swt/k->sw;
769 delta += std::pow(k->z-znew,2.) + std::pow(k->t-tnew,2.);
770 k->z = znew;
771 k->t = tnew;
772 k->Tc = 2*k->swE/k->sw;
773 }else{
774 // cout << " a cluster melted away ? pk=" << k->pk << " sumw=" << k->sw << endl;
775 k->Tc = 0;
776 }
777
778 }
779
780 // return how much the prototypes moved
781 return delta;
782}
783
784//------------------------------------------------------------------------------
785
786static bool merge(vector<vertex_t> &y)
787{
788 // merge clusters that collapsed or never separated, return true if vertices were merged, false otherwise
789
790 if(y.size()<2) return false;
791
792 for(vector<vertex_t>::iterator k=y.begin(); (k+1)!=y.end(); k++){
793 if( std::abs( (k+1)->z - k->z ) < 1.e-3 &&
794 std::abs( (k+1)->t - k->t ) < 1.e-3 ){ // with fabs if only called after freeze-out (splitAll() at highter T)
795 double rho = k->pk + (k+1)->pk;
796 if(rho>0){
797 k->z = ( k->pk * k->z + (k+1)->z * (k+1)->pk)/rho;
798 k->t = ( k->pk * k->t + (k+1)->t * (k+1)->pk)/rho;
799 }else{
800 k->z = 0.5*(k->z + (k+1)->z);
801 k->t = 0.5*(k->t + (k+1)->t);
802 }
803 k->pk = rho;
804
805 y.erase(k+1);
806 return true;
807 }
808 }
809
810 return false;
811}
812
813//------------------------------------------------------------------------------
814
815static bool merge(vector<vertex_t> &y, double &beta)
816{
817 // merge clusters that collapsed or never separated,
818 // only merge if the estimated critical temperature of the merged vertex is below the current temperature
819 // return true if vertices were merged, false otherwise
820 if(y.size()<2) return false;
821
822 for(vector<vertex_t>::iterator k=y.begin(); (k+1)!=y.end(); k++){
823 if ( std::abs((k+1)->z - k->z) < 2.e-3 &&
824 std::abs((k+1)->t - k->t) < 2.e-3 ) {
825 double rho=k->pk + (k+1)->pk;
826 double swE=k->swE+(k+1)->swE - k->pk * (k+1)->pk / rho * ( std::pow((k+1)->z - k->z,2.) +
827 std::pow((k+1)->t - k->t,2.) );
828 double Tc=2*swE/(k->sw+(k+1)->sw);
829
830 if(Tc*beta<1){
831 if(rho>0){
832 k->z = ( k->pk * k->z + (k+1)->z * (k+1)->pk)/rho;
833 k->t = ( k->pk * k->t + (k+1)->t * (k+1)->pk)/rho;
834 }else{
835 k->z = 0.5*(k->z + (k+1)->z);
836 k->t = 0.5*(k->t + (k+1)->t);
837 }
838 k->pk = rho;
839 k->sw += (k+1)->sw;
840 k->swE = swE;
841 k->Tc = Tc;
842 y.erase(k+1);
843 return true;
844 }
845 }
846 }
847
848 return false;
849}
850
851//------------------------------------------------------------------------------
852
853static bool purge(vector<vertex_t> &y, vector<track_t> &tks, double & rho0, const double beta, const double dzCutOff)
854{
855 // eliminate clusters with only one significant/unique track
856 if(y.size()<2) return false;
857
858 unsigned int nt=tks.size();
859 double sumpmin=nt;
860 vector<vertex_t>::iterator k0=y.end();
861 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
862 int nUnique=0;
863 double sump=0;
864 double pmax=k->pk/(k->pk+rho0*exp(-beta*dzCutOff*dzCutOff));
865 for(unsigned int i=0; i<nt; i++){
866 if(tks[i].Z > 0){
867 double p = k->pk * std::exp(-beta*Eik(tks[i],*k)) / tks[i].Z ;
868 sump+=p;
869 if( (p > 0.9*pmax) && (tks[i].pi>0) ){ nUnique++; }
870 }
871 }
872
873 if((nUnique<2)&&(sump<sumpmin)){
874 sumpmin=sump;
875 k0=k;
876 }
877 }
878
879 if(k0!=y.end()){
880 //cout << "eliminating prototype at " << k0->z << "," << k0->t << " with sump=" << sumpmin << endl;
881 //rho0+=k0->pk;
882 y.erase(k0);
883 return true;
884 }else{
885 return false;
886 }
887}
888
889//------------------------------------------------------------------------------
890
891static double beta0(double betamax, vector<track_t> &tks, vector<vertex_t> &y, const double coolingFactor)
892{
893
894 double T0=0; // max Tc for beta=0
895 // estimate critical temperature from beta=0 (T=inf)
896 unsigned int nt=tks.size();
897
898 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
899
900 // vertex fit at T=inf
901 double sumwz=0.;
902 double sumwt=0.;
903 double sumw=0.;
904 for(unsigned int i=0; i<nt; i++){
905 double w = tks[i].pi/(tks[i].dz2 * tks[i].dt2);
906 sumwz += w*tks[i].z;
907 sumwt += w*tks[i].t;
908 sumw += w;
909 }
910 k->z = sumwz/sumw;
911 k->t = sumwt/sumw;
912
913 // estimate Tcrit, eventually do this in the same loop
914 double a=0, b=0;
915 for(unsigned int i=0; i<nt; i++){
916 double dx = tks[i].z-(k->z);
917 double dt = tks[i].t-(k->t);
918 double w = tks[i].pi/(tks[i].dz2 * tks[i].dt2);
919 a += w*(std::pow(dx,2.)/tks[i].dz2 + std::pow(dt,2.)/tks[i].dt2);
920 b += w;
921 }
922 double Tc= 2.*a/b; // the critical temperature of this vertex
923 if(Tc>T0) T0=Tc;
924 }// vertex loop (normally there should be only one vertex at beta=0)
925
926 if (T0>1./betamax){
927 return betamax/pow(coolingFactor, int(std::log(T0*betamax)/std::log(coolingFactor))-1 );
928 }else{
929 // ensure at least one annealing step
930 return betamax/coolingFactor;
931 }
932}
933
934//------------------------------------------------------------------------------
935
936static bool split(double beta, vector<track_t> &tks, vector<vertex_t> &y)
937{
938 // split only critical vertices (Tc >~ T=1/beta <==> beta*Tc>~1)
939 // an update must have been made just before doing this (same beta, no merging)
940 // returns true if at least one cluster was split
941
942 const double epsilon = 1e-3; // split all single vertices by 10 um
943 bool split = false;
944
945 // avoid left-right biases by splitting highest Tc first
946
947 std::vector<std::pair<double, unsigned int> > critical;
948 for(unsigned int ik=0; ik<y.size(); ik++){
949 if (beta*y[ik].Tc > 1.){
950 critical.push_back( make_pair(y[ik].Tc, ik));
951 }
952 }
953 std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int> >() );
954
955 for(unsigned int ic=0; ic<critical.size(); ic++){
956 unsigned int ik=critical[ic].second;
957 // estimate subcluster positions and weight
958 double p1=0, z1=0, t1=0, w1=0;
959 double p2=0, z2=0, t2=0, w2=0;
960 //double sumpi=0;
961 for(unsigned int i=0; i<tks.size(); i++){
962 if(tks[i].Z>0){
963 //sumpi+=tks[i].pi;
964 double p=y[ik].pk * exp(-beta*Eik(tks[i],y[ik])) / tks[i].Z*tks[i].pi;
965 double w=p/(tks[i].dz2 * tks[i].dt2);
966 if(tks[i].z < y[ik].z){
967 p1+=p; z1+=w*tks[i].z; t1+=w*tks[i].t; w1+=w;
968 }else{
969 p2+=p; z2+=w*tks[i].z; t2+=w*tks[i].t; w2+=w;
970 }
971 }
972 }
973 if(w1>0){ z1=z1/w1; t1=t1/w1;} else{ z1=y[ik].z-epsilon; t1=y[ik].t-epsilon; }
974 if(w2>0){ z2=z2/w2; t2=t2/w2;} else{ z2=y[ik].z+epsilon; t2=y[ik].t+epsilon;}
975
976 // reduce split size if there is not enough room
977 if( ( ik > 0 ) && ( y[ik-1].z>=z1 ) ){ z1=0.5*(y[ik].z+y[ik-1].z); t1=0.5*(y[ik].t+y[ik-1].t); }
978 if( ( ik+1 < y.size()) && ( y[ik+1].z<=z2 ) ){ z2=0.5*(y[ik].z+y[ik+1].z); t2=0.5*(y[ik].t+y[ik+1].t); }
979
980 // split if the new subclusters are significantly separated
981 if( (z2-z1)>epsilon || std::abs(t2-t1) > epsilon){
982 split=true;
983 vertex_t vnew;
984 vnew.pk = p1*y[ik].pk/(p1+p2);
985 y[ik].pk= p2*y[ik].pk/(p1+p2);
986 vnew.z = z1;
987 vnew.t = t1;
988 y[ik].z = z2;
989 y[ik].t = t2;
990 y.insert(y.begin()+ik, vnew);
991
992 // adjust remaining pointers
993 for(unsigned int jc=ic; jc<critical.size(); jc++){
994 if (critical[jc].second>ik) {critical[jc].second++;}
995 }
996 }
997 }
998
999 // stable_sort(y.begin(), y.end(), clusterLessZ);
1000 return split;
1001}
1002
1003//------------------------------------------------------------------------------
1004
1005void splitAll(vector<vertex_t> &y)
1006{
1007
1008
1009 const double epsilon=1e-3; // split all single vertices by 10 um
1010 const double zsep=2*epsilon; // split vertices that are isolated by at least zsep (vertices that haven't collapsed)
1011 const double tsep=2*epsilon; // check t as well
1012
1013 vector<vertex_t> y1;
1014
1015 for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
1016 if ( ( (k==y.begin())|| (k-1)->z < k->z - zsep) && (((k+1)==y.end() )|| (k+1)->z > k->z + zsep)) {
1017 // isolated prototype, split
1018 vertex_t vnew;
1019 vnew.z = k->z - epsilon;
1020 vnew.t = k->t - epsilon;
1021 (*k).z = k->z + epsilon;
1022 (*k).t = k->t + epsilon;
1023 vnew.pk= 0.5* (*k).pk;
1024 (*k).pk= 0.5* (*k).pk;
1025 y1.push_back(vnew);
1026 y1.push_back(*k);
1027
1028 }else if( y1.empty() || (y1.back().z < k->z -zsep) || (y1.back().t < k->t - tsep) ){
1029 y1.push_back(*k);
1030 }else{
1031 y1.back().z -= epsilon;
1032 y1.back().t -= epsilon;
1033 k->z += epsilon;
1034 k->t += epsilon;
1035 y1.push_back(*k);
1036 }
1037 }// vertex loop
1038
1039 y=y1;
1040}
1041
1042
1043
Note: See TracBrowser for help on using the repository browser.