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  // get the beam-spot
30  reco::BeamSpot beamspot=(it->stateAtBeamLine()).beamSpot();
31  t.dz2= pow((*it).track().dzError(),2) // track errror
32  + (pow(beamspot.BeamWidthX(),2)+pow(beamspot.BeamWidthY(),2))/pow(tantheta,2) // beam-width induced
33  + pow(vertexSize_,2); // intrinsic vertex size, safer for outliers and short lived decays
34  Measurement1D IP=(*it).stateAtBeamLine().transverseImpactParameter();// error constains beamspot
35  t.pi=1./(1.+exp(pow(IP.value()/IP.error(),2)-pow(3.,2))); // reduce weight for high ip tracks
36  t.tt=&(*it);
37  t.Z=1.;
38  tks.push_back(t);
39  }
40  return tks;
41 }
42 
43 
44 
45 
46 // double DAClusterizerInZ::pik(const double beta, const track_t & t, const vertex_t &k )const{
47 // //note: t.Z = sum_k exp(-beta*Eik ) is assumed to be valid, this is done in updateAndFit
48 // // then we have sum_k pik(beta, t, k) = 1
49 // // the last call of updateAndFit must have been made with the same temperature !
50 // // at low T, the vertex position must be quite accurate (iterated)
51 // double Eik=pow(t.z-k.z,2)/t.dz2;
52 // double o= k.pk*exp(-beta*Eik);
53 // if(isnan(o)){ cout << "NaN in pik " << k.pk << " " << beta << " " << Eik << " Z=" << t.Z << " beta Eik=" << beta*Eik << endl;}
54 // if (t.Z>0) return k.pk*exp(-beta*Eik)/t.Z;
55 // return 0;
56 // }
57 
58 
59 double DAClusterizerInZ::Eik(const track_t & t, const vertex_t &k )const{
60  return pow(t.z-k.z,2)/t.dz2;
61 }
62 
63 
64 
65 
67  double beta,
68  vector<track_t> & tks,
69  vector<vertex_t> & y
70  )const{
71  //update weights and vertex positions
72  // mass constrained annealing without noise
73  // returns the squared sum of changes of vertex positions
74 
75  //cout << endl << endl << "DAClusterizerInZ::update" << endl;
76  unsigned int nt=tks.size();
77 
78  //initialize sums
79  double sumpi=0;
80  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
81  k->se=0; k->sw=0; k->swz=0;
82  }
83 
84 
85  // loop over tracks
86  for(unsigned int i=0; i<nt; i++){
87 
88  // update pik and Zi
89  double Zi=0.;
90  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
91  k->ei=exp(-beta*Eik(tks[i],*k));// cache exponential for one track at a time
92  Zi += k->pk * k->ei;
93  }
94  tks[i].Z=Zi;
95 
96  // normalization for pk
97  if (tks[i].Z>0){
98  sumpi += tks[i].pi;
99  // accumulate weighted z and weights for vertex update
100  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
101  k->se += tks[i].pi* k->ei / Zi;
102  double w = k->pk * tks[i].pi* k->ei / Zi / tks[i].dz2;
103  k->sw += w;
104  k->swz += w * tks[i].z;
105  }
106  }else{
107  sumpi += tks[i].pi;
108  }
109 
110 
111  } // end of track loop
112 
113 
114  // now update z and pk
115  double delta=0;
116  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
117  if ( k->sw > 0){
118  double znew=k->swz/k->sw;
119  delta+=pow(k->z-znew,2);
120  k->z=znew;
121  }else{
122  edm::LogInfo("sumw") << "invalid sum of weights in fit: " << k->sw << endl;
123  if(verbose_){cout << " a cluster melted away ? pk=" << k->pk << " sumw=" << k->sw << endl;}
124  }
125 
126  k->pk = k->pk * k->se / sumpi;
127  }
128 
129  // return how much the prototypes moved
130  return delta;
131 }
132 
133 
134 
135 
136 
138  double beta,
139  vector<track_t> & tks,
140  vector<vertex_t> & y,
141  double & rho0
142  )const{
143  // MVF style, no more vertex weights, update tracks weights and vertex positions, with noise
144  // returns the squared sum of changes of vertex positions
145 
146  unsigned int nt=tks.size();
147 
148  //initialize sums
149  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
150  k->se=0; k->sw=0; k->swz=0;
151  }
152 
153 
154  // loop over tracks
155  for(unsigned int i=0; i<nt; i++){
156 
157  // update pik and Zi
158  double Zi=rho0*exp(-beta*mu0_*mu0_);// cut-off
159  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
160  k->ei=exp(-beta*Eik(tks[i],*k));// cache exponential for one track at a time
161  Zi += k->pk * k->ei;
162  }
163  tks[i].Z=Zi;
164 
165  // normalization
166  if (tks[i].Z>0){
167  //double pi0=exp(-beta*mu0_*mu0_)/tks[i].Z;
168  //sumpi += tks[i].pi*(1.-pi0*rho0); // exclude rho0 from the normalization of the cluster "mass"
169  // accumulate weighted z and weights for vertex update
170  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
171  k->se += tks[i].pi* k->ei / Zi;
172  double w = k->pk * tks[i].pi* k->ei / Zi / tks[i].dz2;
173  k->sw += w;
174  k->swz += w * tks[i].z;
175  }
176  }
177  //else{ sumpi += tks[i].pi; }
178 
179 
180  } // end of track loop
181 
182 
183  // now update z and pk
184  double delta=0;
185  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
186  if ( k->sw > 0){
187  double znew=k->swz/k->sw;
188  delta+=pow(k->z-znew,2);
189  k->z=znew;
190  }else{
191  edm::LogInfo("sumw") << "invalid sum of weights in fit: " << k->sw << endl;
192  if(verbose_){cout << " a cluster melted away ? pk=" << k->pk << " sumw=" << k->sw << endl;}
193  }
194 
195  //k->pk *= k->se / sumpi;
196  }
197 
198  // return how much the prototypes moved
199  return delta;
200 }
201 
202 
203 
204 
205 
206 bool DAClusterizerInZ::merge(vector<vertex_t> & y, int nt)const{
207  // merge clusters that collapsed or never separated, return true if vertices were merged, false otherwise
208 
209  if(y.size()<2) return false;
210 
211  for(vector<vertex_t>::iterator k=y.begin(); (k+1)!=y.end(); k++){
212  //if ((k+1)->z - k->z<1.e-2){ // note, no fabs here, maintains z-ordering (with split()+merge() at every temperature)
213  if (fabs((k+1)->z - k->z)<1.e-2){ // with fabs if only called after freeze-out (splitAll() at highter T)
214  k->pk += (k+1)->pk;
215  k->z=0.5*(k->z+(k+1)->z);
216  y.erase(k+1);
217  return true;
218  }
219  }
220 
221  return false;
222 }
223 
224 
225 
226 bool DAClusterizerInZ::purge(vector<vertex_t> & y, vector<track_t> & tks, double & rho0, const double beta)const{
227  // eliminate clusters with only one significant/unique track
228  if(y.size()<2) return false;
229 
230  unsigned int nt=tks.size();
231  double sumpmin=nt;
232  vector<vertex_t>::iterator k0=y.end();
233  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
234  int nUnique=0;
235  double sump=0;
236  double pmax=k->pk/(k->pk+rho0*exp(-beta*mu0_*mu0_));
237  for(unsigned int i=0; i<nt; i++){
238  if(tks[i].Z>0){
239  //double p=pik(beta,tks[i],*k);
240  double p=k->pk * exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
241  sump+=p;
242  if( (p > 0.9*pmax) && (tks[i].pi>0) ){ nUnique++; }
243  }
244  }
245 
246  if((nUnique<2)&&(sump<sumpmin)){
247  sumpmin=sump;
248  k0=k;
249  }
250  }
251 
252  if(k0!=y.end()){
253  if(verbose_){cout << "eliminating prototype at " << k0->z << " with sump=" << sumpmin << endl;}
254  //rho0+=k0->pk;
255  y.erase(k0);
256  return true;
257  }else{
258  return false;
259  }
260 }
261 
262 
263 
264 
266  double betamax,
267  vector<track_t> & tks,
268  vector<vertex_t> & y
269  )const{
270 
271  double T0=0; // max Tc for beta=0
272  // estimate critical temperature from beta=0 (T=inf)
273  unsigned int nt=tks.size();
274 
275  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
276 
277  // vertex fit at T=inf
278  double sumwz=0;
279  double sumw=0;
280  for(unsigned int i=0; i<nt; i++){
281  double w=tks[i].pi/tks[i].dz2;
282  sumwz+=w*tks[i].z;
283  sumw +=w;
284  }
285  k->z=sumwz/sumw;
286 
287  // estimate Tcrit, eventually do this in the same loop
288  double a=0, b=0;
289  for(unsigned int i=0; i<nt; i++){
290  double dx=tks[i].z-(k->z);
291  double w=tks[i].pi/tks[i].dz2;
292  a+=w*pow(dx,2)/tks[i].dz2;
293  b+=w;
294  }
295  double Tc= 2.*a/b; // the critical temperature of this vertex
296  if(Tc>T0) T0=Tc;
297  }// vertex loop (normally there should be only one vertex at beta=0)
298 
299  if (T0>1./betamax){
300  return betamax/pow(coolingFactor_, int(log(T0*betamax)/log(coolingFactor_))-1 );
301  }else{
302  // ensure at least one annealing step
303  return betamax/coolingFactor_;
304  }
305 }
306 
307 
308 
309 
311  vector<track_t> & tks,
312  vector<vertex_t> & y
313  )const{
314 
315 
316  double epsilon=1e-3; // split all single vertices by 10 um
317  double zsep=2*epsilon; // split vertices that are isolated by at least zsep (vertices that haven't collapsed)
318  vector<vertex_t> y1;
319 
320  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
321  if ( ( (k==y.begin())|| (k-1)->z < k->z - zsep) && (((k+1)==y.end() )|| (k+1)->z > k->z + zsep)) {
322  // isolated prototype, split
323  vertex_t vnew;
324  vnew.z = k->z -epsilon;
325  (*k).z = k->z+epsilon;
326  vnew.pk= 0.5* (*k).pk;
327  (*k).pk= 0.5* (*k).pk;
328  y1.push_back(vnew);
329  y1.push_back(*k);
330 
331  }else if( y1.empty() || (y1.back().z < k->z -zsep)){
332  y1.push_back(*k);
333  }else{
334  y1.back().z -=epsilon;
335  k->z+=epsilon;
336  y1.push_back(*k);
337  }
338  }// vertex loop
339 
340  y=y1;
341 }
342 
343 
345 {
346  // some defaults to avoid uninitialized variables
347  verbose_= conf.getUntrackedParameter<bool>("verbose", false);
348  betamax_=0.1;
349  betastop_ =1.0;
350  coolingFactor_=0.8;
351  maxIterations_=100;
352  vertexSize_=0.05; // 0.5 mm
353  mu0_=4.0; // Adaptive Fitter uses 3.0 but that appears to be a bit tight here sometimes
354 
355  // configure
356 
357  double Tmin = conf.getParameter<double>("Tmin");
358  vertexSize_ = conf.getParameter<double>("vertexSize");
359  coolingFactor_ = conf.getParameter<double>("coolingFactor");
360  maxIterations_=100;
361  if (Tmin==0){
362  cout << "DAClusterizerInZ: invalid Tmin" << Tmin << " reset do default " << 1./betamax_ << endl;
363  }else{
364  betamax_ = 1./Tmin;
365  }
366 
367 }
368 
369 
370 void DAClusterizerInZ::dump(const double beta, const vector<vertex_t> & y, const vector<track_t> & tks0, int verbosity)const{
371 
372  // copy and sort for nicer printout
373  vector<track_t> tks;
374  for(vector<track_t>::const_iterator t=tks0.begin(); t!=tks0.end(); t++){tks.push_back(*t); }
375  stable_sort(tks.begin(), tks.end(), recTrackLessZ1);
376 
377  cout << "-----DAClusterizerInZ::dump ----" << endl;
378  cout << "beta=" << beta << " betamax= " << betamax_ << endl;
379  cout << " z= ";
380  cout.precision(4);
381  for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
382  cout << setw(8) << fixed << k->z ;
383  }
384  cout << endl << "T=" << setw(15) << 1./beta <<" Tc= ";
385  cout << endl << " pk=";
386  double sumpk=0;
387  for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
388  cout << setw(8) << setprecision(3) << fixed << k->pk;
389  sumpk+=k->pk;
390  }
391  cout << endl;
392 
393  if(verbosity>0){
394  double E=0, F=0;
395  cout << endl;
396  cout << "---- z +/- dz ip +/-dip pt phi eta weights ----" << endl;
397  cout.precision(4);
398  for(unsigned int i=0; i<tks.size(); i++){
399  if (tks[i].Z>0){ F-=log(tks[i].Z)/beta;}
400  double tz= tks[i].z;
401  cout << setw (3)<< i << ")" << setw (8) << fixed << setprecision(4)<< tz << " +/-" << setw (6)<< sqrt(tks[i].dz2);
402 
403  if(tks[i].tt->track().quality(reco::TrackBase::highPurity)){ cout << " *";}else{cout <<" ";}
404  if(tks[i].tt->track().hitPattern().hasValidHitInFirstPixelBarrel()){cout <<"+";}else{cout << "-";}
405  cout << setw(1) << tks[i].tt->track().hitPattern().pixelBarrelLayersWithMeasurement(); // see DataFormats/TrackReco/interface/HitPattern.h
406  cout << setw(1) << tks[i].tt->track().hitPattern().pixelEndcapLayersWithMeasurement();
407  cout << setw(1) << hex << tks[i].tt->track().hitPattern().trackerLayersWithMeasurement()-tks[i].tt->track().hitPattern().pixelLayersWithMeasurement() <<dec;
408  cout << "=" << setw(1)<<hex <<tks[i].tt->track().trackerExpectedHitsOuter().numberOfHits() << dec;
409 
410  Measurement1D IP=tks[i].tt->stateAtBeamLine().transverseImpactParameter();
411  cout << setw (8) << IP.value() << "+/-" << setw (6) << IP.error();
412  cout << " " << setw(6) << setprecision(2) << tks[i].tt->track().pt()*tks[i].tt->track().charge();
413  cout << " " << setw(5) << setprecision(2) << tks[i].tt->track().phi()
414  << " " << setw(5) << setprecision(2) << tks[i].tt->track().eta() ;
415 
416  double sump=0.;
417  for(vector<vertex_t>::const_iterator k=y.begin(); k!=y.end(); k++){
418  if((tks[i].pi>0)&&(tks[i].Z>0)){
419  //double p=pik(beta,tks[i],*k);
420  double p=k->pk * exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
421  if( p > 0.0001){
422  cout << setw (8) << setprecision(3) << p;
423  }else{
424  cout << " . ";
425  }
426  E+=p*Eik(tks[i],*k);
427  sump+=p;
428  }else{
429  cout << " ";
430  }
431  }
432  cout << endl;
433  }
434  cout << endl << "T=" << 1/beta << " E=" << E << " n="<< y.size() << " F= " << F << endl << "----------" << endl;
435  }
436 }
437 
438 
439 
440 
441 
442 vector< TransientVertex >
443 DAClusterizerInZ::vertices(const vector<reco::TransientTrack> & tracks, const int verbosity)
444 const
445 {
446 
447  vector<track_t> tks=fill(tracks);
448  unsigned int nt=tracks.size();
449  double rho0=0.0; // start with no outlier rejection
450 
451  vector< TransientVertex > clusters;
452  if (tks.empty()) return clusters;
453 
454  vector<vertex_t> y; // the vertex prototypes
455 
456  // initialize:single vertex at infinite temperature
457  vertex_t vstart;
458  vstart.z=0.;
459  vstart.pk=1.;
460  y.push_back(vstart);
461  int niter=0; // number of iterations
462 
463 
464  // estimate first critical temperature
465  double beta=beta0(betamax_, tks, y);
466  niter=0; while((update(beta, tks,y)>1.e-6) && (niter++ < maxIterations_)){ }
467 
468 
469  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
470  while(beta<betamax_){
471 
472  beta=beta/coolingFactor_;
473  splitAll(tks,y);
474 
475  // make sure we are not too far from equilibrium before cooling further
476  niter=0; while((update(beta, tks,y)>1.e-6) && (niter++ < maxIterations_)){ }
477 
478  }
479 
480 
481  // merge collapsed clusters
482  while(merge(y,tks.size())){}
483  if(verbose_ ){ cout << "dump after 1st merging " << endl; dump(beta,y,tks,2);}
484 
485 
486 
487  // switch on outlier rejection
488  //rho0=exp(beta*mu0_*mu0_)/nt; if(rho0>0.1) rho0=0.1;
489  rho0=1./nt; for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){ k->pk =1.; } // democratic
490  niter=0; while((update(beta, tks,y,rho0) > 1.e-8) && (niter++ < maxIterations_)){ }
491  if(verbose_ ){ cout << "rho0=" << rho0 << " niter=" << niter << endl; dump(beta,y,tks,2);}
492 
493 
494  // merge again (some cluster split by outliers collapse here)
495  while(merge(y,tks.size())){}
496  if(verbose_ ){ cout << "dump after 2nd merging " << endl; dump(beta,y,tks,2);}
497 
498 
499  // continue from freeze-out to Tstop (=1) without splitting, eliminate insignificant vertices
500  while(beta<=betastop_){
501  while(purge(y,tks,rho0, beta)){
502  niter=0; while((update(beta, tks,y,rho0) > 1.e-6) && (niter++ < maxIterations_)){ }
503  }
504  beta/=coolingFactor_;
505  niter=0; while((update(beta, tks,y,rho0) > 1.e-6) && (niter++ < maxIterations_)){ }
506  }
507 
508  if(verbose_){
509  cout << "Final result, rho0=" << rho0 << endl;
510  dump(beta,y,tks,2);
511  }
512 
513 
514  // select significant tracks and use a TransientVertex as a container
515  GlobalError dummyError;
516  for(vector<vertex_t>::iterator k=y.begin(); k!=y.end(); k++){
517  GlobalPoint pos(0, 0, k->z);
518  vector< reco::TransientTrack > vertexTracks;
519  for(unsigned int i=0; i<nt; i++){
520  if(tks[i].Z>0){
521  double p=k->pk * exp(-beta*Eik(tks[i],*k)) / tks[i].Z;
522  if( (tks[i].pi>0) && ( p > 0.5 ) ){ vertexTracks.push_back(*(tks[i].tt));}
523  }
524  }
525  TransientVertex v(pos, dummyError, vertexTracks, 0);
526  clusters.push_back(v);
527  }
528 
529  return clusters;
530 
531 }
532 
533 
534 
535 
536 
537 vector< vector<reco::TransientTrack> >
538 DAClusterizerInZ::clusterize(const vector<reco::TransientTrack> & tracks)
539  const
540 {
541  if(verbose_) {
542  cout << "###################################################" << endl;
543  cout << "# DAClusterizerInZ::clusterize nt="<<tracks.size() << endl;
544  cout << "###################################################" << endl;
545  }
546 
547  vector< vector<reco::TransientTrack> > clusters;
548  vector< TransientVertex > pv=vertices(tracks);
549 
550  if(verbose_){ cout << "# DAClusterizerInZ::clusterize pv.size="<<pv.size() << endl; }
551  if (pv.size()==0){ return clusters;}
552 
553 
554  // fill into clusters and merge
555  vector< reco::TransientTrack> aCluster=pv.begin()->originalTracks();
556 
557  for(vector<TransientVertex>::iterator k=pv.begin()+1; k!=pv.end(); k++){
558  if ( k->position().z() - (k-1)->position().z()> (2*vertexSize_) ){
559  // close a cluster
560  clusters.push_back(aCluster);
561  aCluster.clear();
562  }
563  for(unsigned int i=0; i<k->originalTracks().size(); i++){ aCluster.push_back( k->originalTracks().at(i)); }
564 
565  }
566  clusters.push_back(aCluster);
567 
568 
569  return clusters;
570 
571 }
572 
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
double error() const
Definition: Measurement1D.h:30
std::vector< track_t > fill(const std::vector< reco::TransientTrack > &tracks) const
Exp< T >::type exp(const T &t)
Definition: Exp.h:22
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
bool merge(std::vector< vertex_t > &, int) const
void splitAll(std::vector< track_t > &tks, std::vector< vertex_t > &y) const
Definition: DDAxes.h:10
T sqrt(T t)
Definition: SSEVec.h:28
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
tuple conf
Definition: dbtoconf.py:185
int nt
Definition: AMPTWrapper.h:32
int k[5][pyjets_maxn]
const int verbosity
Log< T >::type log(const T &t)
Definition: Log.h:22
DAClusterizerInZ(const edm::ParameterSet &conf)
tuple tracks
Definition: testEve_cfg.py:39
bool merge(LuminosityBlockRange &lh, LuminosityBlockRange &rh)
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
double update(double beta, std::vector< track_t > &tks, std::vector< vertex_t > &y) const
#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:41
const double epsilon
double pi
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
mathSSE::Vec4< T > v
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