CMS 3D CMS Logo

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