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 
15 //#define DEBUG
16 #ifdef DEBUG
17 #define DEBUGLEVEL 0
18 #endif
19 
21  // hardcoded parameters
22  maxIterations_ = 1000;
23  mintrkweight_ = 0.5;
24 
25  // configurable debug output
26 #ifdef DEBUG
27  zdumpcenter_ = conf.getUntrackedParameter<double>("zdumpcenter", 0.);
28  zdumpwidth_ = conf.getUntrackedParameter<double>("zdumpwidth", 20.);
29 #endif
30 
31  // configurable parameters
32  double Tmin = conf.getParameter<double>("Tmin");
33  double Tpurge = conf.getParameter<double>("Tpurge");
34  double Tstop = conf.getParameter<double>("Tstop");
35  vertexSize_ = conf.getParameter<double>("vertexSize");
36  coolingFactor_ = conf.getParameter<double>("coolingFactor");
37  d0CutOff_ = conf.getParameter<double>("d0CutOff");
38  dzCutOff_ = conf.getParameter<double>("dzCutOff");
39  uniquetrkweight_ = conf.getParameter<double>("uniquetrkweight");
40  uniquetrkminp_ = conf.getParameter<double>("uniquetrkminp");
41  zmerge_ = conf.getParameter<double>("zmerge");
42  sel_zrange_ = conf.getParameter<double>("zrange");
43  convergence_mode_ = conf.getParameter<int>("convergence_mode");
44  delta_lowT_ = conf.getParameter<double>("delta_lowT");
45  delta_highT_ = conf.getParameter<double>("delta_highT");
46  runInBlocks_ = conf.getParameter<bool>("runInBlocks");
47  block_size_ = conf.getParameter<unsigned int>("block_size");
48  overlap_frac_ = conf.getParameter<double>("overlap_frac");
49 
50 #ifdef DEBUG
51  std::cout << "DAClusterizerinZ_vect: mintrkweight = " << mintrkweight_ << std::endl;
52  std::cout << "DAClusterizerinZ_vect: uniquetrkweight = " << uniquetrkweight_ << std::endl;
53  std::cout << "DAClusterizerInZ_vect: uniquetrkminp = " << uniquetrkminp_ << std::endl;
54  std::cout << "DAClusterizerinZ_vect: zmerge = " << zmerge_ << std::endl;
55  std::cout << "DAClusterizerinZ_vect: Tmin = " << Tmin << std::endl;
56  std::cout << "DAClusterizerinZ_vect: Tpurge = " << Tpurge << std::endl;
57  std::cout << "DAClusterizerinZ_vect: Tstop = " << Tstop << std::endl;
58  std::cout << "DAClusterizerinZ_vect: vertexSize = " << vertexSize_ << std::endl;
59  std::cout << "DAClusterizerinZ_vect: coolingFactor = " << coolingFactor_ << std::endl;
60  std::cout << "DAClusterizerinZ_vect: d0CutOff = " << d0CutOff_ << std::endl;
61  std::cout << "DAClusterizerinZ_vect: dzCutOff = " << dzCutOff_ << std::endl;
62  std::cout << "DAClusterizerInZ_vect: zrange = " << sel_zrange_ << std::endl;
63  std::cout << "DAClusterizerinZ_vect: convergence mode = " << convergence_mode_ << std::endl;
64  std::cout << "DAClusterizerinZ_vect: delta_highT = " << delta_highT_ << std::endl;
65  std::cout << "DAClusterizerinZ_vect: delta_lowT = " << delta_lowT_ << std::endl;
66 
67  std::cout << "DAClusterizerinZ_vect: run in blocks = " << runInBlocks_ << std::endl;
68  std::cout << "DAClusterizerinZ_vect: block_size = " << block_size_ << std::endl;
69  std::cout << "DAClusterizerinZ_vect: overlap_fraction = " << overlap_frac_ << std::endl;
70  std::cout << "DAClusterizerinZ_vect: DEBUGLEVEL " << DEBUGLEVEL << std::endl;
71 #endif
72 
73  if (convergence_mode_ > 1) {
74  edm::LogWarning("DAClusterizerinZ_vect")
75  << "DAClusterizerInZ_vect: invalid convergence_mode " << convergence_mode_ << " reset to default " << 0;
76  convergence_mode_ = 0;
77  }
78 
79  if (Tmin == 0) {
80  betamax_ = 1.0;
81  edm::LogWarning("DAClusterizerinZ_vect")
82  << "DAClusterizerInZ_vect: invalid Tmin " << Tmin << " reset to default " << 1. / betamax_;
83  } else {
84  betamax_ = 1. / Tmin;
85  }
86 
87  if ((Tpurge > Tmin) || (Tpurge == 0)) {
88  edm::LogWarning("DAClusterizerinZ_vect")
89  << "DAClusterizerInZ_vect: invalid Tpurge " << Tpurge << " set to " << Tmin;
90  Tpurge = Tmin;
91  }
92  betapurge_ = 1. / Tpurge;
93 
94  if ((Tstop > Tpurge) || (Tstop == 0)) {
95  edm::LogWarning("DAClusterizerinZ_vect")
96  << "DAClusterizerInZ_vect: invalid Tstop " << Tstop << " set to " << max(1., Tpurge);
97  Tstop = max(1., Tpurge);
98  }
99  betastop_ = 1. / Tstop;
100 }
101 
102 namespace {
103  inline double local_exp(double const& inp) { return vdt::fast_exp(inp); }
104 
105  inline void local_exp_list_range(double const* __restrict__ arg_inp,
106  double* __restrict__ arg_out,
107  const int kmin,
108  const int kmax) {
109  for (auto i = kmin; i != kmax; ++i)
110  arg_out[i] = vdt::fast_exp(arg_inp[i]);
111  }
112 
113 } // namespace
114 
115 void DAClusterizerInZ_vect::verify(const vertex_t& v, const track_t& tks, unsigned int nv, unsigned int nt) const {
116  if (!(nv == 999999)) {
117  assert(nv == v.getSize());
118  } else {
119  nv = v.getSize();
120  }
121 
122  if (!(nt == 999999)) {
123  assert(nt == tks.getSize());
124  } else {
125  nt = tks.getSize();
126  }
127 
128  assert(v.zvtx_vec.size() == nv);
129  assert(v.rho_vec.size() == nv);
130  assert(v.swz_vec.size() == nv);
131  assert(v.exp_arg_vec.size() == nv);
132  assert(v.exp_vec.size() == nv);
133  assert(v.se_vec.size() == nv);
134  assert(v.swz_vec.size() == nv);
135  assert(v.swE_vec.size() == nv);
136 
137  assert(v.zvtx == &v.zvtx_vec.front());
138  assert(v.rho == &v.rho_vec.front());
139  assert(v.exp_arg == &v.exp_arg_vec.front());
140  assert(v.sw == &v.sw_vec.front());
141  assert(v.swz == &v.swz_vec.front());
142  assert(v.se == &v.se_vec.front());
143  assert(v.swE == &v.swE_vec.front());
144 
145  for (unsigned int k = 0; k < nv - 1; k++) {
146  if (v.zvtx_vec[k] <= v.zvtx_vec[k + 1])
147  continue;
148  cout << " Z, cluster z-ordering assertion failure z[" << k << "] =" << v.zvtx_vec[k] << " z[" << k + 1
149  << "] =" << v.zvtx_vec[k + 1] << endl;
150  }
151 
152  assert(nt == tks.zpca_vec.size());
153  assert(nt == tks.dz2_vec.size());
154  assert(nt == tks.tt.size());
155  assert(nt == tks.tkwt_vec.size());
156  assert(nt == tks.sum_Z_vec.size());
157  assert(nt == tks.kmin.size());
158  assert(nt == tks.kmax.size());
159 
160  assert(tks.zpca == &tks.zpca_vec.front());
161  assert(tks.dz2 == &tks.dz2_vec.front());
162  assert(tks.tkwt == &tks.tkwt_vec.front());
163  assert(tks.sum_Z == &tks.sum_Z_vec.front());
164 
165  for (unsigned int i = 0; i < nt; i++) {
166  if ((tks.kmin[i] < tks.kmax[i]) && (tks.kmax[i] <= nv))
167  continue;
168  cout << "track vertex range assertion failure" << i << "/" << nt << " kmin,kmax=" << tks.kmin[i] << ", "
169  << tks.kmax[i] << " nv=" << nv << endl;
170  }
171 
172  for (unsigned int i = 0; i < nt; i++) {
173  assert((tks.kmin[i] < tks.kmax[i]) && (tks.kmax[i] <= nv));
174  }
175 }
176 
177 DAClusterizerInZ_vect::track_t DAClusterizerInZ_vect::fill(const vector<reco::TransientTrack>& tracks) const {
178  // prepare track data for clustering
179  track_t tks;
180  double sumtkwt = 0.;
181  for (auto it = tracks.begin(); it != tracks.end(); it++) {
182  if (!(*it).isValid())
183  continue;
184  double t_tkwt = 1.;
185  double t_z = ((*it).stateAtBeamLine().trackStateAtPCA()).position().z();
186  if (std::fabs(t_z) > 1000.)
187  continue;
188  auto const& t_mom = (*it).stateAtBeamLine().trackStateAtPCA().momentum();
189  // get the beam-spot
190  reco::BeamSpot beamspot = (it->stateAtBeamLine()).beamSpot();
191  double t_dz2 = std::pow((*it).track().dzError(), 2) // track errror
192  + (std::pow(beamspot.BeamWidthX() * t_mom.x(), 2) + std::pow(beamspot.BeamWidthY() * t_mom.y(), 2)) *
193  std::pow(t_mom.z(), 2) / std::pow(t_mom.perp2(), 2) // beam spot width
194  + std::pow(vertexSize_, 2); // intrinsic vertex size, safer for outliers and short lived decays
195  t_dz2 = 1. / t_dz2;
196  if (edm::isNotFinite(t_dz2) || t_dz2 < std::numeric_limits<double>::min())
197  continue;
198  if (d0CutOff_ > 0) {
199  Measurement1D atIP = (*it).stateAtBeamLine().transverseImpactParameter(); // error contains beamspot
200  t_tkwt = 1. / (1. + local_exp(std::pow(atIP.value() / atIP.error(), 2) -
201  std::pow(d0CutOff_, 2))); // reduce weight for high ip tracks
203  edm::LogWarning("DAClusterizerinZ_vect") << "rejected track t_tkwt " << t_tkwt;
204  continue; // usually is > 0.99
205  }
206  }
207  tks.addItemSorted(t_z, t_dz2, &(*it), t_tkwt);
208  sumtkwt += t_tkwt;
209  }
210 
211  if (sumtkwt > 0) {
212  tks.extractRaw();
213  tks.osumtkwt = 1. / sumtkwt;
214  } else {
215  tks.osumtkwt = 0.;
216  }
217 
218 #ifdef DEBUG
219  if (DEBUGLEVEL > 0) {
220  std::cout << "Track count (Z) " << tks.getSize() << std::endl;
221  }
222 #endif
223  return tks;
224 }
225 
226 namespace {
227  inline double Eik(double t_z, double k_z, double t_dz2) { return std::pow(t_z - k_z, 2) * t_dz2; }
228 } // namespace
229 
230 void DAClusterizerInZ_vect::set_vtx_range(double beta, track_t& gtracks, vertex_t& gvertices) const {
231  const unsigned int nv = gvertices.getSize();
232  const unsigned int nt = gtracks.getSize();
233 
234  if (nv == 0) {
235  edm::LogWarning("DAClusterizerinZ_vect") << "empty cluster list in set_vtx_range";
236  return;
237  }
238 
239  for (auto itrack = 0U; itrack < nt; ++itrack) {
240  double zrange = max(sel_zrange_ / sqrt(beta * gtracks.dz2[itrack]), zrange_min_);
241 
242  double zmin = gtracks.zpca[itrack] - zrange;
243  unsigned int kmin = min(nv - 1, gtracks.kmin[itrack]);
244  // find the smallest vertex_z that is larger than zmin
245  if (gvertices.zvtx[kmin] > zmin) {
246  while ((kmin > 0) && (gvertices.zvtx[kmin - 1] > zmin)) {
247  kmin--;
248  }
249  } else {
250  while ((kmin < (nv - 1)) && (gvertices.zvtx[kmin] < zmin)) {
251  kmin++;
252  }
253  }
254 
255  double zmax = gtracks.zpca[itrack] + zrange;
256  unsigned int kmax = min(nv - 1, gtracks.kmax[itrack] - 1);
257  // note: kmax points to the last vertex in the range, while gtracks.kmax points to the entry BEHIND the last vertex
258  // find the largest vertex_z that is smaller than zmax
259  if (gvertices.zvtx[kmax] < zmax) {
260  while ((kmax < (nv - 1)) && (gvertices.zvtx[kmax + 1] < zmax)) {
261  kmax++;
262  }
263  } else {
264  while ((kmax > 0) && (gvertices.zvtx[kmax] > zmax)) {
265  kmax--;
266  }
267  }
268 
269  if (kmin <= kmax) {
270  gtracks.kmin[itrack] = kmin;
271  gtracks.kmax[itrack] = kmax + 1;
272  } else {
273  gtracks.kmin[itrack] = max(0U, min(kmin, kmax));
274  gtracks.kmax[itrack] = min(nv, max(kmin, kmax) + 1);
275  }
276  }
277 }
278 
279 void DAClusterizerInZ_vect::clear_vtx_range(track_t& gtracks, vertex_t& gvertices) const {
280  const unsigned int nt = gtracks.getSize();
281  const unsigned int nv = gvertices.getSize();
282  for (auto itrack = 0U; itrack < nt; ++itrack) {
283  gtracks.kmin[itrack] = 0;
284  gtracks.kmax[itrack] = nv;
285  }
286 }
287 
289  double beta, track_t& gtracks, vertex_t& gvertices, const double rho0, const bool updateTc) const {
290  // update weights and vertex positions
291  // returns the maximum of changes of vertex positions
292  // sums needed for Tc are only updated if updateTC == true
293 
294  const unsigned int nt = gtracks.getSize();
295  const unsigned int nv = gvertices.getSize();
296  auto osumtkwt = gtracks.osumtkwt;
297 
298  double Z_init = 0;
299  if (rho0 > 0) {
300  Z_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_); // cut-off
301  }
302 
303  // define kernels
304  auto kernel_calc_exp_arg_range = [beta](const unsigned int itrack,
305  track_t const& tracks,
306  vertex_t const& vertices,
307  const unsigned int kmin,
308  const unsigned int kmax) {
309  const double track_z = tracks.zpca[itrack];
310  const double botrack_dz2 = -beta * tracks.dz2[itrack];
311 
312  // auto-vectorized
313  for (unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
314  auto mult_res = track_z - vertices.zvtx[ivertex];
315  vertices.exp_arg[ivertex] = botrack_dz2 * (mult_res * mult_res);
316  }
317  };
318 
319  auto kernel_add_Z_range = [Z_init](
320  vertex_t const& vertices, const unsigned int kmin, const unsigned int kmax) -> double {
321  double ZTemp = Z_init;
322  for (unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
323  ZTemp += vertices.rho[ivertex] * vertices.exp[ivertex];
324  }
325  return ZTemp;
326  };
327 
328  auto kernel_calc_normalization_range = [updateTc](const unsigned int track_num,
329  track_t& tracks,
331  const unsigned int kmin,
332  const unsigned int kmax) {
333  auto o_trk_sum_Z = tracks.tkwt[track_num] / tracks.sum_Z[track_num];
334  auto o_trk_dz2 = tracks.dz2[track_num];
335  auto tmp_trk_z = tracks.zpca[track_num];
336 
337  // auto-vectorized
338  if (updateTc) {
339 #pragma GCC ivdep
340  for (unsigned int k = kmin; k < kmax; ++k) {
341  vertices.se[k] += vertices.exp[k] * o_trk_sum_Z;
342  auto w = vertices.rho[k] * vertices.exp[k] * (o_trk_sum_Z * o_trk_dz2);
343  vertices.sw[k] += w;
344  vertices.swz[k] += w * tmp_trk_z;
345  vertices.swE[k] += w * vertices.exp_arg[k];
346  }
347  } else {
348  // same loop but without updating sWE
349 #pragma GCC ivdep
350  for (unsigned int k = kmin; k < kmax; ++k) {
351  vertices.se[k] += vertices.exp[k] * o_trk_sum_Z;
352  auto w = vertices.rho[k] * vertices.exp[k] * (o_trk_sum_Z * o_trk_dz2);
353  vertices.sw[k] += w;
354  vertices.swz[k] += w * tmp_trk_z;
355  }
356  }
357  };
358 
359  if (updateTc) {
360  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
361  gvertices.se[ivertex] = 0.0;
362  gvertices.sw[ivertex] = 0.0;
363  gvertices.swz[ivertex] = 0.0;
364  gvertices.swE[ivertex] = 0.0;
365  }
366  } else {
367  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
368  gvertices.se[ivertex] = 0.0;
369  gvertices.sw[ivertex] = 0.0;
370  gvertices.swz[ivertex] = 0.0;
371  }
372  }
373 
374  // loop over tracks
375  for (auto itrack = 0U; itrack < nt; ++itrack) {
376  const unsigned int kmin = gtracks.kmin[itrack];
377  const unsigned int kmax = gtracks.kmax[itrack];
378 
379  kernel_calc_exp_arg_range(itrack, gtracks, gvertices, kmin, kmax);
380  local_exp_list_range(gvertices.exp_arg, gvertices.exp, kmin, kmax);
381  gtracks.sum_Z[itrack] = kernel_add_Z_range(gvertices, kmin, kmax);
382 
383  if (edm::isNotFinite(gtracks.sum_Z[itrack]))
384  gtracks.sum_Z[itrack] = 0.0;
385 
386  if (gtracks.sum_Z[itrack] > 1.e-100) {
387  kernel_calc_normalization_range(itrack, gtracks, gvertices, kmin, kmax);
388  }
389  }
390 
391  // (un-)apply the factor -beta which is needed in exp_arg, but not in swE
392  if (updateTc) {
393  auto obeta = -1. / beta;
394  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
395  gvertices.swE[ivertex] *= obeta;
396  }
397  }
398 
399  // now update z and rho
400  auto kernel_calc_z = [osumtkwt, nv](vertex_t& vertices) -> double {
401  double delta = 0;
402  // does not vectorize
403  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
404  if (vertices.sw[ivertex] > 0) {
405  auto znew = vertices.swz[ivertex] / vertices.sw[ivertex];
406  // prevents from vectorizing if
407  delta = max(std::abs(vertices.zvtx[ivertex] - znew), delta);
408  vertices.zvtx[ivertex] = znew;
409  }
410  }
411 
412  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex)
413  vertices.rho[ivertex] = vertices.rho[ivertex] * vertices.se[ivertex] * osumtkwt;
414 
415  return delta;
416  };
417 
418  double delta = kernel_calc_z(gvertices);
419 
420  // return how much the prototypes moved
421  return delta;
422 }
423 
425  double beta, track_t& tks, vertex_t& v, const double delta_max0, const double rho0) const {
426  unsigned int niter = 0;
427  double delta = 0;
428  double delta_max = delta_lowT_;
429 
430  if (convergence_mode_ == 0) {
431  delta_max = delta_max0;
432  } else if (convergence_mode_ == 1) {
433  delta_max = delta_lowT_ / sqrt(std::max(beta, 1.0));
434  }
435 
436  set_vtx_range(beta, tks, v);
437  double delta_sum_range = 0; // accumulate max(|delta-z|) as a lower bound
438  std::vector<double> z0 = v.zvtx_vec;
439 
440  while (niter++ < maxIterations_) {
441  delta = update(beta, tks, v, rho0);
442  delta_sum_range += delta;
443 
444  if (delta_sum_range > zrange_min_) {
445  for (unsigned int k = 0; k < v.getSize(); k++) {
446  if (std::abs(v.zvtx_vec[k] - z0[k]) > zrange_min_) {
447  set_vtx_range(beta, tks, v);
448  delta_sum_range = 0;
449  z0 = v.zvtx_vec;
450  break;
451  }
452  }
453  }
454 
455  if (delta < delta_max) {
456  break;
457  }
458  }
459 
460 #ifdef DEBUG
461  if (DEBUGLEVEL > 0) {
462  std::cout << "DAClusterizerInZ_vect.thermalize niter = " << niter << " at T = " << 1 / beta
463  << " nv = " << v.getSize() << std::endl;
464  if (DEBUGLEVEL > 2)
465  dump(beta, v, tks, 0, rho0);
466  }
467 #endif
468 
469  return niter;
470 }
471 
472 bool DAClusterizerInZ_vect::merge(vertex_t& y, track_t& tks, double& beta) const {
473  // merge clusters that collapsed or never separated,
474  // only merge if the estimated critical temperature of the merged vertex is below the current temperature
475  // return true if vertices were merged, false otherwise
476  const unsigned int nv = y.getSize();
477 
478  if (nv < 2)
479  return false;
480 
481  // merge the smallest distance clusters first
482  std::vector<std::pair<double, unsigned int>> critical;
483  for (unsigned int k = 0; (k + 1) < nv; k++) {
484  if (std::fabs(y.zvtx[k + 1] - y.zvtx[k]) < zmerge_) {
485  critical.push_back(make_pair(std::fabs(y.zvtx[k + 1] - y.zvtx[k]), k));
486  }
487  }
488  if (critical.empty())
489  return false;
490 
491  std::stable_sort(critical.begin(), critical.end(), std::less<std::pair<double, unsigned int>>());
492 
493  for (unsigned int ik = 0; ik < critical.size(); ik++) {
494  unsigned int k = critical[ik].second;
495  double rho = y.rho[k] + y.rho[k + 1];
496 
497 #ifdef DEBUG
498  assert((k + 1) < nv);
499  if (DEBUGLEVEL > 1) {
500  std::cout << "merging " << fixed << setprecision(4) << y.zvtx[k + 1] << " and " << y.zvtx[k]
501  << " sw = " << y.sw[k] + y.sw[k + 1] << std::endl;
502  }
503 #endif
504 
505  if (rho > 0) {
506  y.zvtx[k] = (y.rho[k] * y.zvtx[k] + y.rho[k + 1] * y.zvtx[k + 1]) / rho;
507  if (not edm::isFinite(y.zvtx[k]))
508  y.zvtx[k] = 0.5 * (y.zvtx[k] + y.zvtx[k + 1]);
509  } else {
510  y.zvtx[k] = 0.5 * (y.zvtx[k] + y.zvtx[k + 1]);
511  }
512  y.rho[k] = rho;
513  y.sw[k] += y.sw[k + 1];
514  y.removeItem(k + 1, tks);
515  set_vtx_range(beta, tks, y);
516  y.extractRaw();
517  return true;
518  }
519 
520  return false;
521 }
522 
523 bool DAClusterizerInZ_vect::purge(vertex_t& y, track_t& tks, double& rho0, const double beta) const {
524  constexpr double eps = 1.e-100;
525  // eliminate clusters with only one significant/unique track
526  const unsigned int nv = y.getSize();
527  const unsigned int nt = tks.getSize();
528 
529  if (nv < 2)
530  return false;
531 
532  std::vector<double> sump_v(nv), arg_cache_v(nv), exp_cache_v(nv), pcut_cache_v(nv);
533  std::vector<int> nUnique_v(nv);
534  double* __restrict__ parg_cache;
535  double* __restrict__ pexp_cache;
536  double* __restrict__ ppcut_cache;
537  double* __restrict__ psump;
538  int* __restrict__ pnUnique;
539  int constexpr nunique_min_ = 2;
540 
541  set_vtx_range(beta, tks, y);
542 
543  parg_cache = arg_cache_v.data();
544  pexp_cache = exp_cache_v.data();
545  ppcut_cache = pcut_cache_v.data();
546  psump = sump_v.data();
547  pnUnique = nUnique_v.data();
548 
549  const auto rhoconst = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
550  for (unsigned int k = 0; k < nv; k++) {
551  const double pmax = y.rho[k] / (y.rho[k] + rhoconst);
552  ppcut_cache[k] = uniquetrkweight_ * pmax;
553  }
554 
555  for (unsigned int i = 0; i < nt; i++) {
556  const auto invZ = ((tks.sum_Z[i] > eps) && (tks.tkwt[i] > uniquetrkminp_)) ? 1. / tks.sum_Z[i] : 0.;
557  const auto track_z = tks.zpca[i];
558  const auto botrack_dz2 = -beta * tks.dz2[i];
559  const auto kmin = tks.kmin[i];
560  const auto kmax = tks.kmax[i];
561 
562  for (unsigned int k = kmin; k < kmax; k++) {
563  const auto mult_resz = track_z - y.zvtx[k];
564  parg_cache[k] = botrack_dz2 * (mult_resz * mult_resz);
565  }
566 
567  local_exp_list_range(parg_cache, pexp_cache, kmin, kmax);
568 
569  for (unsigned int k = kmin; k < kmax; k++) {
570  const double p = y.rho[k] * pexp_cache[k] * invZ;
571  psump[k] += p;
572  pnUnique[k] += (p > ppcut_cache[k]) ? 1 : 0;
573  }
574  }
575 
576  double sumpmin = nt;
577  unsigned int k0 = nv;
578  for (unsigned k = 0; k < nv; k++) {
579  if ((pnUnique[k] < nunique_min_) && (psump[k] < sumpmin)) {
580  sumpmin = psump[k];
581  k0 = k;
582  }
583  }
584 
585  if (k0 != nv) {
586 #ifdef DEBUG
587  assert(k0 < y.getSize());
588  if (DEBUGLEVEL > 1) {
589  std::cout << "eliminating prototype at " << std::setw(10) << std::setprecision(4) << y.zvtx[k0]
590  << " with sump=" << sumpmin << " rho*nt =" << y.rho[k0] * nt << " pnUnique=" << pnUnique[k0] << endl;
591  }
592 #endif
593 
594  y.removeItem(k0, tks);
595  set_vtx_range(beta, tks, y);
596  return true;
597  } else {
598  return false;
599  }
600 }
601 
602 double DAClusterizerInZ_vect::beta0(double betamax, track_t const& tks, vertex_t const& y) const {
603  double T0 = 0; // max Tc for beta=0
604  // estimate critical temperature from beta=0 (T=inf)
605  const unsigned int nt = tks.getSize();
606  const unsigned int nv = y.getSize();
607 
608  for (unsigned int k = 0; k < nv; k++) {
609  // vertex fit at T=inf
610  double sumwz = 0;
611  double sumw = 0;
612  for (unsigned int i = 0; i < nt; i++) {
613  double w = tks.tkwt[i] * tks.dz2[i];
614  sumwz += w * tks.zpca[i];
615  sumw += w;
616  }
617 
618  y.zvtx[k] = sumwz / sumw;
619 
620  // estimate Tcrit
621  double a = 0, b = 0;
622  for (unsigned int i = 0; i < nt; i++) {
623  double dx = tks.zpca[i] - y.zvtx[k];
624  double w = tks.tkwt[i] * tks.dz2[i];
625  a += w * std::pow(dx, 2) * tks.dz2[i];
626  b += w;
627  }
628  double Tc = 2. * a / b; // the critical temperature of this vertex
629 
630  if (Tc > T0)
631  T0 = Tc;
632 
633  } // vertex loop (normally there should be only one vertex at beta=0)
634 
635 #ifdef DEBUG
636  if (DEBUGLEVEL > 0) {
637  std::cout << "DAClusterizerInZ_vect.beta0: Tc = " << T0 << std::endl;
638  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
639  std::cout << "DAClusterizerInZ_vect.beta0: nstep = " << coolingsteps << std::endl;
640  }
641 #endif
642 
643  if (T0 > 1. / betamax) {
644  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
645 
646  return betamax * std::pow(coolingFactor_, coolingsteps);
647  } else {
648  // ensure at least one annealing step
649  return betamax * coolingFactor_;
650  }
651 }
652 
653 bool DAClusterizerInZ_vect::split(const double beta, track_t& tks, vertex_t& y, double threshold) const {
654  // split only critical vertices (Tc >~ T=1/beta <==> beta*Tc>~1)
655  // returns true if at least one cluster was split
656 
657  // update the sums needed for Tc, rho0 is never non-zero while splitting is still active
658  update(beta, tks, y, 0., true); // the "true" flag enables Tc evaluation
659 
660  constexpr double epsilon = 1e-3; // minimum split size
661  unsigned int nv = y.getSize();
662 
663  // avoid left-right biases by splitting highest Tc first
664 
665  std::vector<std::pair<double, unsigned int>> critical;
666  for (unsigned int k = 0; k < nv; k++) {
667  double Tc = 2 * y.swE[k] / y.sw[k];
668  if (edm::isFinite(Tc) and beta * Tc > threshold) {
669  critical.push_back(make_pair(Tc, k));
670  }
671  }
672  if (critical.empty())
673  return false;
674 
675  std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int>>());
676 
677  bool split = false;
678  const unsigned int nt = tks.getSize();
679 
680  for (unsigned int ic = 0; ic < critical.size(); ic++) {
681  unsigned int k = critical[ic].second;
682 
683  // estimate subcluster positions and weight
684  double p1 = 0, z1 = 0, w1 = 0;
685  double p2 = 0, z2 = 0, w2 = 0;
686  for (unsigned int i = 0; i < nt; i++) {
687  if (tks.sum_Z[i] > 1.e-100) {
688  // winner-takes-all, usually overestimates splitting
689  double tl = tks.zpca[i] < y.zvtx[k] ? 1. : 0.;
690  double tr = 1. - tl;
691 
692  // soften it, especially at low T
693  double arg = (tks.zpca[i] - y.zvtx[k]) * sqrt(beta * tks.dz2[i]);
694  if (std::fabs(arg) < 20) {
695  double t = local_exp(-arg);
696  tl = t / (t + 1.);
697  tr = 1 / (t + 1.);
698  }
699 
700  double p = y.rho[k] * tks.tkwt[i] * local_exp(-beta * Eik(tks.zpca[i], y.zvtx[k], tks.dz2[i])) / tks.sum_Z[i];
701  double w = p * tks.dz2[i];
702  p1 += p * tl;
703  z1 += w * tl * tks.zpca[i];
704  w1 += w * tl;
705  p2 += p * tr;
706  z2 += w * tr * tks.zpca[i];
707  w2 += w * tr;
708  }
709  }
710 
711  if (w1 > 0) {
712  z1 = z1 / w1;
713  } else {
714  z1 = y.zvtx[k] - epsilon;
715  }
716  if (w2 > 0) {
717  z2 = z2 / w2;
718  } else {
719  z2 = y.zvtx[k] + epsilon;
720  }
721 
722  // reduce split size if there is not enough room
723  if ((k > 0) && (z1 < (0.6 * y.zvtx[k] + 0.4 * y.zvtx[k - 1]))) {
724  z1 = 0.6 * y.zvtx[k] + 0.4 * y.zvtx[k - 1];
725  }
726  if ((k + 1 < nv) && (z2 > (0.6 * y.zvtx[k] + 0.4 * y.zvtx[k + 1]))) {
727  z2 = 0.6 * y.zvtx[k] + 0.4 * y.zvtx[k + 1];
728  }
729 
730 #ifdef DEBUG
731  assert(k < nv);
732  if (DEBUGLEVEL > 1) {
733  if (std::fabs(y.zvtx[k] - zdumpcenter_) < zdumpwidth_) {
734  std::cout << " T= " << std::setw(8) << 1. / beta << " Tc= " << critical[ic].first << " splitting "
735  << std::fixed << std::setprecision(4) << y.zvtx[k] << " --> " << z1 << "," << z2 << " [" << p1
736  << "," << p2 << "]";
737  if (std::fabs(z2 - z1) > epsilon) {
738  std::cout << std::endl;
739  } else {
740  std::cout << " rejected " << std::endl;
741  }
742  }
743  }
744 #endif
745 
746  // split if the new subclusters are significantly separated
747  if ((z2 - z1) > epsilon) {
748  split = true;
749  double pk1 = p1 * y.rho[k] / (p1 + p2);
750  double pk2 = p2 * y.rho[k] / (p1 + p2);
751 
752  if (not(edm::isFinite(pk1) and edm::isFinite(pk2)))
753  continue;
754 
755  y.zvtx[k] = z2;
756  y.rho[k] = pk2;
757  y.insertItem(k, z1, pk1, tks);
758  if (k == 0)
759  y.extractRaw();
760 
761  nv++;
762 
763  // adjust remaining pointers
764  for (unsigned int jc = ic; jc < critical.size(); jc++) {
765  if (critical[jc].second >= k) {
766  critical[jc].second++;
767  }
768  }
769  } else {
770 #ifdef DEBUG
771  std::cout << "warning ! split rejected, too small." << endl;
772 #endif
773  }
774  }
775 
776  return split;
777 }
778 
779 vector<TransientVertex> DAClusterizerInZ_vect::vertices_no_blocks(const vector<reco::TransientTrack>& tracks) const {
780  track_t&& tks = fill(tracks);
781  vector<TransientVertex> clusters;
782  if (tks.getSize() == 0)
783  return clusters;
784  tks.extractRaw();
785 
786  double rho0 = 0.0; // start with no outlier rejection
787 
788  vertex_t y; // the vertex prototypes
789 
790  // initialize:single vertex at infinite temperature
791  y.addItem(0, 1.0);
792  clear_vtx_range(tks, y);
793 
794  // estimate first critical temperature
795  double beta = beta0(betamax_, tks, y);
796 #ifdef DEBUG
797  if (DEBUGLEVEL > 0)
798  std::cout << "Beta0 is " << beta << std::endl;
799 #endif
800 
801  thermalize(beta, tks, y, delta_highT_);
802 
803  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
804 
805  double betafreeze = betamax_ * sqrt(coolingFactor_);
806 
807  while (beta < betafreeze) {
808  while (merge(y, tks, beta)) {
809  update(beta, tks, y, rho0, false);
810  }
811  split(beta, tks, y);
812 
813  beta = beta / coolingFactor_;
814  thermalize(beta, tks, y, delta_highT_);
815  }
816 
817 #ifdef DEBUG
818  verify(y, tks);
819 
820  if (DEBUGLEVEL > 0) {
821  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
822  << "last round of splitting" << std::endl;
823  }
824 #endif
825 
826  set_vtx_range(beta, tks, y);
827  update(beta, tks, y, rho0, false);
828 
829  while (merge(y, tks, beta)) {
830  set_vtx_range(beta, tks, y);
831  update(beta, tks, y, rho0, false);
832  }
833 
834  unsigned int ntry = 0;
835  double threshold = 1.0;
836  while (split(beta, tks, y, threshold) && (ntry++ < 10)) {
837  thermalize(beta, tks, y, delta_highT_, rho0); // rho0 = 0. here
838  while (merge(y, tks, beta)) {
839  update(beta, tks, y, rho0, false);
840  }
841 
842  // relax splitting a bit to reduce multiple split-merge cycles of the same cluster
843  threshold *= 1.1;
844  }
845 
846 #ifdef DEBUG
847  verify(y, tks);
848  if (DEBUGLEVEL > 0) {
849  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
850  << "turning on outlier rejection at T=" << 1 / beta << std::endl;
851  }
852 #endif
853 
854  // switch on outlier rejection at T=Tmin
855  if (dzCutOff_ > 0) {
856  rho0 = y.getSize() > 1 ? 1. / y.getSize() : 1.;
857  for (unsigned int a = 0; a < 5; a++) {
858  update(beta, tks, y, a * rho0 / 5.); // adiabatic turn-on
859  }
860  }
861 
862  thermalize(beta, tks, y, delta_lowT_, rho0);
863 
864 #ifdef DEBUG
865  verify(y, tks);
866  if (DEBUGLEVEL > 0) {
867  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
868  << "merging with outlier rejection at T=" << 1 / beta << std::endl;
869  }
870  if (DEBUGLEVEL > 2)
871  dump(beta, y, tks, 2, rho0);
872 #endif
873 
874  // merge again (some cluster split by outliers collapse here)
875  while (merge(y, tks, beta)) {
876  set_vtx_range(beta, tks, y);
877  update(beta, tks, y, rho0, false);
878  }
879 
880 #ifdef DEBUG
881  verify(y, tks);
882  if (DEBUGLEVEL > 0) {
883  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
884  << "after merging with outlier rejection at T=" << 1 / beta << std::endl;
885  }
886  if (DEBUGLEVEL > 2)
887  dump(beta, y, tks, 2, rho0);
888 #endif
889 
890  // go down to the purging temperature (if it is lower than tmin)
891  while (beta < betapurge_) {
892  beta = min(beta / coolingFactor_, betapurge_);
893  thermalize(beta, tks, y, delta_lowT_, rho0);
894  }
895 
896 #ifdef DEBUG
897  verify(y, tks);
898  if (DEBUGLEVEL > 0) {
899  std::cout << "DAClusterizerInZ_vect::vertices :"
900  << "purging at T=" << 1 / beta << std::endl;
901  }
902 #endif
903 
904  // eliminate insignificant vertices, this is more restrictive at higher T
905  while (purge(y, tks, rho0, beta)) {
906  thermalize(beta, tks, y, delta_lowT_, rho0);
907  }
908 
909 #ifdef DEBUG
910  verify(y, tks);
911  if (DEBUGLEVEL > 0) {
912  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
913  << "last cooling T=" << 1 / beta << std::endl;
914  }
915 #endif
916 
917  // optionally cool some more without doing anything, to make the assignment harder
918  while (beta < betastop_) {
919  beta = min(beta / coolingFactor_, betastop_);
920  thermalize(beta, tks, y, delta_lowT_, rho0);
921  }
922 
923 #ifdef DEBUG
924  verify(y, tks);
925  if (DEBUGLEVEL > 0) {
926  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
927  << "stop cooling at T=" << 1 / beta << std::endl;
928  }
929  if (DEBUGLEVEL > 2)
930  dump(beta, y, tks, 2, rho0);
931 #endif
932 
933  // assign tracks and fill into transient vertices
934  return fill_vertices(beta, rho0, tks, y);
935 }
936 
937 vector<TransientVertex> DAClusterizerInZ_vect::vertices_in_blocks(const vector<reco::TransientTrack>& tracks) const {
938  vector<reco::TransientTrack> sorted_tracks;
939  vector<pair<float, float>> vertices_tot; // z, rho for each vertex
940  for (unsigned int i = 0; i < tracks.size(); i++) {
941  sorted_tracks.push_back(tracks[i]);
942  }
943  double rho0, beta;
944  std::sort(sorted_tracks.begin(),
945  sorted_tracks.end(),
946  [](const reco::TransientTrack& a, const reco::TransientTrack& b) -> bool {
947  return (a.stateAtBeamLine().trackStateAtPCA()).position().z() <
948  (b.stateAtBeamLine().trackStateAtPCA()).position().z();
949  });
950 
951  unsigned int nBlocks = (unsigned int)std::floor(sorted_tracks.size() / (block_size_ * (1 - overlap_frac_)));
952  if (nBlocks < 1) {
953  nBlocks = 1;
954  edm::LogWarning("DAClusterizerinZ_vect")
955  << "Warning nBlocks was 0 with ntracks = " << sorted_tracks.size() << " block_size = " << block_size_
956  << " and overlap fraction = " << overlap_frac_ << ". Setting nBlocks = 1";
957  }
958  for (unsigned int block = 0; block < nBlocks; block++) {
959  vector<reco::TransientTrack> block_tracks;
960  unsigned int begin = (unsigned int)(block * block_size_ * (1 - overlap_frac_));
961  unsigned int end = (unsigned int)std::min(begin + block_size_, (unsigned int)sorted_tracks.size());
962  for (unsigned int i = begin; i < end; i++) {
963  block_tracks.push_back(sorted_tracks[i]);
964  }
965  if (block_tracks.empty()) {
966  continue;
967  }
968 
969 #ifdef DEBUG
970  std::cout << "Running vertices_in_blocks on" << std::endl;
971  std::cout << "- block no." << block << " on " << nBlocks << " blocks " << std::endl;
972  std::cout << "- block track size: " << sorted_tracks.size() << " - block size: " << block_size_ << std::endl;
973 #endif
974  track_t&& tks = fill(block_tracks);
975  tks.extractRaw();
976 
977  rho0 = 0.0; // start with no outlier rejection
978 
979  vertex_t y; // the vertex prototypes
980 
981  // initialize:single vertex at infinite temperature
982  y.addItem(0, 1.0);
983  clear_vtx_range(tks, y);
984 
985  // estimate first critical temperature
986  beta = beta0(betamax_, tks, y);
987 #ifdef DEBUG
988  if (DEBUGLEVEL > 0)
989  std::cout << "Beta0 is " << beta << std::endl;
990 #endif
991 
992  thermalize(beta, tks, y, delta_highT_);
993 
994  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
995 
996  double betafreeze = betamax_ * sqrt(coolingFactor_);
997  while (beta < betafreeze) {
998  while (merge(y, tks, beta)) {
999  update(beta, tks, y, rho0, false);
1000  }
1001  split(beta, tks, y);
1002 
1003  beta = beta / coolingFactor_;
1004  thermalize(beta, tks, y, delta_highT_);
1005  }
1006 
1007 #ifdef DEBUG
1008  verify(y, tks);
1009 
1010  if (DEBUGLEVEL > 0) {
1011  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1012  << "last round of splitting" << std::endl;
1013  }
1014 #endif
1015 
1016  set_vtx_range(beta, tks, y);
1017  update(beta, tks, y, rho0, false);
1018 
1019  while (merge(y, tks, beta)) {
1020  set_vtx_range(beta, tks, y);
1021  update(beta, tks, y, rho0, false);
1022  }
1023 
1024  unsigned int ntry = 0;
1025  double threshold = 1.0;
1026  while (split(beta, tks, y, threshold) && (ntry++ < 10)) {
1027  thermalize(beta, tks, y, delta_highT_, rho0); // rho0 = 0. here
1028  while (merge(y, tks, beta)) {
1029  update(beta, tks, y, rho0, false);
1030  }
1031 
1032  // relax splitting a bit to reduce multiple split-merge cycles of the same cluster
1033  threshold *= 1.1;
1034  }
1035 
1036 #ifdef DEBUG
1037  verify(y, tks);
1038  if (DEBUGLEVEL > 0) {
1039  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1040  << "turning on outlier rejection at T=" << 1 / beta << std::endl;
1041  }
1042 #endif
1043 
1044  // switch on outlier rejection at T=Tmin
1045  if (dzCutOff_ > 0) {
1046  rho0 = y.getSize() > 1 ? 1. / y.getSize() : 1.;
1047  for (unsigned int a = 0; a < 5; a++) {
1048  update(beta, tks, y, a * rho0 / 5.); // adiabatic turn-on
1049  }
1050  }
1051 
1052  thermalize(beta, tks, y, delta_lowT_, rho0);
1053 
1054 #ifdef DEBUG
1055  verify(y, tks);
1056  if (DEBUGLEVEL > 0) {
1057  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1058  << "merging with outlier rejection at T=" << 1 / beta << std::endl;
1059  }
1060  if (DEBUGLEVEL > 2)
1061  dump(beta, y, tks, 2, rho0);
1062 #endif
1063 
1064  // merge again (some cluster split by outliers collapse here)
1065  while (merge(y, tks, beta)) {
1066  set_vtx_range(beta, tks, y);
1067  update(beta, tks, y, rho0, false);
1068  }
1069 
1070 #ifdef DEBUG
1071  verify(y, tks);
1072  if (DEBUGLEVEL > 0) {
1073  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1074  << "after merging with outlier rejection at T=" << 1 / beta << std::endl;
1075  }
1076  if (DEBUGLEVEL > 2)
1077  dump(beta, y, tks, 2, rho0);
1078 #endif
1079 
1080  // go down to the purging temperature (if it is lower than tmin)
1081  while (beta < betapurge_) {
1082  beta = min(beta / coolingFactor_, betapurge_);
1083  thermalize(beta, tks, y, delta_lowT_, rho0);
1084  }
1085 
1086 #ifdef DEBUG
1087  verify(y, tks);
1088  if (DEBUGLEVEL > 0) {
1089  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1090  << "purging at T=" << 1 / beta << std::endl;
1091  }
1092 #endif
1093 
1094  // eliminate insignificant vertices, this is more restrictive at higher T
1095  while (purge(y, tks, rho0, beta)) {
1096  thermalize(beta, tks, y, delta_lowT_, rho0);
1097  }
1098 
1099 #ifdef DEBUG
1100  verify(y, tks);
1101  if (DEBUGLEVEL > 0) {
1102  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1103  << "last cooling T=" << 1 / beta << std::endl;
1104  }
1105 #endif
1106 
1107  // optionally cool some more without doing anything, to make the assignment harder
1108  while (beta < betastop_) {
1109  beta = min(beta / coolingFactor_, betastop_);
1110  thermalize(beta, tks, y, delta_lowT_, rho0);
1111  }
1112 
1113 #ifdef DEBUG
1114  verify(y, tks);
1115  if (DEBUGLEVEL > 0) {
1116  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1117  << "stop cooling at T=" << 1 / beta << std::endl;
1118  }
1119  if (DEBUGLEVEL > 2)
1120  dump(beta, y, tks, 2, rho0);
1121 #endif
1122 
1123  for (unsigned int ivertex = 0; ivertex < y.getSize(); ivertex++) {
1124  if (y.zvtx_vec[ivertex] != 0 && y.rho_vec[ivertex] != 0) {
1125  vertices_tot.push_back(pair(y.zvtx_vec[ivertex], y.rho_vec[ivertex]));
1126 #ifdef DEBUG
1127  std::cout << "Found new vertex " << y.zvtx_vec[ivertex] << " , " << y.rho_vec[ivertex] << std::endl;
1128 #endif
1129  }
1130  }
1131  }
1132 
1133  std::sort(vertices_tot.begin(),
1134  vertices_tot.end(),
1135  [](const pair<float, float>& a, const pair<float, float>& b) -> bool { return a.first < b.first; });
1136 
1137  // reassign tracks to vertices
1138  track_t&& tracks_tot = fill(tracks);
1139  const unsigned int nv = vertices_tot.size();
1140  const unsigned int nt = tracks_tot.getSize();
1141 
1142  for (auto itrack = 0U; itrack < nt; ++itrack) {
1143  double zrange = max(sel_zrange_ / sqrt(beta * tracks_tot.dz2[itrack]), zrange_min_);
1144 
1145  double zmin = tracks_tot.zpca[itrack] - zrange;
1146  unsigned int kmin = min(nv - 1, tracks_tot.kmin[itrack]);
1147  // find the smallest vertex_z that is larger than zmin
1148  if (vertices_tot[kmin].first > zmin) {
1149  while ((kmin > 0) && (vertices_tot[kmin - 1].first > zmin)) {
1150  kmin--;
1151  }
1152  } else {
1153  while ((kmin < (nv - 1)) && (vertices_tot[kmin].first < zmin)) {
1154  kmin++;
1155  }
1156  }
1157 
1158  double zmax = tracks_tot.zpca[itrack] + zrange;
1159  unsigned int kmax = min(nv - 1, tracks_tot.kmax[itrack] - 1);
1160  // note: kmax points to the last vertex in the range, while gtracks.kmax points to the entry BEHIND the last vertex
1161  // find the largest vertex_z that is smaller than zmax
1162  if (vertices_tot[kmax].first < zmax) {
1163  while ((kmax < (nv - 1)) && (vertices_tot[kmax + 1].first < zmax)) {
1164  kmax++;
1165  }
1166  } else {
1167  while ((kmax > 0) && (vertices_tot[kmax].first > zmax)) {
1168  kmax--;
1169  }
1170  }
1171 
1172  if (kmin <= kmax) {
1173  tracks_tot.kmin[itrack] = kmin;
1174  tracks_tot.kmax[itrack] = kmax + 1;
1175  } else {
1176  tracks_tot.kmin[itrack] = max(0U, min(kmin, kmax));
1177  tracks_tot.kmax[itrack] = min(nv, max(kmin, kmax) + 1);
1178  }
1179  }
1180 
1181  rho0 = nv > 1 ? 1. / nv : 1.;
1182  const auto z_sum_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
1183 
1184  std::vector<std::vector<unsigned int>> vtx_track_indices(nv);
1185  for (unsigned int i = 0; i < nt; i++) {
1186  const auto kmin = tracks_tot.kmin[i];
1187  const auto kmax = tracks_tot.kmax[i];
1188  double p_max = -1;
1189  unsigned int iMax = 10000; //just a "big" number w.r.t. number of vertices
1190  float sum_Z = z_sum_init;
1191  for (auto k = kmin; k < kmax; k++) {
1192  float v_exp = local_exp(-beta * Eik(tracks_tot.zpca[i], vertices_tot[k].first, tracks_tot.dz2[i]));
1193  sum_Z += vertices_tot[k].second * v_exp;
1194  }
1195  double invZ = sum_Z > 1e-100 ? 1. / sum_Z : 0.0;
1196  for (auto k = kmin; k < kmax && invZ != 0.0; k++) {
1197  float v_exp = local_exp(-beta * Eik(tracks_tot.zpca[i], vertices_tot[k].first, tracks_tot.dz2[i]));
1198  double p = vertices_tot[k].second * v_exp * invZ;
1199  if (p > p_max && p > mintrkweight_) {
1200  p_max = p;
1201  iMax = k;
1202  }
1203  }
1204  if (iMax < vtx_track_indices.size()) {
1205  vtx_track_indices[iMax].push_back(i);
1206  }
1207  }
1208 #ifdef DEBUG
1209  for (auto itrack = 0U; itrack < nt; ++itrack) {
1210  std::cout << "itrack " << itrack << " , " << tracks_tot.kmin[itrack] << " , " << tracks_tot.kmax[itrack]
1211  << std::endl;
1212  }
1213 #endif
1214 
1215  vector<TransientVertex> clusters;
1216  if (nv == 0) {
1217  return clusters;
1218  }
1219 
1220  GlobalError dummyError(0.01, 0, 0.01, 0., 0., 0.01);
1221  vector<reco::TransientTrack> vertexTracks;
1222 
1223  for (unsigned int k = 0; k < nv; k++) {
1224  if (!vtx_track_indices[k].empty()) {
1225  for (auto i : vtx_track_indices[k]) {
1226  vertexTracks.push_back(*(tracks_tot.tt[i]));
1227 #ifdef DEBUG
1228  std::cout << vertices_tot[k].first << ","
1229  << (*(tracks_tot.tt[i])).stateAtBeamLine().trackStateAtPCA().position().z() << std::endl;
1230 #endif
1231  }
1232  }
1233 
1234  // implement what clusterize() did before : merge left-to-right if distance < 2 * vertexSize_
1235  if ((k + 1 == nv) || (abs(vertices_tot[k + 1].first - vertices_tot[k].first) > (2 * vertexSize_))) {
1236  // close a cluster
1237  if (vertexTracks.size() > 1) {
1238  GlobalPoint pos(0, 0, vertices_tot[k].first); // only usable with subsequent fit
1239  TransientVertex v(pos, dummyError, vertexTracks, 0);
1240  clusters.push_back(v);
1241  }
1242  vertexTracks.clear();
1243  }
1244  }
1245 
1246  return clusters;
1247 } // end of vertices_in_blocks
1248 
1249 vector<TransientVertex> DAClusterizerInZ_vect::fill_vertices(double beta, double rho0, track_t& tks, vertex_t& y) const {
1250  // select significant tracks and use a TransientVertex as a container
1251 
1252  set_vtx_range(beta, tks, y);
1253  const unsigned int nv = y.getSize();
1254  for (unsigned int k = 0; k < nv; k++) {
1255  if (edm::isNotFinite(y.rho[k]) || edm::isNotFinite(y.zvtx[k])) {
1256  y.rho[k] = 0;
1257  y.zvtx[k] = 0;
1258  }
1259  }
1260 
1261  // ensure consistent assignment probabillities and make a hard assignment
1262  const unsigned int nt = tks.getSize();
1263  const auto z_sum_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
1264  std::vector<std::vector<unsigned int>> vtx_track_indices(nv);
1265  std::vector<std::vector<float>> vtx_track_weights(nv);
1266  for (unsigned int i = 0; i < nt; i++) {
1267  const auto kmin = tks.kmin[i];
1268  const auto kmax = tks.kmax[i];
1269  for (auto k = kmin; k < kmax; k++) {
1270  y.exp_arg[k] = -beta * Eik(tks.zpca[i], y.zvtx[k], tks.dz2[i]);
1271  }
1272 
1273  local_exp_list_range(y.exp_arg, y.exp, kmin, kmax);
1274 
1275  tks.sum_Z[i] = z_sum_init;
1276  for (auto k = kmin; k < kmax; k++) {
1277  tks.sum_Z[i] += y.rho[k] * y.exp[k];
1278  }
1279  const double invZ = tks.sum_Z[i] > 1e-100 ? 1. / tks.sum_Z[i] : 0.0;
1280 
1281  double pmax = -1;
1282  unsigned int k_pmax = 0;
1283  for (auto k = kmin; k < kmax; k++) {
1284  double p = y.rho[k] * y.exp[k] * invZ;
1285  if (p > pmax) {
1286  pmax = p;
1287  k_pmax = k;
1288  }
1289  }
1290 
1291  if (pmax > mintrkweight_) {
1292  // assign to the cluster with the highest assignment weight, if it is at least mintrkweight_
1293  vtx_track_indices[k_pmax].push_back(i);
1294  vtx_track_weights[k_pmax].push_back(pmax);
1295  }
1296  }
1297 
1298  // fill transient vertices
1299  // the position is normally not used, probably not optimal when Tstop <> 2, anyway
1300  vector<TransientVertex> clusters;
1301  for (unsigned int k = 0; k < nv; k++) {
1302  double sump = 0;
1303  double sumw = 0;
1304  double sumwp = 0, sumwz = 0;
1305  if (!vtx_track_indices[k].empty()) {
1306  vector<reco::TransientTrack> vertexTracks;
1308  unsigned int j = 0;
1309  for (auto i : vtx_track_indices[k]) {
1310  auto p = vtx_track_weights[k][j];
1311  vertexTracks.push_back(*(tks.tt[i]));
1312  trkWeightMap[vertexTracks[j]] = p;
1313  auto w = p * tks.dz2[i];
1314  sump += p;
1315  sumw += w;
1316  sumwp += w * p;
1317  sumwz += w * tks.zpca[i];
1318  j++;
1319  }
1320  float zerror_squared = 1.; //
1321  if ((sumw > 0) && (sumwp > 0)) {
1322  zerror_squared = sumwp / (sumw * sumw);
1323  y.zvtx[k] = sumwz / sumw;
1324  }
1325 
1326  reco::BeamSpot bs = vertexTracks[0].stateAtBeamLine().beamSpot();
1327  GlobalPoint pos(bs.x(y.zvtx[k]), bs.y(y.zvtx[k]), y.zvtx[k]);
1328  const float xerror_squared = pow(bs.BeamWidthX(), 2);
1329  const float yerror_squared = pow(bs.BeamWidthY(), 2);
1330  GlobalError err(xerror_squared, 0, yerror_squared, 0., 0., zerror_squared);
1331  TransientVertex v(pos, err, vertexTracks, 0, 2 * sump - 3.);
1332  v.weightMap(trkWeightMap);
1333  clusters.push_back(v);
1334  }
1335  }
1336 
1337  return clusters;
1338 }
1339 
1340 vector<TransientVertex> DAClusterizerInZ_vect::vertices(const vector<reco::TransientTrack>& tracks) const {
1341  if (runInBlocks_ and (block_size_ < tracks.size())) //doesn't bother if low number of tracks
1342  return vertices_in_blocks(tracks);
1343  else
1344  return vertices_no_blocks(tracks);
1345 }
1346 
1347 vector<vector<reco::TransientTrack>> DAClusterizerInZ_vect::clusterize( // OBSOLETE
1348  const vector<reco::TransientTrack>& tracks) const {
1349  vector<vector<reco::TransientTrack>> clusters;
1350  vector<TransientVertex>&& pv = vertices(tracks);
1351 
1352 #ifdef DEBUG
1353  if (DEBUGLEVEL > 0) {
1354  std::cout << "###################################################" << endl;
1355  std::cout << "# vectorized DAClusterizerInZ_vect::clusterize nt=" << tracks.size() << endl;
1356  std::cout << "# DAClusterizerInZ_vect::clusterize pv.size=" << pv.size() << endl;
1357  std::cout << "###################################################" << endl;
1358  }
1359 #endif
1360 
1361  if (pv.empty()) {
1362  return clusters;
1363  }
1364 
1365  // fill into clusters and merge
1366  vector<reco::TransientTrack> aCluster = pv.begin()->originalTracks();
1367 
1368  for (auto k = pv.begin() + 1; k != pv.end(); k++) {
1369  if (std::abs(k->position().z() - (k - 1)->position().z()) > (2 * vertexSize_)) {
1370  // close a cluster
1371  if (aCluster.size() > 1) {
1372  clusters.push_back(aCluster);
1373  }
1374 #ifdef DEBUG
1375  else {
1376  std::cout << " one track cluster at " << k->position().z() << " suppressed" << std::endl;
1377  }
1378 #endif
1379  aCluster.clear();
1380  }
1381  for (unsigned int i = 0; i < k->originalTracks().size(); i++) {
1382  aCluster.push_back(k->originalTracks()[i]);
1383  }
1384  }
1385  clusters.emplace_back(std::move(aCluster));
1386 
1387  return clusters;
1388 }
1389 
1391  const double beta, const vertex_t& y, const track_t& tks, const int verbosity, const double rho0) const {
1392 #ifdef DEBUG
1393  const unsigned int nv = y.getSize();
1394  const unsigned int nt = tks.getSize();
1395  // make a local copies of clusters and tracks to update Tc [ = y_local.swE / y_local.sw ]
1396  vertex_t y_local = y;
1397  track_t tks_local = tks;
1398  update(beta, tks_local, y_local, rho0, true);
1399 
1400  std::vector<unsigned int> iz;
1401  for (unsigned int j = 0; j < nt; j++) {
1402  iz.push_back(j);
1403  }
1404  std::sort(iz.begin(), iz.end(), [tks](unsigned int a, unsigned int b) { return tks.zpca[a] < tks.zpca[b]; });
1405  std::cout << std::endl;
1406  std::cout << "-----DAClusterizerInZ::dump ----" << nv << " clusters " << std::endl;
1407  std::cout << " ";
1408  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1409  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1410  std::cout << " " << setw(3) << ivertex << " ";
1411  }
1412  }
1413  std::cout << endl;
1414  std::cout << " z= ";
1415  std::cout << setprecision(4);
1416  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1417  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1418  std::cout << setw(8) << fixed << y.zvtx[ivertex];
1419  }
1420  }
1421  std::cout << endl
1422  << "T=" << setw(15) << 1. / beta << " Tmin =" << setw(10) << 1. / betamax_
1423  << " Tc= ";
1424  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1425  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1426  double Tc = 2 * y_local.swE[ivertex] / y_local.sw[ivertex];
1427  std::cout << setw(8) << fixed << setprecision(1) << Tc;
1428  }
1429  }
1430  std::cout << endl;
1431 
1432  std::cout << " pk= ";
1433  double sumpk = 0;
1434  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1435  sumpk += y.rho[ivertex];
1436  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1437  continue;
1438  std::cout << setw(8) << setprecision(4) << fixed << y.rho[ivertex];
1439  }
1440  std::cout << endl;
1441 
1442  std::cout << " nt= ";
1443  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1444  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1445  continue;
1446  std::cout << setw(8) << setprecision(1) << fixed << y.rho[ivertex] * nt;
1447  }
1448  std::cout << endl;
1449 
1450  if (verbosity > 0) {
1451  double E = 0, F = 0;
1452  std::cout << endl;
1453  std::cout << "---- z +/- dz ip +/-dip pt phi eta weights ----" << endl;
1454  std::cout << setprecision(4);
1455  for (unsigned int i0 = 0; i0 < nt; i0++) {
1456  unsigned int i = iz[i0];
1457  if (tks.sum_Z[i] > 0) {
1458  F -= std::log(tks.sum_Z[i]) / beta;
1459  }
1460  double tz = tks.zpca[i];
1461 
1462  if (std::fabs(tz - zdumpcenter_) > zdumpwidth_)
1463  continue;
1464  std::cout << setw(4) << i << ")" << setw(8) << fixed << setprecision(4) << tz << " +/-" << setw(6)
1465  << sqrt(1. / tks.dz2[i]);
1466  if ((tks.tt[i] == nullptr)) {
1467  std::cout << " effective track ";
1468  } else {
1469  if (tks.tt[i]->track().quality(reco::TrackBase::highPurity)) {
1470  std::cout << " *";
1471  } else {
1472  std::cout << " ";
1473  }
1474  if (tks.tt[i]->track().hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1)) {
1475  std::cout << "+";
1476  } else {
1477  std::cout << "-";
1478  }
1479  std::cout << setw(1)
1480  << tks.tt[i]
1481  ->track()
1482  .hitPattern()
1483  .pixelBarrelLayersWithMeasurement(); // see DataFormats/TrackReco/interface/HitPattern.h
1484  std::cout << setw(1) << tks.tt[i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
1485  std::cout << setw(1) << hex
1486  << tks.tt[i]->track().hitPattern().trackerLayersWithMeasurement() -
1487  tks.tt[i]->track().hitPattern().pixelLayersWithMeasurement()
1488  << dec;
1489  std::cout << "=" << setw(1) << hex
1490  << tks.tt[i]->track().hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS) << dec;
1491 
1492  Measurement1D IP = tks.tt[i]->stateAtBeamLine().transverseImpactParameter();
1493  std::cout << setw(8) << IP.value() << "+/-" << setw(6) << IP.error();
1494  std::cout << " " << setw(6) << setprecision(2) << tks.tt[i]->track().pt() * tks.tt[i]->track().charge();
1495  std::cout << " " << setw(5) << setprecision(2) << tks.tt[i]->track().phi() << " " << setw(5) << setprecision(2)
1496  << tks.tt[i]->track().eta();
1497  } // not a dummy track
1498 
1499  double sump = 0.;
1500  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1501  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1502  continue;
1503 
1504  if ((tks.tkwt[i] > 0) && (tks.sum_Z[i] > 0)) {
1505  //double p=pik(beta,tks[i],*k);
1506  double p = y.rho[ivertex] * local_exp(-beta * Eik(tks.zpca[i], y.zvtx[ivertex], tks.dz2[i])) / tks.sum_Z[i];
1507  if (p > 0.0001) {
1508  std::cout << setw(8) << setprecision(3) << p;
1509  } else {
1510  std::cout << " . ";
1511  }
1512  E += p * Eik(tks.zpca[i], y.zvtx[ivertex], tks.dz2[i]);
1513  sump += p;
1514  } else {
1515  std::cout << " ";
1516  }
1517  }
1518  std::cout << " ( " << std::setw(3) << tks.kmin[i] << "," << std::setw(3) << tks.kmax[i] - 1 << " ) ";
1519  std::cout << endl;
1520  }
1521  std::cout << " ";
1522  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1523  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1524  std::cout << " " << setw(3) << ivertex << " ";
1525  }
1526  }
1527  std::cout << endl;
1528  std::cout << " z= ";
1529  std::cout << setprecision(4);
1530  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1531  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1532  std::cout << setw(8) << fixed << y.zvtx[ivertex];
1533  }
1534  }
1535  std::cout << endl;
1536  std::cout << endl
1537  << "T=" << 1 / beta << " E=" << E << " n=" << y.getSize() << " F= " << F << endl
1538  << "----------" << endl;
1539  }
1540 #endif
1541 }
1542 
1544  desc.addUntracked<double>("zdumpcenter", 0.);
1545  desc.addUntracked<double>("zdumpwidth", 20.);
1546  desc.add<double>("d0CutOff", 3.0);
1547  desc.add<double>("Tmin", 2.0);
1548  desc.add<double>("delta_lowT", 0.001);
1549  desc.add<double>("zmerge", 0.01);
1550  desc.add<double>("dzCutOff", 3.0);
1551  desc.add<double>("Tpurge", 2.0);
1552  desc.add<int>("convergence_mode", 0);
1553  desc.add<double>("delta_highT", 0.01);
1554  desc.add<double>("Tstop", 0.5);
1555  desc.add<double>("coolingFactor", 0.6);
1556  desc.add<double>("vertexSize", 0.006);
1557  desc.add<double>("uniquetrkweight", 0.8);
1558  desc.add<double>("uniquetrkminp", 0.0);
1559  desc.add<double>("zrange", 4.0);
1560  desc.add<bool>("runInBlocks", false);
1561  desc.add<unsigned int>("block_size", 10000);
1562  desc.add<double>("overlap_frac", 0.0);
1563 }
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
void addItemSorted(double new_zpca, double new_dz2, const reco::TransientTrack *new_tt, double new_tkwt)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool merge(vertex_t &y, track_t &tks, double &beta) const
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
void clear_vtx_range(track_t &gtracks, vertex_t &gvertices) const
int merge(int argc, char *argv[])
Definition: DiMuonVmerge.cc:28
T w() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
bool purge(vertex_t &, track_t &, double &, const double) const
std::map< reco::TransientTrack, float > TransientTrackToFloatMap
static void fillPSetDescription(edm::ParameterSetDescription &desc)
std::vector< const reco::TransientTrack * > tt
assert(be >=bs)
constexpr bool isFinite(T x)
A arg
Definition: Factorize.h:31
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const override
T getUntrackedParameter(std::string const &, T const &) const
U second(std::pair< T, U > const &p)
std::vector< TransientVertex > fill_vertices(double beta, double rho0, track_t &tracks, vertex_t &vertices) const
std::vector< unsigned int > kmin
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
std::vector< unsigned int > kmax
T sqrt(T t)
Definition: SSEVec.h:23
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0, const double rho0=0.) const
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< TransientVertex > vertices_no_blocks(const std::vector< reco::TransientTrack > &tracks) const
std::vector< TransientVertex > vertices_in_blocks(const std::vector< reco::TransientTrack > &tracks) const
int nt
Definition: AMPTWrapper.h:42
const int verbosity
DAClusterizerInZ_vect(const edm::ParameterSet &conf)
unsigned int thermalize(double beta, track_t &gtracks, vertex_t &gvertices, const double delta_max, const double rho0=0.) const
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const
double b
Definition: hdecay.h:120
#define update(a, b)
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
Log< level::Warning, false > LogWarning
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
void verify(const vertex_t &v, const track_t &tks, unsigned int nv=999999, unsigned int nt=999999) const
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
def move(src, dest)
Definition: eostools.py:511