CMS 3D CMS Logo

DAClusterizerInZ_vect.cc
Go to the documentation of this file.
5 
6 #include <cmath>
7 #include <cassert>
8 #include <limits>
9 #include <iomanip>
11 #include "vdt/vdtMath.h"
12 
13 using namespace std;
14 
16 
17  // hardcoded parameters
18  maxIterations_ = 100;
19  mintrkweight_ = 0.5; // conf.getParameter<double>("mintrkweight");
20 
21 
22  // configurable debug outptut debug output
23  verbose_ = conf.getUntrackedParameter<bool> ("verbose", false);
24  zdumpcenter_ = conf.getUntrackedParameter<double> ("zdumpcenter", 0.);
25  zdumpwidth_ = conf.getUntrackedParameter<double> ("zdumpwidth", 20.);
26 
27  // configurable parameters
28  double Tmin = conf.getParameter<double> ("Tmin");
29  double Tpurge = conf.getParameter<double> ("Tpurge");
30  double Tstop = conf.getParameter<double> ("Tstop");
31  vertexSize_ = conf.getParameter<double> ("vertexSize");
32  coolingFactor_ = conf.getParameter<double> ("coolingFactor");
33  useTc_=true;
34  if(coolingFactor_<0){
35  coolingFactor_=-coolingFactor_;
36  useTc_=false;
37  }
38  d0CutOff_ = conf.getParameter<double> ("d0CutOff");
39  dzCutOff_ = conf.getParameter<double> ("dzCutOff");
40  uniquetrkweight_ = conf.getParameter<double>("uniquetrkweight");
41  zmerge_ = conf.getParameter<double>("zmerge");
42 
43  if(verbose_){
44  std::cout << "DAClusterizerinZ_vect: mintrkweight = " << mintrkweight_ << std::endl;
45  std::cout << "DAClusterizerinZ_vect: uniquetrkweight = " << uniquetrkweight_ << std::endl;
46  std::cout << "DAClusterizerinZ_vect: zmerge = " << zmerge_ << std::endl;
47  std::cout << "DAClusterizerinZ_vect: Tmin = " << Tmin << std::endl;
48  std::cout << "DAClusterizerinZ_vect: Tpurge = " << Tpurge << std::endl;
49  std::cout << "DAClusterizerinZ_vect: Tstop = " << Tstop << std::endl;
50  std::cout << "DAClusterizerinZ_vect: vertexSize = " << vertexSize_ << std::endl;
51  std::cout << "DAClusterizerinZ_vect: coolingFactor = " << coolingFactor_ << std::endl;
52  std::cout << "DAClusterizerinZ_vect: d0CutOff = " << d0CutOff_ << std::endl;
53  std::cout << "DAClusterizerinZ_vect: dzCutOff = " << dzCutOff_ << std::endl;
54  }
55 
56 
57  if (Tmin == 0) {
58  edm::LogWarning("DAClusterizerinZ_vectorized") << "DAClusterizerInZ: invalid Tmin" << Tmin
59  << " reset do default " << 1. / betamax_;
60  } else {
61  betamax_ = 1. / Tmin;
62  }
63 
64 
65  if ((Tpurge > Tmin) || (Tpurge == 0)) {
66  edm::LogWarning("DAClusterizerinZ_vectorized") << "DAClusterizerInZ: invalid Tpurge" << Tpurge
67  << " set to " << Tmin;
68  Tpurge = Tmin;
69  }
70  betapurge_ = 1./Tpurge;
71 
72 
73  if ((Tstop > Tpurge) || (Tstop == 0)) {
74  edm::LogWarning("DAClusterizerinZ_vectorized") << "DAClusterizerInZ: invalid Tstop" << Tstop
75  << " set to " << max(1., Tpurge);
76  Tstop = max(1., Tpurge) ;
77  }
78  betastop_ = 1./Tstop;
79 
80 }
81 
82 
83 namespace {
84  inline double local_exp( double const& inp) {
85  return vdt::fast_exp( inp );
86  }
87 
88  inline void local_exp_list( double const * __restrict__ arg_inp, double * __restrict__ arg_out, const int arg_arr_size) {
89  for(int i=0; i!=arg_arr_size; ++i ) arg_out[i]=vdt::fast_exp(arg_inp[i]);
90  }
91 
92 }
93 
94 //todo: use r-value possibility of c++11 here
96 DAClusterizerInZ_vect::fill(const vector<reco::TransientTrack> & tracks) const {
97 
98  // prepare track data for clustering
99  track_t tks;
100  for (auto it = tracks.begin(); it!= tracks.end(); it++){
101  if (!(*it).isValid()) continue;
102  double t_pi=1.;
103  double t_z = ((*it).stateAtBeamLine().trackStateAtPCA()).position().z();
104  if (std::fabs(t_z) > 1000.) continue;
105  auto const & t_mom = (*it).stateAtBeamLine().trackStateAtPCA().momentum();
106  // get the beam-spot
107  reco::BeamSpot beamspot = (it->stateAtBeamLine()).beamSpot();
108  double t_dz2 =
109  std::pow((*it).track().dzError(), 2) // track errror
110  + (std::pow(beamspot.BeamWidthX()*t_mom.x(),2)+std::pow(beamspot.BeamWidthY()*t_mom.y(),2))*std::pow(t_mom.z(),2)/std::pow(t_mom.perp2(),2) // beam spot width
111  + std::pow(vertexSize_, 2); // intrinsic vertex size, safer for outliers and short lived decays
112  t_dz2 = 1./ t_dz2;
113  if (edm::isNotFinite(t_dz2) || t_dz2 < std::numeric_limits<double>::min() ) continue;
114  if (d0CutOff_ > 0) {
116  (*it).stateAtBeamLine().transverseImpactParameter();// error contains beamspot
117  t_pi = 1. / (1. + local_exp(std::pow(atIP.value() / atIP.error(), 2) - std::pow(d0CutOff_, 2))); // reduce weight for high ip tracks
118  if (edm::isNotFinite(t_pi) || t_pi < std::numeric_limits<double>::epsilon()) continue; // usually is > 0.99
119  }
120  LogTrace("DAClusterizerinZ_vectorized") << t_z <<' '<< t_dz2 <<' '<< t_pi;
121  tks.AddItem(t_z, t_dz2, &(*it), t_pi);
122  }
123  tks.ExtractRaw();
124 
125  if (verbose_) {
126  std::cout << "Track count " << tks.GetSize() << std::endl;
127  }
128 
129  return tks;
130 }
131 
132 
133 namespace {
134  inline
135  double Eik(double t_z, double k_z, double t_dz2) {
136  return std::pow(t_z - k_z, 2) * t_dz2;
137  }
138 }
139 
140 double DAClusterizerInZ_vect::update(double beta, track_t & gtracks,
141  vertex_t & gvertices, bool useRho0, const double & rho0) const {
142 
143  //update weights and vertex positions
144  // mass constrained annealing without noise
145  // returns the squared sum of changes of vertex positions
146 
147  const unsigned int nt = gtracks.GetSize();
148  const unsigned int nv = gvertices.GetSize();
149 
150  //initialize sums
151  double sumpi = 0;
152 
153  // to return how much the prototype moved
154  double delta = 0;
155 
156 
157  // intial value of a sum
158  double Z_init = 0;
159  // independpent of loop
160  if ( useRho0 )
161  {
162  Z_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_); // cut-off
163  }
164 
165  // define kernels
166  auto kernel_calc_exp_arg = [ beta, nv ] ( const unsigned int itrack,
167  track_t const& tracks,
168  vertex_t const& vertices ) {
169  const double track_z = tracks._z[itrack];
170  const double botrack_dz2 = -beta*tracks._dz2[itrack];
171 
172  // auto-vectorized
173  for ( unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
174  auto mult_res = track_z - vertices._z[ivertex];
175  vertices._ei_cache[ivertex] = botrack_dz2 * ( mult_res * mult_res );
176  }
177  };
178 
179  auto kernel_add_Z = [ nv, Z_init ] (vertex_t const& vertices) -> double
180  {
181  double ZTemp = Z_init;
182  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
183  ZTemp += vertices._pk[ivertex] * vertices._ei[ivertex];
184  }
185  return ZTemp;
186  };
187 
188  auto kernel_calc_normalization = [ beta, nv ] (const unsigned int track_num,
189  track_t & tks_vec,
190  vertex_t & y_vec ) {
191  auto tmp_trk_pi = tks_vec._pi[track_num];
192  auto o_trk_Z_sum = 1./tks_vec._Z_sum[track_num];
193  auto o_trk_dz2 = tks_vec._dz2[track_num];
194  auto tmp_trk_z = tks_vec._z[track_num];
195  auto obeta = -1./beta;
196 
197  // auto-vectorized
198  for (unsigned int k = 0; k < nv; ++k) {
199  y_vec._se[k] += y_vec._ei[k] * (tmp_trk_pi* o_trk_Z_sum);
200  auto w = y_vec._pk[k] * y_vec._ei[k] * (tmp_trk_pi*o_trk_Z_sum *o_trk_dz2);
201  y_vec._sw[k] += w;
202  y_vec._swz[k] += w * tmp_trk_z;
203  y_vec._swE[k] += w * y_vec._ei_cache[k]*obeta;
204  }
205  };
206 
207 
208  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
209  gvertices._se[ivertex] = 0.0;
210  gvertices._sw[ivertex] = 0.0;
211  gvertices._swz[ivertex] = 0.0;
212  gvertices._swE[ivertex] = 0.0;
213  }
214 
215 
216 
217  // loop over tracks
218  for (auto itrack = 0U; itrack < nt; ++itrack) {
219  kernel_calc_exp_arg(itrack, gtracks, gvertices);
220  local_exp_list(gvertices._ei_cache, gvertices._ei, nv);
221 
222  gtracks._Z_sum[itrack] = kernel_add_Z(gvertices);
223  if (edm::isNotFinite(gtracks._Z_sum[itrack])) gtracks._Z_sum[itrack] = 0.0;
224  // used in the next major loop to follow
225  sumpi += gtracks._pi[itrack];
226 
227  if (gtracks._Z_sum[itrack] > 1.e-100){
228  kernel_calc_normalization(itrack, gtracks, gvertices);
229  }
230  }
231 
232  // now update z and pk
233  auto kernel_calc_z = [ sumpi, nv
234 #ifdef VI_DEBUG
235  , this
236 #endif
237  ] (vertex_t & vertices ) -> double {
238 
239  double delta=0;
240  // does not vectorizes
241  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex ) {
242  if (vertices._sw[ivertex] > 0) {
243  auto znew = vertices._swz[ ivertex ] / vertices._sw[ ivertex ];
244  // prevents from vectorizing if
245  delta += std::pow( vertices._z[ ivertex ] - znew, 2 );
246  vertices._z[ ivertex ] = znew;
247  }
248 #ifdef VI_DEBUG
249  else {
250  edm::LogInfo("sumw") << "invalid sum of weights in fit: " << vertices._sw[ivertex] << endl;
251  if (this->verbose_) {
252  std::cout << " a cluster melted away ? pk=" << vertices._pk[ ivertex ] << " sumw="
253  << vertices._sw[ivertex] << endl;
254  }
255  }
256 #endif
257  }
258 
259  auto osumpi = 1./sumpi;
260  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex )
261  vertices._pk[ ivertex ] = vertices._pk[ ivertex ] * vertices._se[ ivertex ] * osumpi;
262 
263  return delta;
264  };
265 
266  delta += kernel_calc_z(gvertices);
267 
268  // return how much the prototypes moved
269  return delta;
270 }
271 
272 
273 
274 
275 
276 bool DAClusterizerInZ_vect::merge(vertex_t & y, double & beta)const{
277  // merge clusters that collapsed or never separated,
278  // only merge if the estimated critical temperature of the merged vertex is below the current temperature
279  // return true if vertices were merged, false otherwise
280  const unsigned int nv = y.GetSize();
281 
282  if (nv < 2)
283  return false;
284 
285  // merge the smallest distance clusters first
286  std::vector<std::pair<double, unsigned int> > critical;
287  for (unsigned int k = 0; (k + 1) < nv; k++) {
288  if (std::fabs(y._z[k + 1] - y._z[k]) < zmerge_) {
289  critical.push_back( make_pair( std::fabs(y._z[k + 1] - y._z[k]), k) );
290  }
291  }
292  if (critical.empty()) return false;
293 
294  std::stable_sort(critical.begin(), critical.end(), std::less<std::pair<double, unsigned int> >() );
295 
296 
297  for (unsigned int ik=0; ik < critical.size(); ik++){
298  unsigned int k = critical[ik].second;
299  double rho = y._pk[k]+y._pk[k+1];
300  double swE = y._swE[k]+y._swE[k+1]-y._pk[k]*y._pk[k+1] / rho*std::pow(y._z[k+1]-y._z[k],2);
301  double Tc = 2*swE / (y._sw[k]+y._sw[k+1]);
302 
303  if(Tc*beta < 1){
304  if(verbose_){ std::cout << "merging " << y._z[k + 1] << " and " << y._z[k] << " Tc = " << Tc << " sw = " << y._sw[k]+y._sw[k+1] <<std::endl;}
305  if(rho > 0){
306  y._z[k] = (y._pk[k]*y._z[k] + y._pk[k+1]*y._z[k + 1])/rho;
307  }else{
308  y._z[k] = 0.5 * (y._z[k] + y._z[k + 1]);
309  }
310  y._pk[k] = rho;
311  y._sw[k] += y._sw[k+1];
312  y._swE[k] = swE;
313  y.RemoveItem(k+1);
314  return true;
315  }
316  }
317 
318  return false;
319 }
320 
321 
322 
323 
324 bool
325 DAClusterizerInZ_vect::purge(vertex_t & y, track_t & tks, double & rho0, const double beta) const {
326  // eliminate clusters with only one significant/unique track
327  const unsigned int nv = y.GetSize();
328  const unsigned int nt = tks.GetSize();
329 
330  if (nv < 2)
331  return false;
332 
333  double sumpmin = nt;
334  unsigned int k0 = nv;
335 
336  for (unsigned int k = 0; k < nv; k++) {
337 
338  int nUnique = 0;
339  double sump = 0;
340 
341  double pmax = y._pk[k] / (y._pk[k] + rho0 * local_exp(-beta * dzCutOff_* dzCutOff_));
342  for (unsigned int i = 0; i < nt; i++) {
343  if (tks._Z_sum[i] > 1.e-100) {
344  double p = y._pk[k] * local_exp(-beta * Eik(tks._z[i], y._z[k], tks._dz2[i])) / tks._Z_sum[i];
345  sump += p;
346  if ((p > uniquetrkweight_ * pmax) && (tks._pi[i] > 0)) {
347  nUnique++;
348  }
349  }
350  }
351 
352  if ((nUnique < 2) && (sump < sumpmin)) {
353  sumpmin = sump;
354  k0 = k;
355  }
356 
357  }
358 
359  if (k0 != nv) {
360  if (verbose_) {
361  std::cout << "eliminating prototype at " << std::setw(10) << std::setprecision(4) << y._z[k0]
362  << " with sump=" << sumpmin
363  << " rho*nt =" << y._pk[k0]*nt
364  << endl;
365  }
366  y.RemoveItem(k0);
367  return true;
368  } else {
369  return false;
370  }
371 }
372 
373 
374 
375 
376 double
377 DAClusterizerInZ_vect::beta0(double betamax, track_t const & tks, vertex_t const & y) const {
378 
379  double T0 = 0; // max Tc for beta=0
380  // estimate critical temperature from beta=0 (T=inf)
381  const unsigned int nt = tks.GetSize();
382  const unsigned int nv = y.GetSize();
383 
384  for (unsigned int k = 0; k < nv; k++) {
385 
386  // vertex fit at T=inf
387  double sumwz = 0;
388  double sumw = 0;
389  for (unsigned int i = 0; i < nt; i++) {
390  double w = tks._pi[i] * tks._dz2[i];
391  sumwz += w * tks._z[i];
392  sumw += w;
393  }
394  y._z[k] = sumwz / sumw;
395 
396  // estimate Tcrit, eventually do this in the same loop
397  double a = 0, b = 0;
398  for (unsigned int i = 0; i < nt; i++) {
399  double dx = tks._z[i] - y._z[k];
400  double w = tks._pi[i] * tks._dz2[i];
401  a += w * std::pow(dx, 2) * tks._dz2[i];
402  b += w;
403  }
404  double Tc = 2. * a / b; // the critical temperature of this vertex
405  if (Tc > T0) T0 = Tc;
406  }// vertex loop (normally there should be only one vertex at beta=0)
407 
408  if(verbose_){
409  std::cout << "DAClustrizerInZ_vect.beta0: Tc = " << T0 << std::endl;
410  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
411  std::cout << "DAClustrizerInZ_vect.beta0: nstep = " << coolingsteps << std::endl;
412  }
413 
414 
415  if (T0 > 1. / betamax) {
416  return betamax / std::pow(coolingFactor_, int(std::log(T0 * betamax) / std::log(coolingFactor_)) - 1);
417  } else {
418  // ensure at least one annealing step
419  return betamax * coolingFactor_;
420  }
421 }
422 
423 
424 
425 bool
426 DAClusterizerInZ_vect::split(const double beta, track_t &tks, vertex_t & y, double threshold ) const{
427  // split only critical vertices (Tc >~ T=1/beta <==> beta*Tc>~1)
428  // an update must have been made just before doing this (same beta, no merging)
429  // returns true if at least one cluster was split
430 
431  double epsilon=1e-3; // minimum split size
432  unsigned int nv = y.GetSize();
433 
434  // avoid left-right biases by splitting highest Tc first
435 
436  std::vector<std::pair<double, unsigned int> > critical;
437  for(unsigned int k=0; k<nv; k++){
438  double Tc= 2*y._swE[k]/y._sw[k];
439  if (beta*Tc > threshold){
440  critical.push_back( make_pair(Tc, k));
441  }
442  }
443  if (critical.empty()) return false;
444 
445 
446  std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int> >() );
447 
448 
449  bool split=false;
450  const unsigned int nt = tks.GetSize();
451 
452  for(unsigned int ic=0; ic<critical.size(); ic++){
453  unsigned int k=critical[ic].second;
454 
455  // estimate subcluster positions and weight
456  double p1=0, z1=0, w1=0;
457  double p2=0, z2=0, w2=0;
458  for(unsigned int i=0; i<nt; i++){
459  if (tks._Z_sum[i] > 1.e-100) {
460 
461  // winner-takes-all, usually overestimates splitting
462  double tl = tks._z[i] < y._z[k] ? 1.: 0.;
463  double tr = 1. - tl;
464 
465  // soften it, especially at low T
466  double arg = (tks._z[i] - y._z[k]) * sqrt(beta * tks._dz2[i]);
467  if(std::fabs(arg) < 20){
468  double t = local_exp(-arg);
469  tl = t/(t+1.);
470  tr = 1/(t+1.);
471  }
472 
473  double p = y._pk[k] * tks._pi[i] * local_exp(-beta * Eik(tks._z[i], y._z[k], tks._dz2[i])) / tks._Z_sum[i];
474  double w = p*tks._dz2[i];
475  p1 += p*tl ; z1 += w*tl*tks._z[i]; w1 += w*tl;
476  p2 += p*tr; z2 += w*tr*tks._z[i]; w2 += w*tr;
477  }
478  }
479 
480  if(w1>0){z1 = z1/w1;} else {z1=y._z[k]-epsilon;}
481  if(w2>0){z2 = z2/w2;} else {z2=y._z[k]+epsilon;}
482 
483  // reduce split size if there is not enough room
484  if( ( k > 0 ) && ( z1 < (0.6*y._z[k] + 0.4*y._z[k-1])) ){ z1 = 0.6*y._z[k] + 0.4*y._z[k-1]; }
485  if( ( k+1 < nv) && ( z2 > (0.6*y._z[k] + 0.4*y._z[k+1])) ){ z2 = 0.6*y._z[k] + 0.4*y._z[k+1]; }
486 
487  if(verbose_){
488  if (std::fabs(y._z[k] - zdumpcenter_) < zdumpwidth_){
489  std::cout << " T= " << std::setw(8) << 1./beta
490  << " Tc= " << critical[ic].first
491  << " splitting " << std::fixed << std::setprecision(4) << y._z[k]
492  << " --> " << z1 << "," << z2
493  << " [" << p1 << "," << p2 << "]" ;
494  if (std::fabs(z2-z1)>epsilon){
495  std::cout << std::endl;
496  }else{
497  std::cout << " rejected " << std::endl;
498  }
499  }
500  }
501 
502  // split if the new subclusters are significantly separated
503  if( (z2-z1) > epsilon){
504  split = true;
505  double pk1 = p1*y._pk[k]/(p1+p2);
506  double pk2 = p2*y._pk[k]/(p1+p2);
507  y._z[k] = z2;
508  y._pk[k] = pk2;
509  y.InsertItem(k, z1, pk1);
510  nv++;
511 
512  // adjust remaining pointers
513  for(unsigned int jc=ic; jc < critical.size(); jc++){
514  if (critical[jc].second > k) {critical[jc].second++;}
515  }
516  }
517  }
518  return split;
519 }
520 
521 
522 
524 
525  const unsigned int nv = y.GetSize();
526 
527  double epsilon = 1e-3; // split all single vertices by 10 um
528  double zsep = 2 * epsilon; // split vertices that are isolated by at least zsep (vertices that haven't collapsed)
529  vertex_t y1;
530 
531  if (verbose_) {
532  std::cout << "Before Split "<< std::endl;
533  y.DebugOut();
534  }
535 
536  for (unsigned int k = 0; k < nv; k++) {
537  if (
538  ( (k == 0) || ( y._z[k - 1] < (y._z[k] - zsep)) ) &&
539  ( ((k + 1) == nv) || ( y._z[k + 1] > (y._z[k] + zsep)) ) )
540  {
541  // isolated prototype, split
542  double new_z = y.z[k] - epsilon;
543  y.z[k] = y.z[k] + epsilon;
544 
545  double new_pk = 0.5 * y.pk[k];
546  y.pk[k] = 0.5 * y.pk[k];
547 
548  y1.AddItem(new_z, new_pk);
549  y1.AddItem(y._z[k], y._pk[k]);
550  }
551  else if ( (y1.GetSize() == 0 ) ||
552  (y1._z[y1.GetSize() - 1] < (y._z[k] - zsep) ))
553  {
554  y1.AddItem(y._z[k], y._pk[k]);
555  }
556  else
557  {
558  y1._z[y1.GetSize() - 1] = y1._z[y1.GetSize() - 1] - epsilon;
559  y._z[k] = y._z[k] + epsilon;
560  y1.AddItem( y._z[k] , y._pk[k]);
561  }
562  }// vertex loop
563 
564  y = y1;
565  y.ExtractRaw();
566 
567  if (verbose_) {
568  std::cout << "After split " << std::endl;
569  y.DebugOut();
570  }
571 }
572 
573 
574 
575 vector<TransientVertex>
576 DAClusterizerInZ_vect::vertices(const vector<reco::TransientTrack> & tracks, const int verbosity) const {
577  track_t && tks = fill(tracks);
578  tks.ExtractRaw();
579 
580  unsigned int nt = tks.GetSize();
581  double rho0 = 0.0; // start with no outlier rejection
582 
583  vector<TransientVertex> clusters;
584  if (tks.GetSize() == 0) return clusters;
585 
586  vertex_t y; // the vertex prototypes
587 
588  // initialize:single vertex at infinite temperature
589  y.AddItem( 0, 1.0);
590 
591  int niter = 0; // number of iterations
592 
593 
594  // estimate first critical temperature
595  double beta = beta0(betamax_, tks, y);
596  if ( verbose_) std::cout << "Beta0 is " << beta << std::endl;
597 
598  niter = 0;
599  while ((update(beta, tks, y, false, rho0) > 1.e-6) &&
600  (niter++ < maxIterations_)) {}
601 
602  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
603 
604  double betafreeze = betamax_ * sqrt(coolingFactor_);
605 
606  while (beta < betafreeze) {
607  if(useTc_){
608  update(beta, tks,y, false, rho0);
609  while(merge(y, beta)){update(beta, tks, y, false, rho0);}
610  split(beta, tks, y);
611  beta=beta/coolingFactor_;
612  }else{
613  beta=beta/coolingFactor_;
614  splitAll(y);
615  }
616 
617  // make sure we are not too far from equilibrium before cooling further
618  niter = 0;
619  while ((update(beta, tks, y, false, rho0) > 1.e-6) &&
620  (niter++ < maxIterations_)) {}
621 
622  if(verbose_){ dump( beta, y, tks, 0); }
623  }
624 
625 
626  if(useTc_){
627  //last round of splitting, make sure no critical clusters are left
628  if(verbose_){ std::cout << "last spliting at " << 1./beta << std::endl; }
629  update(beta, tks,y, false, rho0);// make sure Tc is up-to-date
630  while(merge(y,beta)){update(beta, tks,y, false, rho0);}
631  unsigned int ntry=0;
632  double threshold = 1.0;
633  while( split(beta, tks, y, threshold) && (ntry++<10) ){
634  niter=0;
635  while((update(beta, tks,y, false, rho0)>1.e-6) && (niter++ < maxIterations_)){}
636  while(merge(y,beta)){update(beta, tks,y, false, rho0);}
637  if(verbose_){
638  std::cout << "after final splitting, try " << ntry << std::endl;
639  dump(beta, y, tks, 2);
640  }
641  // relax splitting a bit to reduce multiple split-merge cycles of the same cluster
642  threshold *= 1.1;
643  }
644  }else{
645  // merge collapsed clusters
646  while(merge(y,beta)){update(beta, tks,y, false, rho0);}
647  }
648 
649  if (verbose_) {
650  update(beta, tks,y, false, rho0);
651  std::cout << "dump after 1st merging " << endl;
652  dump(beta, y, tks, 2);
653  }
654 
655 
656  // switch on outlier rejection at T=Tmin
657  if(dzCutOff_ > 0){
658  rho0 = 1./nt;
659  for(unsigned int a=0; a<10; a++){ update(beta, tks, y, true, a*rho0/10);} // adiabatic turn-on
660  }
661 
662  niter=0;
663  while ((update(beta, tks, y, true, rho0) > 1.e-8) && (niter++ < maxIterations_)) {};
664  if (verbose_) {
665  std::cout << "dump after noise-suppression, rho0=" << rho0 << " niter = " << niter << endl;
666  dump(beta, y, tks, 2);
667  }
668 
669  // merge again (some cluster split by outliers collapse here)
670  while (merge(y, beta)) {update(beta, tks, y, true, rho0); }
671  if (verbose_) {
672  std::cout << "dump after merging " << endl;
673  dump(beta, y, tks, 2);
674  }
675 
676  // go down to the purging temperature (if it is lower than tmin)
677  while( beta < betapurge_ ){
678  beta = min( beta/coolingFactor_, betapurge_);
679  niter = 0;
680  while ((update(beta, tks, y, false, rho0) > 1.e-8) && (niter++ < maxIterations_)) {}
681  }
682 
683 
684  // eliminate insigificant vertices, this is more restrictive at higher T
685  while (purge(y, tks, rho0, beta)) {
686  niter = 0;
687  while ((update(beta, tks, y, true, rho0) > 1.e-6) && (niter++ < maxIterations_)) {}
688  }
689 
690  if (verbose_) {
691  update(beta, tks,y, true, rho0);
692  std::cout << " after purging " << std:: endl;
693  dump(beta, y, tks, 2);
694  }
695 
696  // optionally cool some more without doing anything, to make the assignment harder
697  while( beta < betastop_ ){
698  beta = min( beta/coolingFactor_, betastop_);
699  niter =0;
700  while ((update(beta, tks, y, true, rho0) > 1.e-8) && (niter++ < maxIterations_)) {}
701  }
702 
703 
704  if (verbose_) {
705  std::cout << "Final result, rho0=" << std::scientific << rho0 << endl;
706  dump(beta, y, tks, 2);
707  }
708 
709  // select significant tracks and use a TransientVertex as a container
710  GlobalError dummyError(0.01, 0, 0.01, 0., 0., 0.01);
711 
712  // ensure correct normalization of probabilities, should makes double assignment reasonably impossible
713  const unsigned int nv = y.GetSize();
714  for (unsigned int k = 0; k < nv; k++)
715  if ( edm::isNotFinite(y._pk[k]) || edm::isNotFinite(y._z[k]) ) { y._pk[k]=0; y._z[k]=0;}
716 
717  for (unsigned int i = 0; i < nt; i++) // initialize
718  tks._Z_sum[i] = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
719 
720  // improve vectorization (does not require reduction ....)
721  for (unsigned int k = 0; k < nv; k++) {
722  for (unsigned int i = 0; i < nt; i++)
723  tks._Z_sum[i] += y._pk[k] * local_exp(-beta * Eik(tks._z[i], y._z[k],tks._dz2[i]));
724  }
725 
726 
727  for (unsigned int k = 0; k < nv; k++) {
728  GlobalPoint pos(0, 0, y._z[k]);
729 
730  vector<reco::TransientTrack> vertexTracks;
731  for (unsigned int i = 0; i < nt; i++) {
732  if (tks._Z_sum[i] > 1e-100) {
733 
734  double p = y._pk[k] * local_exp(-beta * Eik(tks._z[i], y._z[k],
735  tks._dz2[i])) / tks._Z_sum[i];
736  if ((tks._pi[i] > 0) && (p > mintrkweight_)) {
737  vertexTracks.push_back(*(tks.tt[i]));
738  tks._Z_sum[i] = 0; // setting Z=0 excludes double assignment
739  }
740  }
741  }
742  TransientVertex v(pos, dummyError, vertexTracks, 0);
743  clusters.push_back(v);
744  }
745 
746  return clusters;
747 
748 }
749 
750 vector<vector<reco::TransientTrack> > DAClusterizerInZ_vect::clusterize(
751  const vector<reco::TransientTrack> & tracks) const {
752 
753  if (verbose_) {
754  std::cout << "###################################################" << endl;
755  std::cout << "# vectorized DAClusterizerInZ_vect::clusterize nt=" << tracks.size() << endl;
756  std::cout << "###################################################" << endl;
757  }
758 
759  vector<vector<reco::TransientTrack> > clusters;
760  vector<TransientVertex> && pv = vertices(tracks);
761 
762  if (verbose_) {
763  std::cout << "# DAClusterizerInZ::clusterize pv.size=" << pv.size()
764  << endl;
765  }
766  if (pv.empty()) {
767  return clusters;
768  }
769 
770  // fill into clusters and merge
771  vector<reco::TransientTrack> aCluster = pv.begin()->originalTracks();
772 
773  for (auto k = pv.begin() + 1; k != pv.end(); k++) {
774  if ( std::abs(k->position().z() - (k - 1)->position().z()) > (2 * vertexSize_)) {
775  // close a cluster
776  if (aCluster.size()>1){
777  clusters.push_back(aCluster);
778  }else{
779  if(verbose_){
780  std::cout << " one track cluster at " << k->position().z() << " suppressed" << std::endl;
781  }
782  }
783  aCluster.clear();
784  }
785  for (unsigned int i = 0; i < k->originalTracks().size(); i++) {
786  aCluster.push_back(k->originalTracks()[i]);
787  }
788 
789  }
790  clusters.emplace_back(std::move(aCluster));
791 
792  return clusters;
793 
794 }
795 
796 
797 
798 void DAClusterizerInZ_vect::dump(const double beta, const vertex_t & y,
799  const track_t & tks, int verbosity) const {
800 
801  const unsigned int nv = y.GetSize();
802  const unsigned int nt = tks.GetSize();
803 
804  std::vector< unsigned int > iz;
805  for(unsigned int j=0; j<nt; j++){ iz.push_back(j); }
806  std::sort(iz.begin(), iz.end(), [tks](unsigned int a, unsigned int b){ return tks._z[a]<tks._z[b];} );
807  std::cout << std::endl;
808  std::cout << "-----DAClusterizerInZ::dump ----" << nv << " clusters " << std::endl;
809  std::cout << " z= ";
810  std::cout << setprecision(4);
811  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
812  if (std::fabs(y._z[ivertex]-zdumpcenter_) < zdumpwidth_){
813  std::cout << setw(8) << fixed << y._z[ivertex];
814  }
815  }
816  std::cout << endl << "T=" << setw(15) << 1. / beta
817  << " Tmin =" << setw(10) << 1./betamax_
818  << " Tc= ";
819  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
820  if (std::fabs(y._z[ivertex]-zdumpcenter_) < zdumpwidth_){
821  double Tc = 2*y._swE[ivertex]/y._sw[ivertex];
822  std::cout << setw(8) << fixed << setprecision(1) << Tc;
823  }
824  }
825  std::cout << endl;
826 
827  std::cout << " pk= ";
828  double sumpk = 0;
829  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
830  sumpk += y._pk[ivertex];
831  if (std::fabs(y._z[ivertex] - zdumpcenter_) > zdumpwidth_) continue;
832  std::cout << setw(8) << setprecision(4) << fixed << y._pk[ivertex];
833  }
834  std::cout << endl;
835 
836  std::cout << " nt= ";
837  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
838  sumpk += y._pk[ivertex];
839  if (std::fabs(y._z[ivertex] - zdumpcenter_) > zdumpwidth_) continue;
840  std::cout << setw(8) << setprecision(1) << fixed << y._pk[ivertex]*nt;
841  }
842  std::cout << endl;
843 
844  if (verbosity > 0) {
845  double E = 0, F = 0;
846  std::cout << endl;
847  std::cout
848  << "---- z +/- dz ip +/-dip pt phi eta weights ----"
849  << endl;
850  std::cout << setprecision(4);
851  for (unsigned int i0 = 0; i0 < nt; i0++) {
852  unsigned int i = iz[i0];
853  if (tks._Z_sum[i] > 0) {
854  F -= std::log(tks._Z_sum[i]) / beta;
855  }
856  double tz = tks._z[i];
857 
858  if( std::fabs(tz - zdumpcenter_) > zdumpwidth_) continue;
859  std::cout << setw(4) << i << ")" << setw(8) << fixed << setprecision(4)
860  << tz << " +/-" << setw(6) << sqrt(1./tks._dz2[i]);
861 
862  if (tks.tt[i]->track().quality(reco::TrackBase::highPurity)) {
863  std::cout << " *";
864  } else {
865  std::cout << " ";
866  }
867  if (tks.tt[i]->track().hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1)) {
868  std::cout << "+";
869  } else {
870  std::cout << "-";
871  }
872  std::cout << setw(1)
873  << tks.tt[i]->track().hitPattern().pixelBarrelLayersWithMeasurement(); // see DataFormats/TrackReco/interface/HitPattern.h
874  std::cout << setw(1)
875  << tks.tt[i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
876  std::cout << setw(1) << hex
877  << tks.tt[i]->track().hitPattern().trackerLayersWithMeasurement()
878  - tks.tt[i]->track().hitPattern().pixelLayersWithMeasurement()
879  << dec;
880  std::cout << "=" << setw(1) << hex
881  << tks.tt[i]->track().hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS)
882  << dec;
883 
884  Measurement1D IP =
885  tks.tt[i]->stateAtBeamLine().transverseImpactParameter();
886  std::cout << setw(8) << IP.value() << "+/-" << setw(6) << IP.error();
887  std::cout << " " << setw(6) << setprecision(2)
888  << tks.tt[i]->track().pt() * tks.tt[i]->track().charge();
889  std::cout << " " << setw(5) << setprecision(2)
890  << tks.tt[i]->track().phi() << " " << setw(5)
891  << setprecision(2) << tks.tt[i]->track().eta();
892 
893  double sump = 0.;
894  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
895  if (std::fabs(y._z[ivertex]-zdumpcenter_) > zdumpwidth_) continue;
896 
897  if ((tks._pi[i] > 0) && (tks._Z_sum[i] > 0)) {
898  //double p=pik(beta,tks[i],*k);
899  double p = y._pk[ivertex] * exp(-beta * Eik(tks._z[i], y._z[ivertex], tks._dz2[i])) / tks._Z_sum[i];
900  if (p > 0.0001) {
901  std::cout << setw(8) << setprecision(3) << p;
902  } else {
903  std::cout << " . ";
904  }
905  E += p * Eik(tks._z[i], y._z[ivertex], tks._dz2[i]);
906  sump += p;
907  } else {
908  std::cout << " ";
909  }
910  }
911  std::cout << endl;
912  }
913  std::cout << endl << "T=" << 1 / beta << " E=" << E << " n=" << y.GetSize()
914  << " F= " << F << endl << "----------" << endl;
915  }
916 }
void AddItem(double new_z, double new_dz2, const reco::TransientTrack *new_tt, double new_pi)
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
const double w
Definition: UKUtility.cc:23
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const
void splitAll(vertex_t &y) const
void AddItem(double new_z, double new_pk)
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
double error() const
Definition: Measurement1D.h:27
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
std::vector< const reco::TransientTrack * > tt
A arg
Definition: Factorize.h:38
U second(std::pair< T, U > const &p)
T sqrt(T t)
Definition: SSEVec.h:18
void InsertItem(unsigned int i, double new_z, double new_pk)
def pv(vc)
Definition: MetAnalyzer.py:7
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const override
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
T min(T a, T b)
Definition: MathUtil.h:58
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
int nt
Definition: AMPTWrapper.h:32
int k[5][pyjets_maxn]
double update(double beta, track_t &gtracks, vertex_t &gvertices, bool useRho0, const double &rho0) const
bool merge(vertex_t &y, double &beta) const
DAClusterizerInZ_vect(const edm::ParameterSet &conf)
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0) const
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:88
double b
Definition: hdecay.h:120
double value() const
Definition: Measurement1D.h:25
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks, const int verbosity=0) const
double p1[4]
Definition: TauolaWrapper.h:89
#define update(a, b)
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:509
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
bool purge(vertex_t &, 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
def move(src, dest)
Definition: eostools.py:511
def merge(dictlist, TELL=False)
Definition: MatrixUtil.py:194