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