CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DAClusterizerInZ.cc
Go to the documentation of this file.
5 
6 
7 using namespace std;
8 
9 
10 namespace {
11 
12  bool recTrackLessZ1(const DAClusterizerInZ::track_t & tk1,
13  const DAClusterizerInZ::track_t & tk2)
14  {
15  return tk1.z < tk2.z;
16  }
17 }
18 
19 
20 vector<DAClusterizerInZ::track_t> DAClusterizerInZ::fill(
21  const vector<reco::TransientTrack> & tracks
22  )const{
23  // prepare track data for clustering
24  vector<track_t> tks;
25  for(vector<reco::TransientTrack>::const_iterator it=tracks.begin(); it!=tracks.end(); it++){
26  track_t t;
27  t.z=((*it).stateAtBeamLine().trackStateAtPCA()).position().z();
28  double tantheta=tan(((*it).stateAtBeamLine().trackStateAtPCA()).momentum().theta());
29  double phi=((*it).stateAtBeamLine().trackStateAtPCA()).momentum().phi();
30  // get the beam-spot
31  reco::BeamSpot beamspot=(it->stateAtBeamLine()).beamSpot();
32  t.dz2= pow((*it).track().dzError(),2) // track errror
33  + (pow(beamspot.BeamWidthX()*cos(phi),2)+pow(beamspot.BeamWidthY()*sin(phi),2))/pow(tantheta,2) // beam-width induced
34  + pow(vertexSize_,2); // intrinsic vertex size, safer for outliers and short lived decays
35  if (d0CutOff_>0){
36  Measurement1D IP=(*it).stateAtBeamLine().transverseImpactParameter();// error constains beamspot
37  t.pi=1./(1.+exp(pow(IP.value()/IP.error(),2)-pow(d0CutOff_ ,2))); // reduce weight for high ip tracks
38  }else{
39  t.pi=1.;
40  }
41  t.tt=&(*it);
42  t.Z=1.;
43  tks.push_back(t);
44  }
45  return tks;
46 }
47 
48 
49 
50 
51 
52 double DAClusterizerInZ::Eik(const track_t & t, const vertex_t &k )const{
53  return pow(t.z-k.z,2)/t.dz2;
54 }
55 
56 
57 
58 
60  double beta,
61  vector<track_t> & tks,
62  vector<vertex_t> & y
63  )const{
64  //update weights and vertex positions
65  // mass constrained annealing without noise
66  // returns the squared sum of changes of vertex positions
67 
68  unsigned int nt=tks.size();
69 
70  //initialize sums
71  double sumpi=0;
72  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
73  k->se=0; k->sw=0; k->swz=0;
74  k->swE=0; k->Tc=0;
75  }
76 
77 
78  // loop over tracks
79  for(unsigned int i=0; i<nt; i++){
80 
81  // update pik and Zi
82  double Zi=0.;
83  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
84  k->ei=exp(-beta*Eik(tks[i],*k));// cache exponential for one track at a time
85  Zi += k->pk * k->ei;
86  }
87  tks[i].Z=Zi;
88 
89  // normalization for pk
90  if (tks[i].Z>0){
91  sumpi += tks[i].pi;
92  // accumulate weighted z and weights for vertex update
93  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
94  k->se += tks[i].pi* k->ei / Zi;
95  double w = k->pk * tks[i].pi* k->ei / Zi / tks[i].dz2;
96  k->sw += w;
97  k->swz += w * tks[i].z;
98  k->swE+= w*Eik(tks[i],*k);
99  }
100  }else{
101  sumpi += tks[i].pi;
102  }
103 
104 
105  } // end of track loop
106 
107 
108  // now update z and pk
109  double delta=0;
110  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
111  if ( k->sw > 0){
112  double znew=k->swz/k->sw;
113  delta+=pow(k->z-znew,2);
114  k->z=znew;
115  k->Tc=2*k->swE/k->sw;
116  }else{
117  edm::LogInfo("sumw") << "invalid sum of weights in fit: " << k->sw << endl;
118  if(verbose_){cout << " a cluster melted away ? pk=" << k->pk << " sumw=" << k->sw << endl;}
119  k->Tc=-1;
120  }
121 
122  k->pk = k->pk * k->se / sumpi;
123  }
124 
125  // return how much the prototypes moved
126  return delta;
127 }
128 
129 
130 
131 
132 
134  double beta,
135  vector<track_t> & tks,
136  vector<vertex_t> & y,
137  double & rho0
138  )const{
139  // MVF style, no more vertex weights, update tracks weights and vertex positions, with noise
140  // returns the squared sum of changes of vertex positions
141 
142  unsigned int nt=tks.size();
143 
144  //initialize sums
145  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
146  k->se=0; k->sw=0; k->swz=0;
147  k->swE=0; k->Tc=0;
148  }
149 
150 
151  // loop over tracks
152  for(unsigned int i=0; i<nt; i++){
153 
154  // update pik and Zi
155  double Zi=rho0*exp(-beta*dzCutOff_*dzCutOff_);// cut-off
156  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
157  k->ei=exp(-beta*Eik(tks[i],*k));// cache exponential for one track at a time
158  Zi += k->pk * k->ei;
159  }
160  tks[i].Z=Zi;
161 
162  // normalization
163  if (tks[i].Z>0){
164  // accumulate weighted z and weights for vertex update
165  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
166  k->se += tks[i].pi* k->ei / Zi;
167  double w = k->pk * tks[i].pi* k->ei / Zi / tks[i].dz2;
168  k->sw += w;
169  k->swz += w * tks[i].z;
170  k->swE+= w*Eik(tks[i],*k);
171  }
172  }
173 
174 
175  } // end of track loop
176 
177 
178  // now update z
179  double delta=0;
180  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
181  if ( k->sw > 0){
182  double znew=k->swz/k->sw;
183  delta+=pow(k->z-znew,2);
184  k->z=znew;
185  k->Tc=2*k->swE/k->sw;
186  }else{
187  edm::LogInfo("sumw") << "invalid sum of weights in fit: " << k->sw << endl;
188  if(verbose_){cout << " a cluster melted away ? pk=" << k->pk << " sumw=" << k->sw << endl;}
189  k->Tc=0;
190  }
191 
192  }
193 
194  // return how much the prototypes moved
195  return delta;
196 }
197 
198 
199 
200 
201 
202 bool DAClusterizerInZ::merge(vector<vertex_t> & y, int nt)const{
203  // merge clusters that collapsed or never separated, return true if vertices were merged, false otherwise
204 
205  if(y.size()<2) return false;
206 
207  for(vector<vertex_t>::iterator k=y.begin(); (k+1)!=y.end(); k++){
208  if( fabs( (k+1)->z - k->z ) < 1.e-3 ){ // with fabs if only called after freeze-out (splitAll() at highter T)
209  double rho=k->pk + (k+1)->pk;
210  if(rho>0){
211  k->z=( k->pk * k->z + (k+1)->z * (k+1)->pk)/rho;
212  }else{
213  k->z=0.5*(k->z + (k+1)->z);
214  }
215  k->pk=rho;
216 
217  y.erase(k+1);
218  return true;
219  }
220  }
221 
222  return false;
223 }
224 
225 
226 
227 
228 bool DAClusterizerInZ::merge(vector<vertex_t> & y, double & beta)const{
229  // merge clusters that collapsed or never separated,
230  // only merge if the estimated critical temperature of the merged vertex is below the current temperature
231  // return true if vertices were merged, false otherwise
232  if(y.size()<2) return false;
233 
234  for(vector<vertex_t>::iterator k=y.begin(); (k+1)!=y.end(); k++){
235  if (fabs((k+1)->z - k->z)<2.e-3){
236  double rho=k->pk + (k+1)->pk;
237  double swE=k->swE+(k+1)->swE - k->pk * (k+1)->pk /rho *pow((k+1)->z - k->z,2);
238  double Tc=2*swE/(k->sw+(k+1)->sw);
239 
240  if(Tc*beta<1){
241  if(rho>0){
242  k->z=( k->pk * k->z + (k+1)->z * (k+1)->pk)/rho;
243  }else{
244  k->z=0.5*(k->z + (k+1)->z);
245  }
246  k->pk=rho;
247  k->sw+=(k+1)->sw;
248  k->swE=swE;
249  k->Tc=Tc;
250  y.erase(k+1);
251  return true;
252  }
253  }
254  }
255 
256  return false;
257 }
258 
259 
260 
261 
262 
263 bool DAClusterizerInZ::purge(vector<vertex_t> & y, vector<track_t> & tks, double & rho0, const double beta)const{
264  // eliminate clusters with only one significant/unique track
265  if(y.size()<2) return false;
266 
267  unsigned int nt=tks.size();
268  double sumpmin=nt;
269  vector<vertex_t>::iterator k0=y.end();
270  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
271  int nUnique=0;
272  double sump=0;
273  double pmax=k->pk/(k->pk+rho0*exp(-beta*dzCutOff_*dzCutOff_));
274  for(unsigned int i=0; i<nt; i++){
275  if(tks[i].Z>0){
276  double p=k->pk * exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
277  sump+=p;
278  if( (p > 0.9*pmax) && (tks[i].pi>0) ){ nUnique++; }
279  }
280  }
281 
282  if((nUnique<2)&&(sump<sumpmin)){
283  sumpmin=sump;
284  k0=k;
285  }
286  }
287 
288  if(k0!=y.end()){
289  if(verbose_){cout << "eliminating prototype at " << k0->z << " with sump=" << sumpmin << endl;}
290  //rho0+=k0->pk;
291  y.erase(k0);
292  return true;
293  }else{
294  return false;
295  }
296 }
297 
298 
299 
300 
302  double betamax,
303  vector<track_t> & tks,
304  vector<vertex_t> & y
305  )const{
306 
307  double T0=0; // max Tc for beta=0
308  // estimate critical temperature from beta=0 (T=inf)
309  unsigned int nt=tks.size();
310 
311  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
312 
313  // vertex fit at T=inf
314  double sumwz=0;
315  double sumw=0;
316  for(unsigned int i=0; i<nt; i++){
317  double w=tks[i].pi/tks[i].dz2;
318  sumwz+=w*tks[i].z;
319  sumw +=w;
320  }
321  k->z=sumwz/sumw;
322 
323  // estimate Tcrit, eventually do this in the same loop
324  double a=0, b=0;
325  for(unsigned int i=0; i<nt; i++){
326  double dx=tks[i].z-(k->z);
327  double w=tks[i].pi/tks[i].dz2;
328  a+=w*pow(dx,2)/tks[i].dz2;
329  b+=w;
330  }
331  double Tc= 2.*a/b; // the critical temperature of this vertex
332  if(Tc>T0) T0=Tc;
333  }// vertex loop (normally there should be only one vertex at beta=0)
334 
335  if (T0>1./betamax){
336  return betamax/pow(coolingFactor_, int(log(T0*betamax)/log(coolingFactor_))-1 );
337  }else{
338  // ensure at least one annealing step
339  return betamax/coolingFactor_;
340  }
341 }
342 
343 
344 
345 
346 
348  double beta,
349  vector<track_t> & tks,
350  vector<vertex_t> & y,
351  double threshold
352  )const{
353  // split only critical vertices (Tc >~ T=1/beta <==> beta*Tc>~1)
354  // an update must have been made just before doing this (same beta, no merging)
355  // returns true if at least one cluster was split
356 
357  double epsilon=1e-3; // split all single vertices by 10 um
358  bool split=false;
359 
360 
361  // avoid left-right biases by splitting highest Tc first
362 
363  std::vector<std::pair<double, unsigned int> > critical;
364  for(unsigned int ik=0; ik<y.size(); ik++){
365  if (beta*y[ik].Tc > 1.){
366  critical.push_back( make_pair(y[ik].Tc, ik));
367  }
368  }
369  stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int> >() );
370 
371 
372  for(unsigned int ic=0; ic<critical.size(); ic++){
373  unsigned int ik=critical[ic].second;
374  // estimate subcluster positions and weight
375  double p1=0, z1=0, w1=0;
376  double p2=0, z2=0, w2=0;
377  //double sumpi=0;
378  for(unsigned int i=0; i<tks.size(); i++){
379  if(tks[i].Z>0){
380  //sumpi+=tks[i].pi;
381  double p=y[ik].pk * exp(-beta*Eik(tks[i],y[ik])) / tks[i].Z*tks[i].pi;
382  double w=p/tks[i].dz2;
383  if(tks[i].z<y[ik].z){
384  p1+=p; z1+=w*tks[i].z; w1+=w;
385  }else{
386  p2+=p; z2+=w*tks[i].z; w2+=w;
387  }
388  }
389  }
390  if(w1>0){ z1=z1/w1;} else{z1=y[ik].z-epsilon;}
391  if(w2>0){ z2=z2/w2;} else{z2=y[ik].z+epsilon;}
392 
393  // reduce split size if there is not enough room
394  if( ( ik > 0 ) && ( y[ik-1].z>=z1 ) ){ z1=0.5*(y[ik].z+y[ik-1].z); }
395  if( ( ik+1 < y.size()) && ( y[ik+1].z<=z2 ) ){ z2=0.5*(y[ik].z+y[ik+1].z); }
396 
397  // split if the new subclusters are significantly separated
398  if( (z2-z1)>epsilon){
399  split=true;
400  vertex_t vnew;
401  vnew.pk = p1*y[ik].pk/(p1+p2);
402  y[ik].pk= p2*y[ik].pk/(p1+p2);
403  vnew.z = z1;
404  y[ik].z = z2;
405  y.insert(y.begin()+ik, vnew);
406 
407  // adjust remaining pointers
408  for(unsigned int jc=ic; jc<critical.size(); jc++){
409  if (critical[jc].second>ik) {critical[jc].second++;}
410  }
411  }
412  }
413 
414  // stable_sort(y.begin(), y.end(), clusterLessZ);
415  return split;
416 }
417 
418 
419 
420 
421 
423  vector<vertex_t> & y
424  )const{
425 
426 
427  double epsilon=1e-3; // split all single vertices by 10 um
428  double zsep=2*epsilon; // split vertices that are isolated by at least zsep (vertices that haven't collapsed)
429  vector<vertex_t> y1;
430 
431  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
432  if ( ( (k==y.begin())|| (k-1)->z < k->z - zsep) && (((k+1)==y.end() )|| (k+1)->z > k->z + zsep)) {
433  // isolated prototype, split
434  vertex_t vnew;
435  vnew.z = k->z -epsilon;
436  (*k).z = k->z+epsilon;
437  vnew.pk= 0.5* (*k).pk;
438  (*k).pk= 0.5* (*k).pk;
439  y1.push_back(vnew);
440  y1.push_back(*k);
441 
442  }else if( y1.empty() || (y1.back().z < k->z -zsep)){
443  y1.push_back(*k);
444  }else{
445  y1.back().z -=epsilon;
446  k->z+=epsilon;
447  y1.push_back(*k);
448  }
449  }// vertex loop
450 
451  y=y1;
452 }
453 
454 
456 {
457  // some defaults to avoid uninitialized variables
458  verbose_= conf.getUntrackedParameter<bool>("verbose", false);
459  useTc_=true;
460  betamax_=0.1;
461  betastop_ =1.0;
462  coolingFactor_=0.8;
463  maxIterations_=100;
464  vertexSize_=0.05; // 0.5 mm
465  dzCutOff_=4.0; // Adaptive Fitter uses 3.0 but that appears to be a bit tight here sometimes
466 
467  // configure
468 
469  double Tmin = conf.getParameter<double>("Tmin");
470  vertexSize_ = conf.getParameter<double>("vertexSize");
471  coolingFactor_ = conf.getParameter<double>("coolingFactor");
472  d0CutOff_ = conf.getParameter<double>("d0CutOff");
473  dzCutOff_ = conf.getParameter<double>("dzCutOff");
474  maxIterations_=100;
475  if (Tmin==0){
476  cout << "DAClusterizerInZ: invalid Tmin" << Tmin << " reset do default " << 1./betamax_ << endl;
477  }else{
478  betamax_ = 1./Tmin;
479  }
480 
481  // for testing, negative cooling factor: revert to old splitting scheme
482  if(coolingFactor_<0){
483  coolingFactor_=-coolingFactor_; useTc_=false;
484  }
485 
486 }
487 
488 
489 void DAClusterizerInZ::dump(const double beta, const vector<vertex_t> & y, const vector<track_t> & tks0, int verbosity)const{
490 
491  // copy and sort for nicer printout
492  vector<track_t> tks;
493  for(vector<track_t>::const_iterator t=tks0.begin(); t!=tks0.end(); t++){tks.push_back(*t); }
494  stable_sort(tks.begin(), tks.end(), recTrackLessZ1);
495 
496  cout << "-----DAClusterizerInZ::dump ----" << endl;
497  cout << "beta=" << beta << " betamax= " << betamax_ << endl;
498  cout << " z= ";
499  cout.precision(4);
500  for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
501  cout << setw(8) << fixed << k->z ;
502  }
503  cout << endl << "T=" << setw(15) << 1./beta <<" Tc= ";
504  for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
505  cout << setw(8) << fixed << k->Tc ;
506  }
507 
508  cout << endl << " pk=";
509  double sumpk=0;
510  for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
511  cout << setw(8) << setprecision(3) << fixed << k->pk;
512  sumpk+=k->pk;
513  }
514  cout << endl;
515 
516  if(verbosity>0){
517  double E=0, F=0;
518  cout << endl;
519  cout << "---- z +/- dz ip +/-dip pt phi eta weights ----" << endl;
520  cout.precision(4);
521  for(unsigned int i=0; i<tks.size(); i++){
522  if (tks[i].Z>0){ F-=log(tks[i].Z)/beta;}
523  double tz= tks[i].z;
524  cout << setw (3)<< i << ")" << setw (8) << fixed << setprecision(4)<< tz << " +/-" << setw (6)<< sqrt(tks[i].dz2);
525 
526  if(tks[i].tt->track().quality(reco::TrackBase::highPurity)){ cout << " *";}else{cout <<" ";}
527  if(tks[i].tt->track().hitPattern().hasValidHitInFirstPixelBarrel()){cout <<"+";}else{cout << "-";}
528  cout << setw(1) << tks[i].tt->track().hitPattern().pixelBarrelLayersWithMeasurement(); // see DataFormats/TrackReco/interface/HitPattern.h
529  cout << setw(1) << tks[i].tt->track().hitPattern().pixelEndcapLayersWithMeasurement();
530  cout << setw(1) << hex << tks[i].tt->track().hitPattern().trackerLayersWithMeasurement()-tks[i].tt->track().hitPattern().pixelLayersWithMeasurement() <<dec;
531  cout << "=" << setw(1)<<hex <<tks[i].tt->track().trackerExpectedHitsOuter().numberOfHits() << dec;
532 
533  Measurement1D IP=tks[i].tt->stateAtBeamLine().transverseImpactParameter();
534  cout << setw (8) << IP.value() << "+/-" << setw (6) << IP.error();
535  cout << " " << setw(6) << setprecision(2) << tks[i].tt->track().pt()*tks[i].tt->track().charge();
536  cout << " " << setw(5) << setprecision(2) << tks[i].tt->track().phi()
537  << " " << setw(5) << setprecision(2) << tks[i].tt->track().eta() ;
538 
539  double sump=0.;
540  for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
541  if((tks[i].pi>0)&&(tks[i].Z>0)){
542  //double p=pik(beta,tks[i],*k);
543  double p=k->pk * exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
544  if( p > 0.0001){
545  cout << setw (8) << setprecision(3) << p;
546  }else{
547  cout << " . ";
548  }
549  E+=p*Eik(tks[i],*k);
550  sump+=p;
551  }else{
552  cout << " ";
553  }
554  }
555  cout << endl;
556  }
557  cout << endl << "T=" << 1/beta << " E=" << E << " n="<< y.size() << " F= " << F << endl << "----------" << endl;
558  }
559 }
560 
561 
562 
563 
564 
565 vector< TransientVertex >
566 DAClusterizerInZ::vertices(const vector<reco::TransientTrack> & tracks, const int verbosity)
567 const
568 {
569 
570  vector<track_t> tks=fill(tracks);
571  unsigned int nt=tracks.size();
572  double rho0=0.0; // start with no outlier rejection
573 
574  vector< TransientVertex > clusters;
575  if (tks.empty()) return clusters;
576 
577  vector<vertex_t> y; // the vertex prototypes
578 
579  // initialize:single vertex at infinite temperature
580  vertex_t vstart;
581  vstart.z=0.;
582  vstart.pk=1.;
583  y.push_back(vstart);
584  int niter=0; // number of iterations
585 
586 
587  // estimate first critical temperature
588  double beta=beta0(betamax_, tks, y);
589  niter=0; while((update(beta, tks,y)>1.e-6) && (niter++ < maxIterations_)){ }
590 
591 
592  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
593  while(beta<betamax_){
594 
595  if(useTc_){
596  update(beta, tks,y);
597  while(merge(y,beta)){update(beta, tks,y);}
598  split(beta, tks,y, 1.);
599  beta=beta/coolingFactor_;
600  }else{
601  beta=beta/coolingFactor_;
602  splitAll(y);
603  }
604 
605 
606  // make sure we are not too far from equilibrium before cooling further
607  niter=0; while((update(beta, tks,y)>1.e-6) && (niter++ < maxIterations_)){ }
608 
609  }
610 
611  if(useTc_){
612  // last round of splitting, make sure no critical clusters are left
613  update(beta, tks,y);
614  while(merge(y,beta)){update(beta, tks,y);}
615  unsigned int ntry=0;
616  while( split(beta, tks,y,1.) && (ntry++<10) ){
617  niter=0;
618  while((update(beta, tks,y)>1.e-6) && (niter++ < maxIterations_)){}
619  merge(y,beta);
620  update(beta, tks,y);
621  }
622  }else{
623  // merge collapsed clusters
624  while(merge(y,beta)){update(beta, tks,y);}
625  if(verbose_ ){ cout << "dump after 1st merging " << endl; dump(beta,y,tks,2);}
626  }
627 
628 
629 
630 
631  // switch on outlier rejection
632  rho0=1./nt; for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){ k->pk =1.; } // democratic
633  niter=0; while((update(beta, tks,y,rho0) > 1.e-8) && (niter++ < maxIterations_)){ }
634  if(verbose_ ){ cout << "rho0=" << rho0 << " niter=" << niter << endl; dump(beta,y,tks,2);}
635 
636 
637  // merge again (some cluster split by outliers collapse here)
638  while(merge(y,tks.size())){}
639  if(verbose_ ){ cout << "dump after 2nd merging " << endl; dump(beta,y,tks,2);}
640 
641 
642  // continue from freeze-out to Tstop (=1) without splitting, eliminate insignificant vertices
643  while(beta<=betastop_){
644  while(purge(y,tks,rho0, beta)){
645  niter=0; while((update(beta, tks,y,rho0) > 1.e-6) && (niter++ < maxIterations_)){ }
646  }
647  beta/=coolingFactor_;
648  niter=0; while((update(beta, tks,y,rho0) > 1.e-6) && (niter++ < maxIterations_)){ }
649  }
650 
651 
652 // // new, one last round of cleaning at T=Tstop
653 // while(purge(y,tks,rho0, beta)){
654 // niter=0; while((update(beta, tks,y,rho0) > 1.e-6) && (niter++ < maxIterations_)){ }
655 // }
656 
657 
658  if(verbose_){
659  cout << "Final result, rho0=" << rho0 << endl;
660  dump(beta,y,tks,2);
661  }
662 
663 
664  // select significant tracks and use a TransientVertex as a container
665  GlobalError dummyError;
666 
667 
668  // ensure correct normalization of probabilities, should make double assginment reasonably impossible
669  for(unsigned int i=0; i<nt; i++){
670  tks[i].Z=rho0*exp(-beta*dzCutOff_*dzCutOff_);
671  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
672  tks[i].Z+=k->pk * exp(-beta*Eik(tks[i],*k));
673  }
674  }
675 
676 
677  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
678  GlobalPoint pos(0, 0, k->z);
679  vector< reco::TransientTrack > vertexTracks;
680  for(unsigned int i=0; i<nt; i++){
681  if(tks[i].Z>0){
682  double p=k->pk * exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
683  if( (tks[i].pi>0) && ( p > 0.5 ) ){ vertexTracks.push_back(*(tks[i].tt)); tks[i].Z=0; } // setting Z=0 excludes double assignment
684  }
685  }
686  TransientVertex v(pos, dummyError, vertexTracks, 5);
687  clusters.push_back(v);
688  }
689 
690 
691  return clusters;
692 
693 }
694 
695 
696 
697 
698 
699 vector< vector<reco::TransientTrack> >
700 DAClusterizerInZ::clusterize(const vector<reco::TransientTrack> & tracks)
701  const
702 {
703  if(verbose_) {
704  cout << "###################################################" << endl;
705  cout << "# DAClusterizerInZ::clusterize nt="<<tracks.size() << endl;
706  cout << "###################################################" << endl;
707  }
708 
709  vector< vector<reco::TransientTrack> > clusters;
710  vector< TransientVertex > pv=vertices(tracks);
711 
712  if(verbose_){ cout << "# DAClusterizerInZ::clusterize pv.size="<<pv.size() << endl; }
713  if (pv.size()==0){ return clusters;}
714 
715 
716  // fill into clusters and merge
717  vector< reco::TransientTrack> aCluster=pv.begin()->originalTracks();
718 
719  for(vector<TransientVertex>::iterator k=pv.begin()+1; k!=pv.end(); k++){
720  if ( fabs(k->position().z() - (k-1)->position().z()) > (2*vertexSize_) ){
721  // close a cluster
722  clusters.push_back(aCluster);
723  aCluster.clear();
724  }
725  for(unsigned int i=0; i<k->originalTracks().size(); i++){ aCluster.push_back( k->originalTracks().at(i)); }
726 
727  }
728  clusters.push_back(aCluster);
729 
730 
731  return clusters;
732 
733 }
734 
const double Z[kNumberCalorimeter]
const double beta
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
string fill
Definition: lumiContext.py:319
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Definition: DDAxes.h:10
double error() const
Definition: Measurement1D.h:30
std::vector< track_t > fill(const std::vector< reco::TransientTrack > &tracks) const
float float float z
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
U second(std::pair< T, U > const &p)
bool merge(std::vector< vertex_t > &, int) const
bool split(double beta, std::vector< track_t > &tks, std::vector< vertex_t > &y, double threshold) const
T sqrt(T t)
Definition: SSEVec.h:48
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
tuple IP
Definition: listHistos.py:76
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:87
double p2[4]
Definition: TauolaWrapper.h:90
tuple conf
Definition: dbtoconf.py:185
int nt
Definition: AMPTWrapper.h:32
int k[5][pyjets_maxn]
DAClusterizerInZ(const edm::ParameterSet &conf)
tuple tracks
Definition: testEve_cfg.py:39
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:89
const reco::TransientTrack * tt
double b
Definition: hdecay.h:120
void dump(const double beta, const std::vector< vertex_t > &y, const std::vector< track_t > &tks, const int verbosity=0) const
double value() const
Definition: Measurement1D.h:28
void splitAll(std::vector< vertex_t > &y) const
double update(double beta, std::vector< track_t > &tks, std::vector< vertex_t > &y) const
double p1[4]
Definition: TauolaWrapper.h:89
#define update(a, b)
double Eik(const track_t &t, const vertex_t &k) const
double a
Definition: hdecay.h:121
tuple cout
Definition: gather_cfg.py:121
T w() const
const double epsilon
double pi
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
double beta0(const double betamax, std::vector< track_t > &tks, std::vector< vertex_t > &y) const
bool purge(std::vector< vertex_t > &, std::vector< track_t > &, double &, const double) const
double split
Definition: MVATrainer.cc:139
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks, const int verbosity=0) const
Definition: DDAxes.h:10