CMS 3D CMS Logo

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