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 
47 #ifdef DEBUG
48  std::cout << "DAClusterizerinZ_vect: mintrkweight = " << mintrkweight_ << std::endl;
49  std::cout << "DAClusterizerinZ_vect: uniquetrkweight = " << uniquetrkweight_ << std::endl;
50  std::cout << "DAClusterizerInZ_vect: uniquetrkminp = " << uniquetrkminp_ << std::endl;
51  std::cout << "DAClusterizerinZ_vect: zmerge = " << zmerge_ << std::endl;
52  std::cout << "DAClusterizerinZ_vect: Tmin = " << Tmin << std::endl;
53  std::cout << "DAClusterizerinZ_vect: Tpurge = " << Tpurge << std::endl;
54  std::cout << "DAClusterizerinZ_vect: Tstop = " << Tstop << std::endl;
55  std::cout << "DAClusterizerinZ_vect: vertexSize = " << vertexSize_ << std::endl;
56  std::cout << "DAClusterizerinZ_vect: coolingFactor = " << coolingFactor_ << std::endl;
57  std::cout << "DAClusterizerinZ_vect: d0CutOff = " << d0CutOff_ << std::endl;
58  std::cout << "DAClusterizerinZ_vect: dzCutOff = " << dzCutOff_ << std::endl;
59  std::cout << "DAClusterizerInZ_vect: zrange = " << sel_zrange_ << std::endl;
60  std::cout << "DAClusterizerinZ_vect: convergence mode = " << convergence_mode_ << std::endl;
61  std::cout << "DAClusterizerinZ_vect: delta_highT = " << delta_highT_ << std::endl;
62  std::cout << "DAClusterizerinZ_vect: delta_lowT = " << delta_lowT_ << std::endl;
63  std::cout << "DAClusterizerinZ_vect: DEBUGLEVEL " << DEBUGLEVEL << std::endl;
64 #endif
65 
66  if (convergence_mode_ > 1) {
67  edm::LogWarning("DAClusterizerinZ_vect")
68  << "DAClusterizerInZ_vect: invalid convergence_mode " << convergence_mode_ << " reset to default " << 0;
69  convergence_mode_ = 0;
70  }
71 
72  if (Tmin == 0) {
73  betamax_ = 1.0;
74  edm::LogWarning("DAClusterizerinZ_vect")
75  << "DAClusterizerInZ_vect: invalid Tmin " << Tmin << " reset to default " << 1. / betamax_;
76  } else {
77  betamax_ = 1. / Tmin;
78  }
79 
80  if ((Tpurge > Tmin) || (Tpurge == 0)) {
81  edm::LogWarning("DAClusterizerinZ_vect")
82  << "DAClusterizerInZ_vect: invalid Tpurge " << Tpurge << " set to " << Tmin;
83  Tpurge = Tmin;
84  }
85  betapurge_ = 1. / Tpurge;
86 
87  if ((Tstop > Tpurge) || (Tstop == 0)) {
88  edm::LogWarning("DAClusterizerinZ_vect")
89  << "DAClusterizerInZ_vect: invalid Tstop " << Tstop << " set to " << max(1., Tpurge);
90  Tstop = max(1., Tpurge);
91  }
92  betastop_ = 1. / Tstop;
93 }
94 
95 namespace {
96  inline double local_exp(double const& inp) { return vdt::fast_exp(inp); }
97 
98  inline void local_exp_list_range(double const* __restrict__ arg_inp,
99  double* __restrict__ arg_out,
100  const int kmin,
101  const int kmax) {
102  for (auto i = kmin; i != kmax; ++i)
103  arg_out[i] = vdt::fast_exp(arg_inp[i]);
104  }
105 
106 } // namespace
107 
108 void DAClusterizerInZ_vect::verify(const vertex_t& v, const track_t& tks, unsigned int nv, unsigned int nt) const {
109  if (!(nv == 999999)) {
110  assert(nv == v.getSize());
111  } else {
112  nv = v.getSize();
113  }
114 
115  if (!(nt == 999999)) {
116  assert(nt == tks.getSize());
117  } else {
118  nt = tks.getSize();
119  }
120 
121  assert(v.zvtx_vec.size() == nv);
122  assert(v.rho_vec.size() == nv);
123  assert(v.swz_vec.size() == nv);
124  assert(v.exp_arg_vec.size() == nv);
125  assert(v.exp_vec.size() == nv);
126  assert(v.se_vec.size() == nv);
127  assert(v.swz_vec.size() == nv);
128  assert(v.swE_vec.size() == nv);
129 
130  assert(v.zvtx == &v.zvtx_vec.front());
131  assert(v.rho == &v.rho_vec.front());
132  assert(v.exp_arg == &v.exp_arg_vec.front());
133  assert(v.sw == &v.sw_vec.front());
134  assert(v.swz == &v.swz_vec.front());
135  assert(v.se == &v.se_vec.front());
136  assert(v.swE == &v.swE_vec.front());
137 
138  for (unsigned int k = 0; k < nv - 1; k++) {
139  if (v.zvtx_vec[k] <= v.zvtx_vec[k + 1])
140  continue;
141  cout << " Z, cluster z-ordering assertion failure z[" << k << "] =" << v.zvtx_vec[k] << " z[" << k + 1
142  << "] =" << v.zvtx_vec[k + 1] << endl;
143  }
144 
145  assert(nt == tks.zpca_vec.size());
146  assert(nt == tks.dz2_vec.size());
147  assert(nt == tks.tt.size());
148  assert(nt == tks.tkwt_vec.size());
149  assert(nt == tks.sum_Z_vec.size());
150  assert(nt == tks.kmin.size());
151  assert(nt == tks.kmax.size());
152 
153  assert(tks.zpca == &tks.zpca_vec.front());
154  assert(tks.dz2 == &tks.dz2_vec.front());
155  assert(tks.tkwt == &tks.tkwt_vec.front());
156  assert(tks.sum_Z == &tks.sum_Z_vec.front());
157 
158  for (unsigned int i = 0; i < nt; i++) {
159  if ((tks.kmin[i] < tks.kmax[i]) && (tks.kmax[i] <= nv))
160  continue;
161  cout << "track vertex range assertion failure" << i << "/" << nt << " kmin,kmax=" << tks.kmin[i] << ", "
162  << tks.kmax[i] << " nv=" << nv << endl;
163  }
164 
165  for (unsigned int i = 0; i < nt; i++) {
166  assert((tks.kmin[i] < tks.kmax[i]) && (tks.kmax[i] <= nv));
167  }
168 }
169 
170 DAClusterizerInZ_vect::track_t DAClusterizerInZ_vect::fill(const vector<reco::TransientTrack>& tracks) const {
171  // prepare track data for clustering
172  track_t tks;
173  double sumtkwt = 0.;
174  for (auto it = tracks.begin(); it != tracks.end(); it++) {
175  if (!(*it).isValid())
176  continue;
177  double t_tkwt = 1.;
178  double t_z = ((*it).stateAtBeamLine().trackStateAtPCA()).position().z();
179  if (std::fabs(t_z) > 1000.)
180  continue;
181  auto const& t_mom = (*it).stateAtBeamLine().trackStateAtPCA().momentum();
182  // get the beam-spot
183  reco::BeamSpot beamspot = (it->stateAtBeamLine()).beamSpot();
184  double t_dz2 = std::pow((*it).track().dzError(), 2) // track errror
185  + (std::pow(beamspot.BeamWidthX() * t_mom.x(), 2) + std::pow(beamspot.BeamWidthY() * t_mom.y(), 2)) *
186  std::pow(t_mom.z(), 2) / std::pow(t_mom.perp2(), 2) // beam spot width
187  + std::pow(vertexSize_, 2); // intrinsic vertex size, safer for outliers and short lived decays
188  t_dz2 = 1. / t_dz2;
189  if (edm::isNotFinite(t_dz2) || t_dz2 < std::numeric_limits<double>::min())
190  continue;
191  if (d0CutOff_ > 0) {
192  Measurement1D atIP = (*it).stateAtBeamLine().transverseImpactParameter(); // error contains beamspot
193  t_tkwt = 1. / (1. + local_exp(std::pow(atIP.value() / atIP.error(), 2) -
194  std::pow(d0CutOff_, 2))); // reduce weight for high ip tracks
196  edm::LogWarning("DAClusterizerinZ_vect") << "rejected track t_tkwt " << t_tkwt;
197  continue; // usually is > 0.99
198  }
199  }
200  tks.addItemSorted(t_z, t_dz2, &(*it), t_tkwt);
201  sumtkwt += t_tkwt;
202  }
203 
204  tks.extractRaw();
205  tks.osumtkwt = sumtkwt > 0 ? 1. / sumtkwt : 0.;
206 
207 #ifdef DEBUG
208  if (DEBUGLEVEL > 0) {
209  std::cout << "Track count (Z) " << tks.getSize() << std::endl;
210  }
211 #endif
212 
213  return tks;
214 }
215 
216 namespace {
217  inline double Eik(double t_z, double k_z, double t_dz2) { return std::pow(t_z - k_z, 2) * t_dz2; }
218 } // namespace
219 
220 void DAClusterizerInZ_vect::set_vtx_range(double beta, track_t& gtracks, vertex_t& gvertices) const {
221  const unsigned int nv = gvertices.getSize();
222  const unsigned int nt = gtracks.getSize();
223 
224  if (nv == 0) {
225  edm::LogWarning("DAClusterizerinZ_vect") << "empty cluster list in set_vtx_range";
226  return;
227  }
228 
229  for (auto itrack = 0U; itrack < nt; ++itrack) {
230  double zrange = max(sel_zrange_ / sqrt(beta * gtracks.dz2[itrack]), zrange_min_);
231 
232  double zmin = gtracks.zpca[itrack] - zrange;
233  unsigned int kmin = min(nv - 1, gtracks.kmin[itrack]);
234  // find the smallest vertex_z that is larger than zmin
235  if (gvertices.zvtx[kmin] > zmin) {
236  while ((kmin > 0) && (gvertices.zvtx[kmin - 1] > zmin)) {
237  kmin--;
238  }
239  } else {
240  while ((kmin < (nv - 1)) && (gvertices.zvtx[kmin] < zmin)) {
241  kmin++;
242  }
243  }
244 
245  double zmax = gtracks.zpca[itrack] + zrange;
246  unsigned int kmax = min(nv - 1, gtracks.kmax[itrack] - 1);
247  // note: kmax points to the last vertex in the range, while gtracks.kmax points to the entry BEHIND the last vertex
248  // find the largest vertex_z that is smaller than zmax
249  if (gvertices.zvtx[kmax] < zmax) {
250  while ((kmax < (nv - 1)) && (gvertices.zvtx[kmax + 1] < zmax)) {
251  kmax++;
252  }
253  } else {
254  while ((kmax > 0) && (gvertices.zvtx[kmax] > zmax)) {
255  kmax--;
256  }
257  }
258 
259  if (kmin <= kmax) {
260  gtracks.kmin[itrack] = kmin;
261  gtracks.kmax[itrack] = kmax + 1;
262  } else {
263  gtracks.kmin[itrack] = max(0U, min(kmin, kmax));
264  gtracks.kmax[itrack] = min(nv, max(kmin, kmax) + 1);
265  }
266  }
267 }
268 
269 void DAClusterizerInZ_vect::clear_vtx_range(track_t& gtracks, vertex_t& gvertices) const {
270  const unsigned int nt = gtracks.getSize();
271  const unsigned int nv = gvertices.getSize();
272  for (auto itrack = 0U; itrack < nt; ++itrack) {
273  gtracks.kmin[itrack] = 0;
274  gtracks.kmax[itrack] = nv;
275  }
276 }
277 
279  double beta, track_t& gtracks, vertex_t& gvertices, const double rho0, const bool updateTc) const {
280  // update weights and vertex positions
281  // returns the maximum of changes of vertex positions
282  // sums needed for Tc are only updated if updateTC == true
283 
284  const unsigned int nt = gtracks.getSize();
285  const unsigned int nv = gvertices.getSize();
286  auto osumtkwt = gtracks.osumtkwt;
287 
288  double Z_init = 0;
289  if (rho0 > 0) {
290  Z_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_); // cut-off
291  }
292 
293  // define kernels
294  auto kernel_calc_exp_arg_range = [beta](const unsigned int itrack,
295  track_t const& tracks,
296  vertex_t const& vertices,
297  const unsigned int kmin,
298  const unsigned int kmax) {
299  const double track_z = tracks.zpca[itrack];
300  const double botrack_dz2 = -beta * tracks.dz2[itrack];
301 
302  // auto-vectorized
303  for (unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
304  auto mult_res = track_z - vertices.zvtx[ivertex];
305  vertices.exp_arg[ivertex] = botrack_dz2 * (mult_res * mult_res);
306  }
307  };
308 
309  auto kernel_add_Z_range = [Z_init](
310  vertex_t const& vertices, const unsigned int kmin, const unsigned int kmax) -> double {
311  double ZTemp = Z_init;
312  for (unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
313  ZTemp += vertices.rho[ivertex] * vertices.exp[ivertex];
314  }
315  return ZTemp;
316  };
317 
318  auto kernel_calc_normalization_range = [updateTc](const unsigned int track_num,
319  track_t& tracks,
321  const unsigned int kmin,
322  const unsigned int kmax) {
323  auto o_trk_sum_Z = tracks.tkwt[track_num] / tracks.sum_Z[track_num];
324  auto o_trk_dz2 = tracks.dz2[track_num];
325  auto tmp_trk_z = tracks.zpca[track_num];
326 
327  // auto-vectorized
328  if (updateTc) {
329 #pragma GCC ivdep
330  for (unsigned int k = kmin; k < kmax; ++k) {
331  vertices.se[k] += vertices.exp[k] * o_trk_sum_Z;
332  auto w = vertices.rho[k] * vertices.exp[k] * (o_trk_sum_Z * o_trk_dz2);
333  vertices.sw[k] += w;
334  vertices.swz[k] += w * tmp_trk_z;
335  vertices.swE[k] += w * vertices.exp_arg[k];
336  }
337  } else {
338  // same loop but without updating sWE
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  }
346  }
347  };
348 
349  if (updateTc) {
350  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
351  gvertices.se[ivertex] = 0.0;
352  gvertices.sw[ivertex] = 0.0;
353  gvertices.swz[ivertex] = 0.0;
354  gvertices.swE[ivertex] = 0.0;
355  }
356  } else {
357  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
358  gvertices.se[ivertex] = 0.0;
359  gvertices.sw[ivertex] = 0.0;
360  gvertices.swz[ivertex] = 0.0;
361  }
362  }
363 
364  // loop over tracks
365  for (auto itrack = 0U; itrack < nt; ++itrack) {
366  const unsigned int kmin = gtracks.kmin[itrack];
367  const unsigned int kmax = gtracks.kmax[itrack];
368 
369  kernel_calc_exp_arg_range(itrack, gtracks, gvertices, kmin, kmax);
370  local_exp_list_range(gvertices.exp_arg, gvertices.exp, kmin, kmax);
371  gtracks.sum_Z[itrack] = kernel_add_Z_range(gvertices, kmin, kmax);
372 
373  if (edm::isNotFinite(gtracks.sum_Z[itrack]))
374  gtracks.sum_Z[itrack] = 0.0;
375 
376  if (gtracks.sum_Z[itrack] > 1.e-100) {
377  kernel_calc_normalization_range(itrack, gtracks, gvertices, kmin, kmax);
378  }
379  }
380 
381  // (un-)apply the factor -beta which is needed in exp_arg, but not in swE
382  if (updateTc) {
383  auto obeta = -1. / beta;
384  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
385  gvertices.swE[ivertex] *= obeta;
386  }
387  }
388 
389  // now update z and rho
390  auto kernel_calc_z = [osumtkwt, nv](vertex_t& vertices) -> double {
391  double delta = 0;
392  // does not vectorize
393  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
394  if (vertices.sw[ivertex] > 0) {
395  auto znew = vertices.swz[ivertex] / vertices.sw[ivertex];
396  // prevents from vectorizing if
397  delta = max(std::abs(vertices.zvtx[ivertex] - znew), delta);
398  vertices.zvtx[ivertex] = znew;
399  }
400  }
401 
402  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex)
403  vertices.rho[ivertex] = vertices.rho[ivertex] * vertices.se[ivertex] * osumtkwt;
404 
405  return delta;
406  };
407 
408  double delta = kernel_calc_z(gvertices);
409 
410  // return how much the prototypes moved
411  return delta;
412 }
413 
415  double beta, track_t& tks, vertex_t& v, const double delta_max0, const double rho0) const {
416  unsigned int niter = 0;
417  double delta = 0;
418  double delta_max = delta_lowT_;
419 
420  if (convergence_mode_ == 0) {
421  delta_max = delta_max0;
422  } else if (convergence_mode_ == 1) {
423  delta_max = delta_lowT_ / sqrt(std::max(beta, 1.0));
424  }
425 
426  set_vtx_range(beta, tks, v);
427  double delta_sum_range = 0; // accumulate max(|delta-z|) as a lower bound
428  std::vector<double> z0 = v.zvtx_vec;
429 
430  while (niter++ < maxIterations_) {
431  delta = update(beta, tks, v, rho0);
432  delta_sum_range += delta;
433 
434  if (delta_sum_range > zrange_min_) {
435  for (unsigned int k = 0; k < v.getSize(); k++) {
436  if (std::abs(v.zvtx_vec[k] - z0[k]) > zrange_min_) {
437  set_vtx_range(beta, tks, v);
438  delta_sum_range = 0;
439  z0 = v.zvtx_vec;
440  break;
441  }
442  }
443  }
444 
445  if (delta < delta_max) {
446  break;
447  }
448  }
449 
450 #ifdef DEBUG
451  if (DEBUGLEVEL > 0) {
452  std::cout << "DAClusterizerInZ_vect.thermalize niter = " << niter << " at T = " << 1 / beta
453  << " nv = " << v.getSize() << std::endl;
454  if (DEBUGLEVEL > 2)
455  dump(beta, v, tks, 0, rho0);
456  }
457 #endif
458 
459  return niter;
460 }
461 
462 bool DAClusterizerInZ_vect::merge(vertex_t& y, track_t& tks, double& beta) const {
463  // merge clusters that collapsed or never separated,
464  // only merge if the estimated critical temperature of the merged vertex is below the current temperature
465  // return true if vertices were merged, false otherwise
466  const unsigned int nv = y.getSize();
467 
468  if (nv < 2)
469  return false;
470 
471  // merge the smallest distance clusters first
472  std::vector<std::pair<double, unsigned int> > critical;
473  for (unsigned int k = 0; (k + 1) < nv; k++) {
474  if (std::fabs(y.zvtx[k + 1] - y.zvtx[k]) < zmerge_) {
475  critical.push_back(make_pair(std::fabs(y.zvtx[k + 1] - y.zvtx[k]), k));
476  }
477  }
478  if (critical.empty())
479  return false;
480 
481  std::stable_sort(critical.begin(), critical.end(), std::less<std::pair<double, unsigned int> >());
482 
483  for (unsigned int ik = 0; ik < critical.size(); ik++) {
484  unsigned int k = critical[ik].second;
485  double rho = y.rho[k] + y.rho[k + 1];
486 
487 #ifdef DEBUG
488  assert((k + 1) < nv);
489  if (DEBUGLEVEL > 1) {
490  std::cout << "merging " << fixed << setprecision(4) << y.zvtx[k + 1] << " and " << y.zvtx[k]
491  << " sw = " << y.sw[k] + y.sw[k + 1] << std::endl;
492  }
493 #endif
494 
495  if (rho > 0) {
496  y.zvtx[k] = (y.rho[k] * y.zvtx[k] + y.rho[k + 1] * y.zvtx[k + 1]) / rho;
497  } else {
498  y.zvtx[k] = 0.5 * (y.zvtx[k] + y.zvtx[k + 1]);
499  }
500  y.rho[k] = rho;
501  y.sw[k] += y.sw[k + 1];
502  y.removeItem(k + 1, tks);
503  set_vtx_range(beta, tks, y);
504  y.extractRaw();
505  return true;
506  }
507 
508  return false;
509 }
510 
511 bool DAClusterizerInZ_vect::purge(vertex_t& y, track_t& tks, double& rho0, const double beta) const {
512  constexpr double eps = 1.e-100;
513  // eliminate clusters with only one significant/unique track
514  const unsigned int nv = y.getSize();
515  const unsigned int nt = tks.getSize();
516 
517  if (nv < 2)
518  return false;
519 
520  std::vector<double> sump_v(nv), arg_cache_v(nv), exp_cache_v(nv), pcut_cache_v(nv);
521  std::vector<int> nUnique_v(nv);
522  double* __restrict__ parg_cache;
523  double* __restrict__ pexp_cache;
524  double* __restrict__ ppcut_cache;
525  double* __restrict__ psump;
526  int* __restrict__ pnUnique;
527  int constexpr nunique_min_ = 2;
528 
529  set_vtx_range(beta, tks, y);
530 
531  parg_cache = arg_cache_v.data();
532  pexp_cache = exp_cache_v.data();
533  ppcut_cache = pcut_cache_v.data();
534  psump = sump_v.data();
535  pnUnique = nUnique_v.data();
536 
537  const auto rhoconst = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
538  for (unsigned int k = 0; k < nv; k++) {
539  const double pmax = y.rho[k] / (y.rho[k] + rhoconst);
540  ppcut_cache[k] = uniquetrkweight_ * pmax;
541  }
542 
543  for (unsigned int i = 0; i < nt; i++) {
544  const auto invZ = ((tks.sum_Z[i] > eps) && (tks.tkwt[i] > uniquetrkminp_)) ? 1. / tks.sum_Z[i] : 0.;
545  const auto track_z = tks.zpca[i];
546  const auto botrack_dz2 = -beta * tks.dz2[i];
547  const auto kmin = tks.kmin[i];
548  const auto kmax = tks.kmax[i];
549 
550  for (unsigned int k = kmin; k < kmax; k++) {
551  const auto mult_resz = track_z - y.zvtx[k];
552  parg_cache[k] = botrack_dz2 * (mult_resz * mult_resz);
553  }
554 
555  local_exp_list_range(parg_cache, pexp_cache, kmin, kmax);
556 
557  for (unsigned int k = kmin; k < kmax; k++) {
558  const double p = y.rho[k] * pexp_cache[k] * invZ;
559  psump[k] += p;
560  pnUnique[k] += (p > ppcut_cache[k]) ? 1 : 0;
561  }
562  }
563 
564  double sumpmin = nt;
565  unsigned int k0 = nv;
566  for (unsigned k = 0; k < nv; k++) {
567  if ((pnUnique[k] < nunique_min_) && (psump[k] < sumpmin)) {
568  sumpmin = psump[k];
569  k0 = k;
570  }
571  }
572 
573  if (k0 != nv) {
574 #ifdef DEBUG
575  assert(k0 < y.getSize());
576  if (DEBUGLEVEL > 1) {
577  std::cout << "eliminating prototype at " << std::setw(10) << std::setprecision(4) << y.zvtx[k0]
578  << " with sump=" << sumpmin << " rho*nt =" << y.rho[k0] * nt << endl;
579  }
580 #endif
581 
582  y.removeItem(k0, tks);
583  set_vtx_range(beta, tks, y);
584  return true;
585  } else {
586  return false;
587  }
588 }
589 
590 double DAClusterizerInZ_vect::beta0(double betamax, track_t const& tks, vertex_t const& y) const {
591  double T0 = 0; // max Tc for beta=0
592  // estimate critical temperature from beta=0 (T=inf)
593  const unsigned int nt = tks.getSize();
594  const unsigned int nv = y.getSize();
595 
596  for (unsigned int k = 0; k < nv; k++) {
597  // vertex fit at T=inf
598  double sumwz = 0;
599  double sumw = 0;
600  for (unsigned int i = 0; i < nt; i++) {
601  double w = tks.tkwt[i] * tks.dz2[i];
602  sumwz += w * tks.zpca[i];
603  sumw += w;
604  }
605 
606  y.zvtx[k] = sumwz / sumw;
607 
608  // estimate Tcrit
609  double a = 0, b = 0;
610  for (unsigned int i = 0; i < nt; i++) {
611  double dx = tks.zpca[i] - y.zvtx[k];
612  double w = tks.tkwt[i] * tks.dz2[i];
613  a += w * std::pow(dx, 2) * tks.dz2[i];
614  b += w;
615  }
616  double Tc = 2. * a / b; // the critical temperature of this vertex
617 
618  if (Tc > T0)
619  T0 = Tc;
620 
621  } // vertex loop (normally there should be only one vertex at beta=0)
622 
623 #ifdef DEBUG
624  if (DEBUGLEVEL > 0) {
625  std::cout << "DAClusterizerInZ_vect.beta0: Tc = " << T0 << std::endl;
626  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
627  std::cout << "DAClusterizerInZ_vect.beta0: nstep = " << coolingsteps << std::endl;
628  }
629 #endif
630 
631  if (T0 > 1. / betamax) {
632  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
633 
634  return betamax * std::pow(coolingFactor_, coolingsteps);
635  } else {
636  // ensure at least one annealing step
637  return betamax * coolingFactor_;
638  }
639 }
640 
641 bool DAClusterizerInZ_vect::split(const double beta, track_t& tks, vertex_t& y, double threshold) const {
642  // split only critical vertices (Tc >~ T=1/beta <==> beta*Tc>~1)
643  // returns true if at least one cluster was split
644 
645  // update the sums needed for Tc, rho0 is never non-zero while splitting is still active
646  update(beta, tks, y, 0., true); // the "true" flag enables Tc evaluation
647 
648  constexpr double epsilon = 1e-3; // minimum split size
649  unsigned int nv = y.getSize();
650 
651  // avoid left-right biases by splitting highest Tc first
652 
653  std::vector<std::pair<double, unsigned int> > critical;
654  for (unsigned int k = 0; k < nv; k++) {
655  double Tc = 2 * y.swE[k] / y.sw[k];
656  if (beta * Tc > threshold) {
657  critical.push_back(make_pair(Tc, k));
658  }
659  }
660  if (critical.empty())
661  return false;
662 
663  std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int> >());
664 
665  bool split = false;
666  const unsigned int nt = tks.getSize();
667 
668  for (unsigned int ic = 0; ic < critical.size(); ic++) {
669  unsigned int k = critical[ic].second;
670 
671  // estimate subcluster positions and weight
672  double p1 = 0, z1 = 0, w1 = 0;
673  double p2 = 0, z2 = 0, w2 = 0;
674  for (unsigned int i = 0; i < nt; i++) {
675  if (tks.sum_Z[i] > 1.e-100) {
676  // winner-takes-all, usually overestimates splitting
677  double tl = tks.zpca[i] < y.zvtx[k] ? 1. : 0.;
678  double tr = 1. - tl;
679 
680  // soften it, especially at low T
681  double arg = (tks.zpca[i] - y.zvtx[k]) * sqrt(beta * tks.dz2[i]);
682  if (std::fabs(arg) < 20) {
683  double t = local_exp(-arg);
684  tl = t / (t + 1.);
685  tr = 1 / (t + 1.);
686  }
687 
688  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];
689  double w = p * tks.dz2[i];
690  p1 += p * tl;
691  z1 += w * tl * tks.zpca[i];
692  w1 += w * tl;
693  p2 += p * tr;
694  z2 += w * tr * tks.zpca[i];
695  w2 += w * tr;
696  }
697  }
698 
699  if (w1 > 0) {
700  z1 = z1 / w1;
701  } else {
702  z1 = y.zvtx[k] - epsilon;
703  }
704  if (w2 > 0) {
705  z2 = z2 / w2;
706  } else {
707  z2 = y.zvtx[k] + epsilon;
708  }
709 
710  // reduce split size if there is not enough room
711  if ((k > 0) && (z1 < (0.6 * y.zvtx[k] + 0.4 * y.zvtx[k - 1]))) {
712  z1 = 0.6 * y.zvtx[k] + 0.4 * y.zvtx[k - 1];
713  }
714  if ((k + 1 < nv) && (z2 > (0.6 * y.zvtx[k] + 0.4 * y.zvtx[k + 1]))) {
715  z2 = 0.6 * y.zvtx[k] + 0.4 * y.zvtx[k + 1];
716  }
717 
718 #ifdef DEBUG
719  assert(k < nv);
720  if (DEBUGLEVEL > 1) {
721  if (std::fabs(y.zvtx[k] - zdumpcenter_) < zdumpwidth_) {
722  std::cout << " T= " << std::setw(8) << 1. / beta << " Tc= " << critical[ic].first << " splitting "
723  << std::fixed << std::setprecision(4) << y.zvtx[k] << " --> " << z1 << "," << z2 << " [" << p1
724  << "," << p2 << "]";
725  if (std::fabs(z2 - z1) > epsilon) {
726  std::cout << std::endl;
727  } else {
728  std::cout << " rejected " << std::endl;
729  }
730  }
731  }
732 #endif
733 
734  // split if the new subclusters are significantly separated
735  if ((z2 - z1) > epsilon) {
736  split = true;
737  double pk1 = p1 * y.rho[k] / (p1 + p2);
738  double pk2 = p2 * y.rho[k] / (p1 + p2);
739  y.zvtx[k] = z2;
740  y.rho[k] = pk2;
741  y.insertItem(k, z1, pk1, tks);
742  if (k == 0)
743  y.extractRaw();
744 
745  nv++;
746 
747  // adjust remaining pointers
748  for (unsigned int jc = ic; jc < critical.size(); jc++) {
749  if (critical[jc].second >= k) {
750  critical[jc].second++;
751  }
752  }
753  } else {
754 #ifdef DEBUG
755  std::cout << "warning ! split rejected, too small." << endl;
756 #endif
757  }
758  }
759 
760  return split;
761 }
762 
763 vector<TransientVertex> DAClusterizerInZ_vect::vertices(const vector<reco::TransientTrack>& tracks) const {
764  track_t&& tks = fill(tracks);
765  tks.extractRaw();
766 
767  unsigned int nt = tks.getSize();
768  double rho0 = 0.0; // start with no outlier rejection
769 
770  vector<TransientVertex> clusters;
771  if (tks.getSize() == 0)
772  return clusters;
773 
774  vertex_t y; // the vertex prototypes
775 
776  // initialize:single vertex at infinite temperature
777  y.addItem(0, 1.0);
778  clear_vtx_range(tks, y);
779 
780  // estimate first critical temperature
781  double beta = beta0(betamax_, tks, y);
782 #ifdef DEBUG
783  if (DEBUGLEVEL > 0)
784  std::cout << "Beta0 is " << beta << std::endl;
785 #endif
786 
787  thermalize(beta, tks, y, delta_highT_);
788 
789  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
790 
791  double betafreeze = betamax_ * sqrt(coolingFactor_);
792 
793  while (beta < betafreeze) {
794  while (merge(y, tks, beta)) {
795  update(beta, tks, y, rho0, false);
796  }
797  split(beta, tks, y);
798 
799  beta = beta / coolingFactor_;
800  thermalize(beta, tks, y, delta_highT_);
801  }
802 
803 #ifdef DEBUG
804  verify(y, tks);
805 
806  if (DEBUGLEVEL > 0) {
807  std::cout << "DAClusterizerInZ_vect::vertices :"
808  << "last round of splitting" << std::endl;
809  }
810 #endif
811 
812  set_vtx_range(beta, tks, y);
813  update(beta, tks, y, rho0, false);
814 
815  while (merge(y, tks, beta)) {
816  set_vtx_range(beta, tks, y);
817  update(beta, tks, y, rho0, false);
818  }
819 
820  unsigned int ntry = 0;
821  double threshold = 1.0;
822  while (split(beta, tks, y, threshold) && (ntry++ < 10)) {
823  thermalize(beta, tks, y, delta_highT_, rho0); // rho0 = 0. here
824  while (merge(y, tks, beta)) {
825  update(beta, tks, y, rho0, false);
826  }
827 
828  // relax splitting a bit to reduce multiple split-merge cycles of the same cluster
829  threshold *= 1.1;
830  }
831 
832 #ifdef DEBUG
833  verify(y, tks);
834  if (DEBUGLEVEL > 0) {
835  std::cout << "DAClusterizerInZ_vect::vertices :"
836  << "turning on outlier rejection at T=" << 1 / beta << std::endl;
837  }
838 #endif
839 
840  // switch on outlier rejection at T=Tmin
841  if (dzCutOff_ > 0) {
842  rho0 = y.getSize() > 1 ? 1. / y.getSize() : 1.;
843  for (unsigned int a = 0; a < 5; a++) {
844  update(beta, tks, y, a * rho0 / 5.); // adiabatic turn-on
845  }
846  }
847 
848  thermalize(beta, tks, y, delta_lowT_, rho0);
849 
850 #ifdef DEBUG
851  verify(y, tks);
852  if (DEBUGLEVEL > 0) {
853  std::cout << "DAClusterizerInZ_vect::vertices :"
854  << "merging with outlier rejection at T=" << 1 / beta << std::endl;
855  }
856  if (DEBUGLEVEL > 2)
857  dump(beta, y, tks, 2, rho0);
858 #endif
859 
860  // merge again (some cluster split by outliers collapse here)
861  while (merge(y, tks, beta)) {
862  set_vtx_range(beta, tks, y);
863  update(beta, tks, y, rho0, false);
864  }
865 
866 #ifdef DEBUG
867  verify(y, tks);
868  if (DEBUGLEVEL > 0) {
869  std::cout << "DAClusterizerInZ_vect::vertices :"
870  << "after merging with outlier rejection at T=" << 1 / beta << std::endl;
871  }
872  if (DEBUGLEVEL > 2)
873  dump(beta, y, tks, 2, rho0);
874 #endif
875 
876  // go down to the purging temperature (if it is lower than tmin)
877  while (beta < betapurge_) {
878  beta = min(beta / coolingFactor_, betapurge_);
879  thermalize(beta, tks, y, delta_lowT_, rho0);
880  }
881 
882 #ifdef DEBUG
883  verify(y, tks);
884  if (DEBUGLEVEL > 0) {
885  std::cout << "DAClusterizerInZ_vect::vertices :"
886  << "purging at T=" << 1 / beta << std::endl;
887  }
888 #endif
889 
890  // eliminate insignificant vertices, this is more restrictive at higher T
891  while (purge(y, tks, rho0, beta)) {
892  thermalize(beta, tks, y, delta_lowT_, rho0);
893  }
894 
895 #ifdef DEBUG
896  verify(y, tks);
897  if (DEBUGLEVEL > 0) {
898  std::cout << "DAClusterizerInZ_vect::vertices :"
899  << "last cooling T=" << 1 / beta << std::endl;
900  }
901 #endif
902 
903  // optionally cool some more without doing anything, to make the assignment harder
904  while (beta < betastop_) {
905  beta = min(beta / coolingFactor_, betastop_);
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 :"
913  << "stop cooling at T=" << 1 / beta << std::endl;
914  }
915  if (DEBUGLEVEL > 2)
916  dump(beta, y, tks, 2, rho0);
917 #endif
918 
919  // select significant tracks and use a TransientVertex as a container
920 
921  set_vtx_range(beta, tks, y);
922  const unsigned int nv = y.getSize();
923  for (unsigned int k = 0; k < nv; k++) {
924  if (edm::isNotFinite(y.rho[k]) || edm::isNotFinite(y.zvtx[k])) {
925  y.rho[k] = 0;
926  y.zvtx[k] = 0;
927  }
928  }
929 
930  const auto z_sum_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
931  std::vector<std::vector<unsigned int> > vtx_track_indices(nv);
932  for (unsigned int i = 0; i < nt; i++) {
933  const auto kmin = tks.kmin[i];
934  const auto kmax = tks.kmax[i];
935  for (auto k = kmin; k < kmax; k++) {
936  y.exp_arg[k] = -beta * Eik(tks.zpca[i], y.zvtx[k], tks.dz2[i]);
937  }
938 
939  local_exp_list_range(y.exp_arg, y.exp, kmin, kmax);
940 
941  tks.sum_Z[i] = z_sum_init;
942  for (auto k = kmin; k < kmax; k++) {
943  tks.sum_Z[i] += y.rho[k] * y.exp[k];
944  }
945  const double invZ = tks.sum_Z[i] > 1e-100 ? 1. / tks.sum_Z[i] : 0.0;
946 
947  for (auto k = kmin; k < kmax; k++) {
948  double p = y.rho[k] * y.exp[k] * invZ;
949  if (p > mintrkweight_) {
950  // assign track i -> vertex k (hard, mintrkweight_ should be >= 0.5 here
951  vtx_track_indices[k].push_back(i);
952  break;
953  }
954  }
955 
956  } // track loop
957 
958  GlobalError dummyError(0.01, 0, 0.01, 0., 0., 0.01);
959  for (unsigned int k = 0; k < nv; k++) {
960  if (!vtx_track_indices[k].empty()) {
961  GlobalPoint pos(0, 0, y.zvtx[k]);
962  vector<reco::TransientTrack> vertexTracks;
963  for (auto i : vtx_track_indices[k]) {
964  vertexTracks.push_back(*(tks.tt[i]));
965  }
966  TransientVertex v(pos, dummyError, vertexTracks, 0);
967  clusters.push_back(v);
968  }
969  }
970 
971  return clusters;
972 }
973 
974 vector<vector<reco::TransientTrack> > DAClusterizerInZ_vect::clusterize(
975  const vector<reco::TransientTrack>& tracks) const {
976  vector<vector<reco::TransientTrack> > clusters;
977  vector<TransientVertex>&& pv = vertices(tracks);
978 
979 #ifdef DEBUG
980  if (DEBUGLEVEL > 0) {
981  std::cout << "###################################################" << endl;
982  std::cout << "# vectorized DAClusterizerInZ_vect::clusterize nt=" << tracks.size() << endl;
983  std::cout << "# DAClusterizerInZ_vect::clusterize pv.size=" << pv.size() << endl;
984  std::cout << "###################################################" << endl;
985  }
986 #endif
987 
988  if (pv.empty()) {
989  return clusters;
990  }
991 
992  // fill into clusters and merge
993  vector<reco::TransientTrack> aCluster = pv.begin()->originalTracks();
994 
995  for (auto k = pv.begin() + 1; k != pv.end(); k++) {
996  if (std::abs(k->position().z() - (k - 1)->position().z()) > (2 * vertexSize_)) {
997  // close a cluster
998  if (aCluster.size() > 1) {
999  clusters.push_back(aCluster);
1000  }
1001 #ifdef DEBUG
1002  else {
1003  std::cout << " one track cluster at " << k->position().z() << " suppressed" << std::endl;
1004  }
1005 #endif
1006  aCluster.clear();
1007  }
1008  for (unsigned int i = 0; i < k->originalTracks().size(); i++) {
1009  aCluster.push_back(k->originalTracks()[i]);
1010  }
1011  }
1012  clusters.emplace_back(std::move(aCluster));
1013 
1014  return clusters;
1015 }
1016 
1018  const double beta, const vertex_t& y, const track_t& tks, const int verbosity, const double rho0) const {
1019 #ifdef DEBUG
1020  const unsigned int nv = y.getSize();
1021  const unsigned int nt = tks.getSize();
1022  // make a local copies of clusters and tracks to update Tc [ = y_local.swE / y_local.sw ]
1023  vertex_t y_local = y;
1024  track_t tks_local = tks;
1025  update(beta, tks_local, y_local, rho0, true);
1026 
1027  std::vector<unsigned int> iz;
1028  for (unsigned int j = 0; j < nt; j++) {
1029  iz.push_back(j);
1030  }
1031  std::sort(iz.begin(), iz.end(), [tks](unsigned int a, unsigned int b) { return tks.zpca[a] < tks.zpca[b]; });
1032  std::cout << std::endl;
1033  std::cout << "-----DAClusterizerInZ::dump ----" << nv << " clusters " << std::endl;
1034  std::cout << " ";
1035  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1036  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1037  std::cout << " " << setw(3) << ivertex << " ";
1038  }
1039  }
1040  std::cout << endl;
1041  std::cout << " z= ";
1042  std::cout << setprecision(4);
1043  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1044  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1045  std::cout << setw(8) << fixed << y.zvtx[ivertex];
1046  }
1047  }
1048  std::cout << endl
1049  << "T=" << setw(15) << 1. / beta << " Tmin =" << setw(10) << 1. / betamax_
1050  << " Tc= ";
1051  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1052  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1053  double Tc = 2 * y_local.swE[ivertex] / y_local.sw[ivertex];
1054  std::cout << setw(8) << fixed << setprecision(1) << Tc;
1055  }
1056  }
1057  std::cout << endl;
1058 
1059  std::cout << " pk= ";
1060  double sumpk = 0;
1061  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1062  sumpk += y.rho[ivertex];
1063  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1064  continue;
1065  std::cout << setw(8) << setprecision(4) << fixed << y.rho[ivertex];
1066  }
1067  std::cout << endl;
1068 
1069  std::cout << " nt= ";
1070  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1071  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1072  continue;
1073  std::cout << setw(8) << setprecision(1) << fixed << y.rho[ivertex] * nt;
1074  }
1075  std::cout << endl;
1076 
1077  if (verbosity > 0) {
1078  double E = 0, F = 0;
1079  std::cout << endl;
1080  std::cout << "---- z +/- dz ip +/-dip pt phi eta weights ----" << endl;
1081  std::cout << setprecision(4);
1082  for (unsigned int i0 = 0; i0 < nt; i0++) {
1083  unsigned int i = iz[i0];
1084  if (tks.sum_Z[i] > 0) {
1085  F -= std::log(tks.sum_Z[i]) / beta;
1086  }
1087  double tz = tks.zpca[i];
1088 
1089  if (std::fabs(tz - zdumpcenter_) > zdumpwidth_)
1090  continue;
1091  std::cout << setw(4) << i << ")" << setw(8) << fixed << setprecision(4) << tz << " +/-" << setw(6)
1092  << sqrt(1. / tks.dz2[i]);
1093  if ((tks.tt[i] == nullptr)) {
1094  std::cout << " effective track ";
1095  } else {
1096  if (tks.tt[i]->track().quality(reco::TrackBase::highPurity)) {
1097  std::cout << " *";
1098  } else {
1099  std::cout << " ";
1100  }
1101  if (tks.tt[i]->track().hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1)) {
1102  std::cout << "+";
1103  } else {
1104  std::cout << "-";
1105  }
1106  std::cout << setw(1)
1107  << tks.tt[i]
1108  ->track()
1109  .hitPattern()
1110  .pixelBarrelLayersWithMeasurement(); // see DataFormats/TrackReco/interface/HitPattern.h
1111  std::cout << setw(1) << tks.tt[i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
1112  std::cout << setw(1) << hex
1113  << tks.tt[i]->track().hitPattern().trackerLayersWithMeasurement() -
1114  tks.tt[i]->track().hitPattern().pixelLayersWithMeasurement()
1115  << dec;
1116  std::cout << "=" << setw(1) << hex
1117  << tks.tt[i]->track().hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS) << dec;
1118 
1119  Measurement1D IP = tks.tt[i]->stateAtBeamLine().transverseImpactParameter();
1120  std::cout << setw(8) << IP.value() << "+/-" << setw(6) << IP.error();
1121  std::cout << " " << setw(6) << setprecision(2) << tks.tt[i]->track().pt() * tks.tt[i]->track().charge();
1122  std::cout << " " << setw(5) << setprecision(2) << tks.tt[i]->track().phi() << " " << setw(5) << setprecision(2)
1123  << tks.tt[i]->track().eta();
1124  } // not a dummy track
1125 
1126  double sump = 0.;
1127  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1128  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1129  continue;
1130 
1131  if ((tks.tkwt[i] > 0) && (tks.sum_Z[i] > 0)) {
1132  //double p=pik(beta,tks[i],*k);
1133  double p = y.rho[ivertex] * local_exp(-beta * Eik(tks.zpca[i], y.zvtx[ivertex], tks.dz2[i])) / tks.sum_Z[i];
1134  if (p > 0.0001) {
1135  std::cout << setw(8) << setprecision(3) << p;
1136  } else {
1137  std::cout << " . ";
1138  }
1139  E += p * Eik(tks.zpca[i], y.zvtx[ivertex], tks.dz2[i]);
1140  sump += p;
1141  } else {
1142  std::cout << " ";
1143  }
1144  }
1145  std::cout << " ( " << std::setw(3) << tks.kmin[i] << "," << std::setw(3) << tks.kmax[i] - 1 << " ) ";
1146  std::cout << endl;
1147  }
1148  std::cout << " ";
1149  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1150  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1151  std::cout << " " << setw(3) << ivertex << " ";
1152  }
1153  }
1154  std::cout << endl;
1155  std::cout << " z= ";
1156  std::cout << setprecision(4);
1157  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1158  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1159  std::cout << setw(8) << fixed << y.zvtx[ivertex];
1160  }
1161  }
1162  std::cout << endl;
1163  std::cout << endl
1164  << "T=" << 1 / beta << " E=" << E << " n=" << y.getSize() << " F= " << F << endl
1165  << "----------" << endl;
1166  }
1167 #endif
1168 }
1169 
1171  desc.addUntracked<double>("zdumpcenter", 0.);
1172  desc.addUntracked<double>("zdumpwidth", 20.);
1173  desc.add<double>("d0CutOff", 3.0);
1174  desc.add<double>("Tmin", 2.0);
1175  desc.add<double>("delta_lowT", 0.001);
1176  desc.add<double>("zmerge", 0.01);
1177  desc.add<double>("dzCutOff", 3.0);
1178  desc.add<double>("Tpurge", 2.0);
1179  desc.add<int>("convergence_mode", 0);
1180  desc.add<double>("delta_highT", 0.01);
1181  desc.add<double>("Tstop", 0.5);
1182  desc.add<double>("coolingFactor", 0.6);
1183  desc.add<double>("vertexSize", 0.006);
1184  desc.add<double>("uniquetrkweight", 0.8);
1185  desc.add<double>("uniquetrkminp", 0.0);
1186  desc.add<double>("zrange", 4.0);
1187 }
HIPAlignmentAlgorithm_cfi.verbosity
verbosity
Definition: HIPAlignmentAlgorithm_cfi.py:7
reco::HitPattern::MISSING_OUTER_HITS
Definition: HitPattern.h:155
alignBH_cfg.fixed
fixed
Definition: alignBH_cfg.py:54
DAClusterizerInZ_vect::fill
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
Definition: DAClusterizerInZ_vect.cc:170
Measurement1D
Definition: Measurement1D.h:11
w2
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
mps_fire.i
i
Definition: mps_fire.py:428
HLT_FULL_cff.atIP
atIP
Definition: HLT_FULL_cff.py:46147
pwdgSkimBPark_cfi.beamSpot
beamSpot
Definition: pwdgSkimBPark_cfi.py:5
MessageLogger.h
DAClusterizerInZ_vect.h
TkClusParameters_cff.zrange
zrange
Definition: TkClusParameters_cff.py:7
reco::ParticleMasses::k0
const double k0
Definition: ParticleMasses.h:7
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
nt
int nt
Definition: AMPTWrapper.h:42
DAClusterizerInZ_vect::track_t::sum_Z
double *__restrict__ sum_Z
Definition: DAClusterizerInZ_vect.h:72
min
T min(T a, T b)
Definition: MathUtil.h:58
HLT_FULL_cff.beta
beta
Definition: HLT_FULL_cff.py:8686
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
DAClusterizerInZ_vect::thermalize
unsigned int thermalize(double beta, track_t &gtracks, vertex_t &gvertices, const double delta_max, const double rho0=0.) const
Definition: DAClusterizerInZ_vect.cc:414
DAClusterizerInZ_vect::track_t::tkwt
double *__restrict__ tkwt
Definition: DAClusterizerInZ_vect.h:71
gpuVertexFinder::eps
WorkSpace int float eps
Definition: gpuClusterTracksDBSCAN.h:18
HLT_FULL_cff.Tmin
Tmin
Definition: HLT_FULL_cff.py:34844
gather_cfg.cout
cout
Definition: gather_cfg.py:144
pos
Definition: PixelAliasList.h:18
AlignmentPI::t_z
Definition: AlignmentPayloadInspectorHelper.h:36
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
DAClusterizerInZ_vect::track_t::addItemSorted
void addItemSorted(double new_zpca, double new_dz2, const reco::TransientTrack *new_tt, double new_tkwt)
Definition: DAClusterizerInZ_vect.h:37
Measurement1D.h
DAClusterizerInZ_vect::DAClusterizerInZ_vect
DAClusterizerInZ_vect(const edm::ParameterSet &conf)
Definition: DAClusterizerInZ_vect.cc:20
cms::cuda::assert
assert(be >=bs)
edm::second
U second(std::pair< T, U > const &p)
Definition: ParameterSet.cc:222
DAClusterizerInZ_vect::track_t::tt
std::vector< const reco::TransientTrack * > tt
Definition: DAClusterizerInZ_vect.h:33
HLT_FULL_cff.Tstop
Tstop
Definition: HLT_FULL_cff.py:34838
edm::ParameterSet::getUntrackedParameter
T getUntrackedParameter(std::string const &, T const &) const
DAClusterizerInZ_vect::set_vtx_range
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const
Definition: DAClusterizerInZ_vect.cc:220
findQualityFiles.v
v
Definition: findQualityFiles.py:179
DAClusterizerInZ_vect::track_t::extractRaw
void extractRaw()
Definition: DAClusterizerInZ_vect.h:61
DAClusterizerInZ_vect::verify
void verify(const vertex_t &v, const track_t &tks, unsigned int nv=999999, unsigned int nt=999999) const
Definition: DAClusterizerInZ_vect.cc:108
DAClusterizerInZ_vect::vertex_t::zvtx
double *__restrict__ zvtx
Definition: DAClusterizerInZ_vect.h:151
SiStripMonitorCluster_cfi.zmin
zmin
Definition: SiStripMonitorCluster_cfi.py:200
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
DAClusterizerInZ_vect::track_t::kmax
std::vector< unsigned int > kmax
Definition: DAClusterizerInZ_vect.h:32
DAClusterizerInZ_vect::track_t::kmin
std::vector< unsigned int > kmin
Definition: DAClusterizerInZ_vect.h:31
geometryDiff.epsilon
int epsilon
Definition: geometryDiff.py:26
F
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163
testProducerWithPsetDescEmpty_cfi.z2
z2
Definition: testProducerWithPsetDescEmpty_cfi.py:41
DAClusterizerInZ_vect::split
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const
Definition: DAClusterizerInZ_vect.cc:641
GeomDetEnumerators::PixelBarrel
Definition: GeomDetEnumerators.h:11
DAClusterizerInZ_vect::vertex_t::se
double *__restrict__ se
Definition: DAClusterizerInZ_vect.h:157
DAClusterizerInZ_vect::track_t::zpca
double *__restrict__ zpca
Definition: DAClusterizerInZ_vect.h:69
SiStripMonitorCluster_cfi.zmax
zmax
Definition: SiStripMonitorCluster_cfi.py:201
submitPVValidationJobs.split
def split(sequence, size)
Definition: submitPVValidationJobs.py:352
w
const double w
Definition: UKUtility.cc:23
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
reco::BeamSpot
Definition: BeamSpot.h:21
p2
double p2[4]
Definition: TauolaWrapper.h:90
DAClusterizerInZ_vect::vertex_t::swE
double *__restrict__ swE
Definition: DAClusterizerInZ_vect.h:158
beamspot
Definition: BeamSpotWrite2Txt.h:8
DAClusterizerInZ_vect::track_t::dz2_vec
std::vector< double > dz2_vec
Definition: DAClusterizerInZ_vect.h:28
HLTMuonOfflineAnalyzer_cfi.z0
z0
Definition: HLTMuonOfflineAnalyzer_cfi.py:98
dqmdumpme.k
k
Definition: dqmdumpme.py:60
DAClusterizerInZ_vect::track_t::sum_Z_vec
std::vector< double > sum_Z_vec
Definition: DAClusterizerInZ_vect.h:29
Point3DBase< float, GlobalTag >
b
double b
Definition: hdecay.h:118
mitigatedMETSequence_cff.U
U
Definition: mitigatedMETSequence_cff.py:36
DAClusterizerInZ_vect::beta0
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
Definition: DAClusterizerInZ_vect.cc:590
DAClusterizerInZ_vect::track_t::dz2
double *__restrict__ dz2
Definition: DAClusterizerInZ_vect.h:70
DAClusterizerInZ_vect::vertex_t::sw
double *__restrict__ sw
Definition: DAClusterizerInZ_vect.h:155
DAClusterizerInZ_vect::track_t
Definition: DAClusterizerInZ_vect.h:26
bsc_activity_cfg.clusters
clusters
Definition: bsc_activity_cfg.py:36
DAClusterizerInZ_vect::track_t::tkwt_vec
std::vector< double > tkwt_vec
Definition: DAClusterizerInZ_vect.h:30
DAClusterizerInZ_vect::merge
bool merge(vertex_t &y, track_t &tks, double &beta) const
Definition: DAClusterizerInZ_vect.cc:462
listHistos.IP
IP
Definition: listHistos.py:76
ntuplemaker.fill
fill
Definition: ntuplemaker.py:304
edm::ParameterSet
Definition: ParameterSet.h:47
a
double a
Definition: hdecay.h:119
tracks
const uint32_t *__restrict__ const HitContainer *__restrict__ TkSoA *__restrict__ tracks
Definition: CAHitNtupletGeneratorKernelsImpl.h:159
SiStripPI::max
Definition: SiStripPayloadInspectorHelper.h:169
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
HLT_FULL_cff.Tpurge
Tpurge
Definition: HLT_FULL_cff.py:34843
dumpMFGeometry_cfg.delta
delta
Definition: dumpMFGeometry_cfg.py:25
DAClusterizerInZ_vect::vertex_t::swz
double *__restrict__ swz
Definition: DAClusterizerInZ_vect.h:156
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
FrontierConditions_GlobalTag_cff.dump
dump
Definition: FrontierConditions_GlobalTag_cff.py:12
createfilelist.int
int
Definition: createfilelist.py:10
DAClusterizerInZ_vect::vertex_t::getSize
unsigned int getSize() const
Definition: DAClusterizerInZ_vect.h:87
MetAnalyzer.pv
def pv(vc)
Definition: MetAnalyzer.py:7
DAClusterizerInZ_vect::track_t::getSize
unsigned int getSize() const
Definition: DAClusterizerInZ_vect.h:58
hlt_jetmet_dqm_QT_fromfile_cfg.critical
critical
Definition: hlt_jetmet_dqm_QT_fromfile_cfg.py:109
p1
double p1[4]
Definition: TauolaWrapper.h:89
GlobalErrorBase< double, ErrorMatrixTag >
DAClusterizerInZ_vect::clear_vtx_range
void clear_vtx_range(track_t &gtracks, vertex_t &gvertices) const
Definition: DAClusterizerInZ_vect.cc:269
DAClusterizerInZ_vect::fillPSetDescription
static void fillPSetDescription(edm::ParameterSetDescription &desc)
Definition: DAClusterizerInZ_vect.cc:1170
TransientVertex
Definition: TransientVertex.h:18
cms::cuda::for
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
Definition: HistoContainer.h:27
DAClusterizerInZ_vect::clusterize
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const override
Definition: DAClusterizerInZ_vect.cc:974
DAClusterizerInZ_vect::vertex_t::exp_arg
double *__restrict__ exp_arg
Definition: DAClusterizerInZ_vect.h:153
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
isFinite.h
DAClusterizerInZ_vect::track_t::zpca_vec
std::vector< double > zpca_vec
Definition: DAClusterizerInZ_vect.h:27
relativeConstraints.empty
bool empty
Definition: relativeConstraints.py:46
DAClusterizerInZ_vect::purge
bool purge(vertex_t &, track_t &, double &, const double) const
Definition: DAClusterizerInZ_vect.cc:511
DAClusterizerInZ_vect::update
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
Definition: DAClusterizerInZ_vect.cc:278
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
DAClusterizerInZ_vect::vertex_t::exp
double *__restrict__ exp
Definition: DAClusterizerInZ_vect.h:154
VertexException.h
dqm-mbProfile.log
log
Definition: dqm-mbProfile.py:17
funct::arg
A arg
Definition: Factorize.h:31
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
DAClusterizerInZ_vect::vertices
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const
Definition: DAClusterizerInZ_vect.cc:763
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
DAClusterizerInZ_vect::track_t::osumtkwt
double osumtkwt
Definition: DAClusterizerInZ_vect.h:35
remoteMonitoring_LED_IterMethod_cfg.threshold
threshold
Definition: remoteMonitoring_LED_IterMethod_cfg.py:430
MatrixUtil.merge
def merge(dictlist, TELL=False)
Definition: MatrixUtil.py:201
submitPVValidationJobs.t
string t
Definition: submitPVValidationJobs.py:644
DAClusterizerInZ_vect::dump
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0, const double rho0=0.) const
Definition: DAClusterizerInZ_vect.cc:1017
TauDecayModes.dec
dec
Definition: TauDecayModes.py:143
PVValHelper::dx
Definition: PVValidationHelpers.h:49
DAClusterizerInZ_vect::vertex_t
Definition: DAClusterizerInZ_vect.h:76
update
#define update(a, b)
Definition: TrackClassifier.cc:10
pwdgSkimBPark_cfi.vertices
vertices
Definition: pwdgSkimBPark_cfi.py:7
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37
reco::TrackBase::highPurity
Definition: TrackBase.h:154