CMS 3D CMS Logo

DAClusterizerInZT_vect.cc
Go to the documentation of this file.
5 
6 #include <cmath>
7 #include <cassert>
8 #include <limits>
9 #include <iomanip>
11 #include "vdt/vdtMath.h"
12 
13 using namespace std;
14 //#define VI_DEBUG
15 
17 
18  // hardcoded parameters
19  maxIterations_ = 100;
20  mintrkweight_ = 0.5; // conf.getParameter<double>("mintrkweight");
21 
22 
23  // configurable debug outptut debug output
24  verbose_ = conf.getUntrackedParameter<bool> ("verbose", false);
25  zdumpcenter_ = conf.getUntrackedParameter<double> ("zdumpcenter", 0.);
26  zdumpwidth_ = conf.getUntrackedParameter<double> ("zdumpwidth", 20.);
27 
28  // configurable parameters
29  double minT = conf.getParameter<double> ("Tmin");
30  double purgeT = conf.getParameter<double> ("Tpurge");
31  double stopT = conf.getParameter<double> ("Tstop");
32  vertexSize_ = conf.getParameter<double> ("vertexSize");
33  vertexSizeTime_ = conf.getParameter<double> ("vertexSizeTime");
34  coolingFactor_ = conf.getParameter<double> ("coolingFactor");
35  useTc_=true;
36  if(coolingFactor_<0){
37  coolingFactor_=-coolingFactor_;
38  useTc_=false;
39  }
40  d0CutOff_ = conf.getParameter<double> ("d0CutOff");
41  dzCutOff_ = conf.getParameter<double> ("dzCutOff");
42  dtCutOff_ = conf.getParameter<double> ("dtCutOff");
43  t0Max_ = conf.getParameter<double> ("t0Max");
44  uniquetrkweight_ = conf.getParameter<double>("uniquetrkweight");
45  zmerge_ = conf.getParameter<double>("zmerge");
46  tmerge_ = conf.getParameter<double>("tmerge");
47 
48 #ifdef VI_DEBUG
49  if(verbose_){
50  std::cout << "DAClusterizerinZT_vect: mintrkweight = " << mintrkweight_ << std::endl;
51  std::cout << "DAClusterizerinZT_vect: uniquetrkweight = " << uniquetrkweight_ << std::endl;
52  std::cout << "DAClusterizerinZT_vect: zmerge = " << zmerge_ << std::endl;
53  std::cout << "DAClusterizerinZT_vect: tmerge = " << tmerge_ << std::endl;
54  std::cout << "DAClusterizerinZT_vect: Tmin = " << minT << std::endl;
55  std::cout << "DAClusterizerinZT_vect: Tpurge = " << purgeT << std::endl;
56  std::cout << "DAClusterizerinZT_vect: Tstop = " << stopT << std::endl;
57  std::cout << "DAClusterizerinZT_vect: vertexSize = " << vertexSize_ << std::endl;
58  std::cout << "DAClusterizerinZT_vect: vertexSizeTime = " << vertexSizeTime_ << std::endl;
59  std::cout << "DAClusterizerinZT_vect: coolingFactor = " << coolingFactor_ << std::endl;
60  std::cout << "DAClusterizerinZT_vect: d0CutOff = " << d0CutOff_ << std::endl;
61  std::cout << "DAClusterizerinZT_vect: dzCutOff = " << dzCutOff_ << std::endl;
62  std::cout << "DAClusterizerinZT_vect: dtCutoff = " << dtCutOff_ << std::endl;
63  }
64 #endif
65 
66 
67  if (minT == 0) {
68  edm::LogWarning("DAClusterizerinZT_vectorized") << "DAClusterizerInZT: invalid Tmin" << minT
69  << " reset do default " << 1. / betamax_;
70  } else {
71  betamax_ = 1. / minT;
72  }
73 
74 
75  if ((purgeT > minT) || (purgeT == 0)) {
76  edm::LogWarning("DAClusterizerinZT_vectorized") << "DAClusterizerInZT: invalid Tpurge" << purgeT
77  << " set to " << minT;
78  purgeT = minT;
79  }
80  betapurge_ = 1./purgeT;
81 
82 
83  if ((stopT > purgeT) || (stopT == 0)) {
84  edm::LogWarning("DAClusterizerinZT_vectorized") << "DAClusterizerInZT: invalid Tstop" << stopT
85  << " set to " << max(1., purgeT);
86  stopT = max(1., purgeT) ;
87  }
88  betastop_ = 1./stopT;
89 
90 }
91 
92 
93 namespace {
94  inline double local_exp( double const& inp) {
95  return vdt::fast_exp( inp );
96  }
97 
98  inline void local_exp_list( double const *arg_inp,
99  double *arg_out,
100  const unsigned arg_arr_size) {
101  for(unsigned i=0; i!=arg_arr_size; ++i ) arg_out[i]=vdt::fast_exp(arg_inp[i]);
102  }
103 
104 }
105 
106 //todo: use r-value possibility of c++11 here
108 DAClusterizerInZT_vect::fill(const vector<reco::TransientTrack> & tracks) const {
109 
110  // prepare track data for clustering
111  track_t tks;
112  for( const auto& tk : tracks ) {
113  if (!tk.isValid()) continue;
114  double t_pi=1.;
115  double t_z = tk.stateAtBeamLine().trackStateAtPCA().position().z();
116  double t_t = tk.timeExt();
117  if (std::fabs(t_z) > 1000.) continue;
118  if (std::abs(t_t) > t0Max_) continue;
119  auto const & t_mom = tk.stateAtBeamLine().trackStateAtPCA().momentum();
120  // get the beam-spot
121  reco::BeamSpot beamspot = tk.stateAtBeamLine().beamSpot();
122  double t_dz2 =
123  std::pow(tk.track().dzError(), 2) // track errror
124  + (std::pow(beamspot.BeamWidthX()*t_mom.x(),2)+std::pow(beamspot.BeamWidthY()*t_mom.y(),2))*std::pow(t_mom.z(),2)/std::pow(t_mom.perp2(),2) // beam spot width
125  + std::pow(vertexSize_, 2); // intrinsic vertex size, safer for outliers and short lived decays
126  t_dz2 = 1./ t_dz2;
127  if (edm::isNotFinite(t_dz2) || t_dz2 < std::numeric_limits<double>::min() ) continue;
128 
129 
130  double t_dt2 =std::pow(tk.dtErrorExt(),2.) + std::pow(vertexSizeTime_,2.); // the ~injected~ timing error, need to add a small minimum vertex size in time
131  if (tk.dtErrorExt()>0.3){
132  t_dt2 = 0; // tracks with no time measurement
133  }else{
134  t_dt2 = 1./t_dt2;
136  }
137 
138 
139  if (d0CutOff_ > 0) {
141  tk.stateAtBeamLine().transverseImpactParameter();// error contains beamspot
142  t_pi = 1. / (1. + local_exp(std::pow(atIP.value() / atIP.error(), 2) - std::pow(d0CutOff_, 2))); // reduce weight for high ip tracks
143  if (edm::isNotFinite(t_pi) || t_pi < std::numeric_limits<double>::epsilon()) continue; // usually is > 0.99
144  }
145  LogTrace("DAClusterizerinZT_vectorized") << t_z << ' ' << t_t <<' '<< t_dz2 << ' ' << t_dt2 <<' '<< t_pi;
146  tks.addItem(t_z, t_t, t_dz2, t_dt2, &tk, t_pi);
147  }
148  tks.extractRaw();
149 
150 #ifdef VI_DEBUG
151  if (verbose_) {
152  std::cout << "Track count " << tks.getSize() << std::endl;
153  }
154 #endif
155 
156  return tks;
157 }
158 
159 
160 namespace {
161  inline
162  double Eik(double t_z, double k_z, double t_dz2, double t_t, double k_t, double t_dt2) {
163  return std::pow(t_z - k_z, 2) * t_dz2 + std::pow(t_t - k_t,2) * t_dt2;
164  }
165 }
166 
167 double DAClusterizerInZT_vect::update(double beta, track_t & gtracks,
168  vertex_t & gvertices, bool useRho0, const double & rho0) const {
169 
170  //update weights and vertex positions
171  // mass constrained annealing without noise
172  // returns the squared sum of changes of vertex positions
173 
174  const unsigned int nt = gtracks.getSize();
175  const unsigned int nv = gvertices.getSize();
176 
177  //initialize sums
178  double sumpi = 0.;
179 
180  // to return how much the prototype moved
181  double delta = 0.;
182 
183 
184  // intial value of a sum
185  double Z_init = 0;
186  // independpent of loop
187  if ( useRho0 )
188  {
189  Z_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_); // cut-off
190  }
191 
192  // define kernels
193  auto kernel_calc_exp_arg = [ beta, nv ] ( const unsigned int itrack,
194  track_t const& tracks,
195  vertex_t const& vertices ) {
196 
197  const auto track_z = tracks.z_[itrack];
198  const auto track_t = tracks.t_[itrack];
199  const auto botrack_dz2 = -beta*tracks.dz2_[itrack];
200  const auto botrack_dt2 = -beta*tracks.dt2_[itrack];
201 
202  // auto-vectorized
203  for ( unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
204  const auto mult_resz = track_z - vertices.z_[ivertex];
205  const auto mult_rest = track_t - vertices.t_[ivertex];
206  vertices.ei_cache_[ivertex] = botrack_dz2 * ( mult_resz * mult_resz ) + botrack_dt2 * ( mult_rest * mult_rest );
207  }
208  };
209 
210  auto kernel_add_Z = [ nv, Z_init ] (vertex_t const& vertices) -> double
211  {
212  double ZTemp = Z_init;
213  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
214  ZTemp += vertices.pk_[ivertex] * vertices.ei_[ivertex];
215  }
216  return ZTemp;
217  };
218 
219  auto kernel_calc_normalization = [ nv ] (const unsigned int track_num,
220  track_t & tks_vec,
221  vertex_t & y_vec ) {
222  auto tmp_trk_pi = tks_vec.pi_[track_num];
223  auto o_trk_Z_sum = 1./tks_vec.Z_sum_[track_num];
224  auto o_trk_err_z = tks_vec.dz2_[track_num];
225  auto o_trk_err_t = tks_vec.dt2_[track_num];
226  auto tmp_trk_z = tks_vec.z_[track_num];
227  auto tmp_trk_t = tks_vec.t_[track_num];
228 
229 
230  // auto-vectorized
231  for (unsigned int k = 0; k < nv; ++k) {
232  // parens are important for numerical stability
233  y_vec.se_[k] += tmp_trk_pi*( y_vec.ei_[k] * o_trk_Z_sum );
234  const auto w = tmp_trk_pi * (y_vec.pk_[k] * y_vec.ei_[k] * o_trk_Z_sum); // p_{ik}
235  const auto wz = w * o_trk_err_z;
236  const auto wt = w * o_trk_err_t;
237  y_vec.nuz_[k] += wz;
238  y_vec.nut_[k] += wt;
239  y_vec.swz_[k] += wz * tmp_trk_z;
240  y_vec.swt_[k] += wt * tmp_trk_t;
241  /* this is really only needed when we want to get Tc too, mayb better to do it elsewhere? */
242  const auto dsz = (tmp_trk_z - y_vec.z[k]) * o_trk_err_z;
243  const auto dst = (tmp_trk_t - y_vec.t[k]) * o_trk_err_t;
244  y_vec.szz_[k] += w * dsz * dsz;
245  y_vec.stt_[k] += w * dst * dst;
246  y_vec.szt_[k] += w * dsz * dst;
247  }
248  };
249 
250 
251  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
252  gvertices.se_[ivertex] = 0.0;
253  gvertices.nuz_[ivertex] = 0.0;
254  gvertices.nut_[ivertex] = 0.0;
255  gvertices.swz_[ivertex] = 0.0;
256  gvertices.swt_[ivertex] = 0.0;
257  gvertices.szz_[ivertex] = 0.0;
258  gvertices.stt_[ivertex] = 0.0;
259  gvertices.szt_[ivertex] = 0.0;
260  }
261 
262 
263  // loop over tracks
264  for (auto itrack = 0U; itrack < nt; ++itrack) {
265  kernel_calc_exp_arg(itrack, gtracks, gvertices);
266  local_exp_list(gvertices.ei_cache_, gvertices.ei_, nv);
267 
268  gtracks.Z_sum_[itrack] = kernel_add_Z(gvertices);
269  if (edm::isNotFinite(gtracks.Z_sum_[itrack])) gtracks.Z_sum_[itrack] = 0.0;
270  // used in the next major loop to follow
271  sumpi += gtracks.pi_[itrack];
272 
273  if (gtracks.Z_sum_[itrack] > 1.e-100){
274  kernel_calc_normalization(itrack, gtracks, gvertices);
275  }
276  }
277 
278  // now update z, t, and pk
279  auto kernel_calc_zt = [ sumpi, nv
280 #ifdef VI_DEBUG
281  , this
282 #endif
283  ] (vertex_t & vertices ) -> double {
284 
285  double delta=0;
286 
287  // does not vectorizes
288  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex ) {
289 
290  if (vertices.nuz_[ivertex] > 0.) {
291  auto znew = vertices.swz_[ ivertex ] / vertices.nuz_[ ivertex ];
292  // prevents from vectorizing if
293  delta += std::pow( vertices.z_[ ivertex ] - znew, 2 );
294  vertices.z_[ ivertex ] = znew;
295  }
296 
297  if(vertices.nut_[ivertex] > 0.) {
298  auto tnew = vertices.swt_[ ivertex ] / vertices.nut_[ ivertex ];
299  // prevents from vectorizing if
300  delta += std::pow( vertices.t_[ ivertex ] - tnew, 2 );
301  vertices.t_[ ivertex ] = tnew;
302  }
303 
304 #ifdef VI_DEBUG
305  else {
306  edm::LogInfo("sumw") << "invalid sum of weights in fit: " << endl;//vertices.sw_[ivertex] << endl;
307  if (this->verbose_) {
308  std::cout << " a cluster melted away ? pk=" << vertices.pk_[ ivertex ] << " sumw="
309  << /*vertices.sw_[ivertex] << */ endl;
310  }
311  }
312 #endif
313  }
314 
315  auto osumpi = 1./sumpi;
316  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex )
317  vertices.pk_[ ivertex ] = vertices.pk_[ ivertex ] * vertices.se_[ ivertex ] * osumpi;
318 
319  return delta;
320  };
321 
322  delta += kernel_calc_zt(gvertices);
323 
324  // return how much the prototypes moved
325  return delta;
326 }
327 
328 
329 
330 
331 
332 
334  const unsigned int nv = y.getSize();
335 
336  if (nv < 2) return;
337 
338 
339  bool reordering = true;
340  bool modified = false;
341 
342  while(reordering){
343  reordering = false;
344  for (unsigned int k = 0; (k + 1) < nv; k++) {
345  if ( y.z[k + 1] < y.z[k] ) {
346  auto ztemp = y.z[k]; y.z[k] = y.z[k+1]; y.z[k+1] = ztemp;
347  auto ttemp = y.t[k]; y.t[k] = y.t[k+1]; y.t[k+1] = ttemp;
348  auto ptemp = y.pk[k]; y.pk[k] = y.pk[k+1]; y.pk[k+1] = ptemp;
349  reordering = true;
350  }
351  }
352  modified |= reordering;
353  }
354 
355  if( modified ) y.extractRaw();
356 }
357 
358 
359 
360 bool DAClusterizerInZT_vect::find_nearest( double z, double t, vertex_t & y, unsigned int & k_min, double dz, double dt)const{
361  // find the cluster nearest to (z,t)
362  // distance measure is delta = (delta_z / dz )**2 + (delta_t/ d_t)**2
363  // assumes that clusters are ordered n z
364  // return value is false if o neighbour with distance < 1 is found
365 
366  unsigned int nv = y.getSize();
367 
368  // find nearest in z, binary search later
369  unsigned int k = 0;
370  for(unsigned int k0 = 1; k0 < nv; k0++){
371  if (std::abs( y.z_[k0]-z ) < std::abs( y.z_[k]-z )){
372  k = k0;
373  }
374  }
375 
376  double delta_min = 1.;
377 
378  //search left
379  unsigned int k1 = k;
380  while( ( k1 > 0) && ((y.z[k] - y.z[--k1]) < dz ) ){
381  auto delta = std::pow( (y.z_[k] - y.z_[k1]) / dz, 2)
382  +std::pow( (y.t_[k] - y.t_[k1]) / dt, 2);
383  if (delta < delta_min){
384  k_min = k1;
385  delta_min = delta;
386  }
387  }
388 
389  //search right
390  k1 = k;
391  while( ((++k1) < nv) && ((y.z[k1] - y.z[k]) < dz ) ){
392  auto delta = std::pow( (y.z_[k1] - y.z_[k]) / dz, 2)
393  +std::pow( (y.t_[k1] - y.t_[k]) / dt, 2);
394  if (delta < delta_min){
395  k_min = k1;
396  delta_min = delta;
397  }
398  }
399 
400  return (delta_min < 1.);
401 }
402 
403 
404 
405 
406 bool DAClusterizerInZT_vect::merge(vertex_t & y, double & beta)const{
407  // merge clusters that collapsed or never separated,
408  // return true if vertices were merged, false otherwise
409  const unsigned int nv = y.getSize();
410 
411  if (nv < 2)
412  return false;
413 
414  // merge the smallest distance clusters first
415  unsigned int k1_min = 0, k2_min = 0;
416  double delta_min = 0;
417 
418  for (unsigned int k1 = 0; (k1 + 1) < nv; k1++) {
419  unsigned int k2 = k1;
420  while( (++k2 < nv) && ( std::fabs(y.z[k2] - y.z_[k1] ) < zmerge_ ) ){
421  auto delta = std::pow( (y.z_[k2] - y.z_[k1]) / zmerge_, 2)
422  +std::pow( (y.t_[k2] - y.t_[k1]) / tmerge_, 2);
423  if ( (delta < delta_min) || (k1_min == k2_min) ){
424  k1_min = k1;
425  k2_min = k2;
426  delta_min = delta;
427  }
428  }
429  }
430 
431  if( (k1_min == k2_min) || (delta_min > 1) ){
432  return false;
433  }
434 
435 
436  double rho = y.pk_[k1_min] + y.pk_[k2_min];
437 #ifdef VI_DEBUG
438  if(verbose_){ std::cout << "merging (" << setw(8) << fixed << setprecision(4) << y.z_[k1_min] << ',' << y.t_[k1_min] << ") and (" << y.z_[k2_min] << ',' << y.t_[k2_min] << ")" << " idx=" << k1_min << "," << k2_min << std::endl;}
439 #endif
440  if(rho > 0){
441  y.z_[k1_min] = (y.pk_[k1_min] * y.z_[k1_min] + y.pk_[k2_min] * y.z_[k2_min])/rho;
442  y.t_[k1_min] = (y.pk_[k1_min] * y.t_[k1_min] + y.pk_[k2_min] * y.t_[k2_min])/rho;
443  }else{
444  y.z_[k1_min] = 0.5 * (y.z_[k1_min] + y.z_[k2_min]);
445  y.t_[k1_min] = 0.5 * (y.t_[k1_min] + y.t_[k2_min]);
446  }
447  y.pk_[k1_min] = rho;
448  y.removeItem(k2_min);
449  y.extractRaw();
450  return true;
451 }
452 
453 bool
454 DAClusterizerInZT_vect::purge(vertex_t & y, track_t & tks, double & rho0, const double beta) const {
455  constexpr double eps = 1.e-100;
456  // eliminate clusters with only one significant/unique track
457  const unsigned int nv = y.getSize();
458  const unsigned int nt = tks.getSize();
459 
460  if (nv < 2)
461  return false;
462 
463  double sumpmin = nt;
464  unsigned int k0 = nv;
465 
466  int nUnique = 0;
467  double sump = 0;
468 
469  std::vector<double> inverse_zsums(nt), arg_cache(nt), eik_cache(nt);
470  double * pinverse_zsums;
471  double * parg_cache;
472  double * peik_cache;
473  pinverse_zsums = inverse_zsums.data();
474  parg_cache = arg_cache.data();
475  peik_cache = eik_cache.data();
476  for(unsigned i = 0; i < nt; ++i) {
477  inverse_zsums[i] = tks.Z_sum_[i] > eps ? 1./tks.Z_sum_[i] : 0.0;
478  }
479 
480  for (unsigned int k = 0; k < nv; ++k) {
481 
482  nUnique = 0;
483  sump = 0;
484 
485  const double pmax = y.pk_[k] / (y.pk_[k] + rho0 * local_exp(-beta * dzCutOff_* dzCutOff_));
486  const double pcut = uniquetrkweight_ * pmax;
487  for(unsigned i = 0; i < nt; ++i) {
488  const auto track_z = tks.z_[i];
489  const auto track_t = tks.t_[i];
490  const auto botrack_dz2 = -beta*tks.dz2_[i];
491  const auto botrack_dt2 = -beta*tks.dt2_[i];
492 
493  const auto mult_resz = track_z - y.z_[k];
494  const auto mult_rest = track_t - y.t_[k];
495  parg_cache[i] = botrack_dz2 * ( mult_resz * mult_resz ) + botrack_dt2 * ( mult_rest * mult_rest );
496  }
497  local_exp_list(parg_cache, peik_cache, nt);
498  for (unsigned int i = 0; i < nt; ++i) {
499  const double p = y.pk_[k] * peik_cache[i] * pinverse_zsums[i];
500  sump += p;
501  nUnique += ( ( p > pcut ) & ( tks.pi_[i] > 0 ) );
502  }
503 
504  if ((nUnique < 2) && (sump < sumpmin)) {
505  sumpmin = sump;
506  k0 = k;
507  }
508 
509  }
510 
511  if (k0 != nv) {
512 #ifdef VI_DEBUG
513  if (verbose_) {
514  std::cout << "eliminating prototype at " << std::setw(10) << std::setprecision(4) << y.z_[k0] << "," << y.t_[k0]
515  << " with sump=" << sumpmin
516  << " rho*nt =" << y.pk_[k0]*nt
517  << endl;
518  }
519 #endif
520  y.removeItem(k0);
521  return true;
522  } else {
523  return false;
524  }
525 }
526 
527 
528 
529 
530 double
531 DAClusterizerInZT_vect::beta0(double betamax, track_t const & tks, vertex_t const & y) const {
532 
533  double T0 = 0; // max Tc for beta=0
534  // estimate critical temperature from beta=0 (T=inf)
535  const unsigned int nt = tks.getSize();
536  const unsigned int nv = y.getSize();
537 
538  for (unsigned int k = 0; k < nv; k++) {
539 
540  // vertex fit at T=inf
541  double sumwz = 0;
542  double sumwt = 0;
543  double sumw_z = 0;
544  double sumw_t = 0;
545  for (unsigned int i = 0; i < nt; i++) {
546  double w_z = tks.pi_[i] * tks.dz2_[i];
547  double w_t = tks.pi_[i] * tks.dt2_[i];
548  sumwz += w_z * tks.z_[i];
549  sumwt += w_t * tks.t_[i];
550  sumw_z += w_z;
551  sumw_t += w_t;
552  }
553  y.z_[k] = sumwz / sumw_z;
554  y.t_[k] = sumwt / sumw_t;
555 
556  // estimate Tc, eventually do this in the same loop
557  double szz = 0, stt = 0, szt = 0;
558  double nuz = 0, nut = 0;
559  for (unsigned int i = 0; i < nt; i++) {
560  double dz = (tks.z_[i] - y.z_[k]) * tks.dz2_[i];
561  double dt = (tks.t_[i] - y.t_[k]) * tks.dt2_[i];
562  double w = tks.pi_[i];
563  szz += w * dz * dz;
564  stt += w * dt * dt;
565  szt += w * dz * dt;
566  nuz += w * tks.dz2_[i];
567  nut += w * tks.dt2_[i];
568  }
569  double Tz = szz/nuz;
570  double Tt = stt/nut;
571  double Tc = Tz + Tt + sqrt( pow(Tz - Tt, 2) + 4 * szt * szt / nuz / nut);
572 
573  if (Tc > T0) T0 = Tc;
574  }// vertex loop (normally there should be only one vertex at beta=0)
575 
576 #ifdef VI_DEBUG
577  if(verbose_){
578  std::cout << "DAClustrizerInZT_vect.beta0: Tc = " << T0 << std::endl;
579  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
580  std::cout << "DAClustrizerInZT_vect.beta0: nstep = " << coolingsteps << std::endl;
581  }
582 #endif
583 
584 
585  if (T0 > 1. / betamax) {
586  return betamax / std::pow(coolingFactor_, int(std::log(T0 * betamax) / std::log(coolingFactor_)) - 1);
587  } else {
588  // ensure at least one annealing step
589  return betamax * coolingFactor_;
590  }
591 }
592 
593 double DAClusterizerInZT_vect::get_Tc(const vertex_t & y, int k)const{
594  double Tz = y.szz_[k]/y.nuz_[k]; // actually 0.5*Tc(z)
595  double Tt = y.stt_[k]/y.nut_[k];
596  double mx = y.szt_[k]/y.nuz_[k]*y.szt_[k]/y.nut_[k];
597  return Tz + Tt + sqrt(pow(Tz - Tt, 2) + 4 * mx);
598 }
599 
600 bool
601 DAClusterizerInZT_vect::split(const double beta, track_t &tks, vertex_t & y, double threshold ) const{
602  // split only critical vertices (Tc >~ T=1/beta <==> beta*Tc>~1)
603  // an update must have been made just before doing this (same beta, no merging)
604  // returns true if at least one cluster was split
605 
606  constexpr double epsilonz=1e-3; // minimum split size z
607  constexpr double epsilont=1e-2; // minimum split size t
608  unsigned int nv = y.getSize();
609  const double twoBeta = 2.0*beta;
610 
611  // avoid left-right biases by splitting highest Tc first
612 
613  std::vector<std::pair<double, unsigned int> > critical;
614  for(unsigned int k=0; k<nv; k++){
615  double Tc= get_Tc(y, k);
616  if (beta*Tc > threshold){
617  critical.push_back( make_pair(Tc, k));
618  }
619  }
620  if (critical.empty()) return false;
621 
622 
623  std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int> >() );
624 
625  bool split=false;
626  const unsigned int nt = tks.getSize();
627 
628  for(unsigned int ic=0; ic<critical.size(); ic++){
629  unsigned int k=critical[ic].second;
630 
631  // split direction in the (z,t)-plane
632 
633  double Mzz = y.nuz_[k] - twoBeta*y.szz_[k];
634  double Mtt = y.nut_[k] - twoBeta*y.stt_[k];
635  double Mzt = - twoBeta*y.szt_[k];
636  const double twoMzt = 2.0*Mzt;
637  double D = sqrt( pow(Mtt-Mzz, 2) + twoMzt*twoMzt);
638  double q1 = atan2(-Mtt+Mzz+D, -twoMzt);
639  double l1 = 0.5*(-Mzz-Mtt+D);
640  double l2 = 0.5*(-Mzz-Mtt-D);
641  if( (std::abs(l1) < 1e-4)&&(std::abs(l2) < 1e-4) ){
642  edm::LogWarning("DAClusterizerInZT_vect") << "warning, bad eigenvalues! idx=" << k << " z= " << y.z_[k] << " Mzz=" << Mzz << " Mtt=" << Mtt << " Mzt=" << Mzt << endl;
643  }
644 
645  double qsplit = q1;
646  double cq = cos(qsplit);
647  double sq = sin(qsplit);
648  if(cq < 0){ cq=-cq; sq=-sq; } // prefer cq>0 to keep z-ordering
649 
650  // estimate subcluster positions and weight
651  double p1=0, z1=0, t1=0, wz1=0, wt1=0;
652  double p2=0, z2=0, t2=0, wz2=0, wt2=0;
653  for(unsigned int i=0; i<nt; ++i){
654  if (tks.Z_sum_[i] > 1.e-100) {
655  double lr = (tks.z_[i]-y.z_[k]) * cq + (tks.t[i] - y.t_[k]) * sq;
656  // winner-takes-all, usually overestimates splitting
657  double tl = lr < 0 ? 1.: 0.;
658  double tr = 1. - tl;
659 
660  // soften it, especially at low T
661  double arg = lr * std::sqrt(beta * ( cq*cq*tks.dz2_[i] + sq*sq*tks.dt2_[i] ) );
662  if(std::abs(arg) < 20){
663  double t = local_exp(-arg);
664  tl = t/(t+1.);
665  tr = 1/(t+1.);
666  }
667 
668  double p = y.pk_[k] * tks.pi_[i] * local_exp(-beta * Eik(tks.z_[i], y.z_[k], tks.dz2_[i],
669  tks.t_[i], y.t_[k], tks.dt2_[i])) / tks.Z_sum_[i];
670  double wz = p*tks.dz2_[i];
671  double wt = p*tks.dt2_[i];
672  p1 += p*tl; z1 += wz*tl*tks.z_[i]; t1 += wt*tl*tks.t_[i]; wz1 += wz*tl; wt1 += wt*tl;
673  p2 += p*tr; z2 += wz*tr*tks.z_[i]; t2 += wt*tr*tks.t_[i]; wz2 += wz*tr; wt2 += wt*tr;
674  }
675  }
676 
677  if(wz1 > 0){z1 /= wz1;} else {z1 = y.z_[k] - epsilonz * cq; edm::LogWarning("DAClusterizerInZT_vect") << "warning, wz1 = " << scientific << wz1 << endl;}
678  if(wt1 > 0){t1 /= wt1;} else {t1 = y.t_[k] - epsilont * sq; edm::LogWarning("DAClusterizerInZT_vect") << "warning, w11 = " << scientific << wt1 << endl;}
679  if(wz2 > 0){z2 /= wz2;} else {z2 = y.z_[k] + epsilonz * cq; edm::LogWarning("DAClusterizerInZT_vect") << "warning, wz2 = " << scientific << wz2 << endl;}
680  if(wt2 > 0){t2 /= wt2;} else {t2 = y.t_[k] + epsilont * sq; edm::LogWarning("DAClusterizerInZT_vect") << "warning, wt2 = " << scientific << wt2 << endl;}
681 
682  unsigned int k_min1 = k, k_min2 = k;
683  constexpr double spliteps = 1e-8;
684  while( ((find_nearest(z1, t1, y, k_min1, epsilonz, epsilont) && (k_min1 != k))
685  || (find_nearest(z2,t2, y, k_min2, epsilonz, epsilont) && (k_min2 !=k))) && (std::abs(z2-z1)>spliteps || std::abs(t2-t1)>spliteps) ){
686  z1 = 0.5*( z1 + y.z_[k]); t1 = 0.5*( t1 + y.t_[k]);
687  z2 = 0.5*( z2 + y.z_[k]); t2 = 0.5*( t2 + y.t_[k]);
688  }
689 
690 #ifdef VI_DEBUG
691  if(verbose_){
692  if (std::fabs(y.z_[k] - zdumpcenter_) < zdumpwidth_){
693  std::cout << " T= " << std::setw(10) << std::setprecision(1) << 1./beta
694  << " Tc= " << critical[ic].first
695  << " direction =" << std::setprecision(4) << qsplit
696  << " splitting (" << std::setw(8) << std::fixed << std::setprecision(4) << y.z_[k] <<"," << y.t_[k] << ")"
697  << " --> (" << z1 << ',' << t1<< "),(" << z2 << ',' << t2
698  << ") [" << p1 << "," << p2 << "]" ;
699  if (std::fabs(z2-z1) > epsilonz || std::fabs(t2-t1) > epsilont){
700  std::cout << std::endl;
701  }else{
702  std::cout << " rejected " << std::endl;
703  }
704  }
705  }
706 #endif
707 
708  if(z1>z2){
709  edm::LogInfo("DAClusterizerInZT") << "warning swapping z in split qsplit=" << qsplit << " cq=" << cq << " sq=" << sq << endl;
710  auto ztemp = z1;
711  auto ttemp = t1;
712  auto ptemp = p1;
713  z1 = z2; t1=t2; p1=p2;
714  z2 = ztemp; t2=ttemp; p2=ptemp;
715  }
716 
717 
718 
719 
720  // split if the new subclusters are significantly separated
721  if( std::fabs(z2-z1) > epsilonz || std::fabs(t2-t1) > epsilont){
722  split = true;
723  double pk1 = p1*y.pk_[k]/(p1+p2);
724  double pk2 = p2*y.pk_[k]/(p1+p2);
725 
726  // replace the original by (z2,t2)
727  y.removeItem(k);
728  unsigned int k2 = y.insertOrdered(z2, t2, pk2);
729  // adjust pointers if necessary
730  if( !(k == k2) ){
731  for(unsigned int jc=ic; jc < critical.size(); jc++){
732  if (critical[jc].second > k ) {critical[jc].second--;}
733  if (critical[jc].second >= k2) {critical[jc].second++;}
734  }
735  }
736 
737  // insert (z1,t1) where it belongs
738  unsigned int k1 = y.insertOrdered(z1, t1, pk1);
739  nv++;
740 
741  // adjust remaining pointers
742  for(unsigned int jc=ic; jc < critical.size(); jc++){
743  if (critical[jc].second >= k1) {critical[jc].second++;} // need to backport the ">-"?
744  }
745 
746  }else{
747 #ifdef VI_DEBUG
748  std::cout << "warning ! split rejected, too small." << endl;
749 #endif
750  }
751 
752  }
753  return split;
754 }
755 
756 
757 
759 
760  const unsigned int nv = y.getSize();
761 
762  constexpr double epsilonz = 1e-3; // split all single vertices by 10 um
763  constexpr double epsilont = 1e-2; // split all single vertices by 10 ps
764  constexpr double zsep = 2 * epsilonz; // split vertices that are isolated by at least zsep (vertices that haven't collapsed)
765  constexpr double tsep = 2 * epsilont;
766  vertex_t y1;
767 
768 #ifdef VI_DEBUG
769  if (verbose_) {
770  std::cout << "Before Split "<< std::endl;
771  y.debugOut();
772  }
773 #endif
774 
775  for (unsigned int k = 0; k < nv; k++) {
776  if (
777  ( ( (k == 0) || ( y.z_[k - 1] < (y.z_[k] - zsep)) ) &&
778  ( ((k + 1) == nv) || ( y.z_[k + 1] > (y.z_[k] + zsep)) ) ) )
779  {
780  // isolated prototype, split
781  double new_z = y.z[k] - epsilonz;
782  double new_t = y.t[k] - epsilont;
783  y.z[k] = y.z[k] + epsilonz;
784  y.t[k] = y.t[k] + epsilont;
785 
786  double new_pk = 0.5 * y.pk[k];
787  y.pk[k] = 0.5 * y.pk[k];
788 
789  y1.addItem(new_z, new_t, new_pk);
790  y1.addItem(y.z_[k], y.t_[k], y.pk_[k]);
791  }
792  else if ( (y1.getSize() == 0 ) ||
793  (y1.z_[y1.getSize() - 1] < (y.z_[k] - zsep) ) ||
794  (y1.t_[y1.getSize() - 1] < (y.t_[k] - tsep) ))
795  {
796  y1.addItem(y.z_[k], y.t_[k], y.pk_[k]);
797  }
798  else
799  {
800  y1.z_[y1.getSize() - 1] = y1.z_[y1.getSize() - 1] - epsilonz;
801  y1.t_[y1.getSize() - 1] = y1.t_[y1.getSize() - 1] - epsilont;
802  y.z_[k] = y.z_[k] + epsilonz;
803  y.t_[k] = y.t_[k] + epsilont;
804  y1.addItem( y.z_[k], y.t_[k] , y.pk_[k]);
805  }
806  }// vertex loop
807 
808  y = y1;
809  y.extractRaw();
810 
811 #ifdef VI_DEBUG
812  if (verbose_) {
813  std::cout << "After split " << std::endl;
814  y.debugOut();
815  }
816 #endif
817 }
818 
819 
820 
821 vector<TransientVertex>
822 DAClusterizerInZT_vect::vertices(const vector<reco::TransientTrack> & tracks, const int verbosity) const {
823  track_t && tks = fill(tracks);
824  tks.extractRaw();
825 
826  unsigned int nt = tks.getSize();
827  double rho0 = 0.0; // start with no outlier rejection
828 
829  vector<TransientVertex> clusters;
830  if (tks.getSize() == 0) return clusters;
831 
832  vertex_t y; // the vertex prototypes
833 
834  // initialize:single vertex at infinite temperature
835  y.addItem( 0, 0, 1.0);
836 
837  int niter = 0; // number of iterations
838 
839 
840  // estimate first critical temperature
841  double beta = beta0(betamax_, tks, y);
842 #ifdef VI_DEBUG
843  if ( verbose_) std::cout << "Beta0 is " << beta << std::endl;
844 #endif
845 
846  niter = 0;
847  while ((update(beta, tks, y, false, rho0) > 1.e-6) &&
848  (niter++ < maxIterations_)) {}
849 
850  // annealing loop, stop when T<minT (i.e. beta>1/minT)
851 
852  double betafreeze = betamax_ * sqrt(coolingFactor_);
853 
854  while (beta < betafreeze) {
855  if(useTc_){
856  update(beta, tks,y, false, rho0);
857  zorder(y);
858  while(merge(y, beta)){update(beta, tks, y, false, rho0);}
859  split(beta, tks, y);
860  beta=beta/coolingFactor_;
861  }else{
862  beta=beta/coolingFactor_;
863  splitAll(y);
864  }
865 
866  // make sure we are not too far from equilibrium before cooling further
867  niter = 0;
868  while ((update(beta, tks, y, false, rho0) > 1.e-6) &&
869  (niter++ < maxIterations_)) {}
870 
871 #ifdef VI_DEBUG
872  if(verbose_){ dump( beta, y, tks, 0); }
873 #endif
874  }
875 
876 
877  if(useTc_){
878  //last round of splitting, make sure no critical clusters are left
879 #ifdef VI_DEBUG
880  if(verbose_){ std::cout << "last spliting at " << 1./beta << std::endl; }
881 #endif
882  update(beta, tks,y, false, rho0);// make sure Tc is up-to-date
883  zorder(y);
884  while(merge(y,beta)){update(beta, tks,y, false, rho0);}
885  unsigned int ntry=0;
886  double threshold = 1.0;
887  while( split(beta, tks, y, threshold) && (ntry++<10) ){
888  niter=0;
889  while((update(beta, tks,y, false, rho0)>1.e-6) && (niter++ < maxIterations_)){}
890  zorder(y);
891 #ifdef VI_DEBUG
892  if(verbose_) dump(beta, y, tks, 2);
893 #endif
894 
895  while(merge(y,beta)){update(beta, tks,y, false, rho0);}
896 #ifdef VI_DEBUG
897  if(verbose_){
898  std::cout << "after final splitting, try " << ntry << std::endl;
899  dump(beta, y, tks, 2);
900  }
901 #endif
902  // relax splitting a bit to reduce multiple split-merge cycles of the same cluster
903  threshold *= 1.1;
904  }
905  }else{
906  // merge collapsed clusters
907  while(merge(y,beta)){update(beta, tks,y, false, rho0);}
908  }
909 
910 #ifdef VI_DEBUG
911  if (verbose_) {
912  update(beta, tks,y, false, rho0);
913  std::cout << "dump after 1st merging " << endl;
914  dump(beta, y, tks, 2);
915  }
916 #endif
917 
918 
919  // switch on outlier rejection at T=minT
920  if(dzCutOff_ > 0){
921  rho0 = 1./nt;
922  for(unsigned int a=0; a<10; a++){ update(beta, tks, y, true, a*rho0/10);} // adiabatic turn-on
923  }
924 
925  niter=0;
926  while ((update(beta, tks, y, true, rho0) > 1.e-8) && (niter++ < maxIterations_)) {};
927 #ifdef VI_DEBUG
928  if (verbose_) {
929  std::cout << "dump after noise-suppression, rho0=" << rho0 << " niter = " << niter << endl;
930  dump(beta, y, tks, 2);
931  }
932 #endif
933 
934  // merge again (some cluster split by outliers collapse here)
935  zorder(y);
936  while (merge(y, beta)) {update(beta, tks, y, true, rho0); }
937 #ifdef VI_DEBUG
938  if (verbose_) {
939  std::cout << "dump after merging " << endl;
940  dump(beta, y, tks, 2);
941  }
942 #endif
943 
944  // go down to the purging temperature (if it is lower than tmin)
945  while( beta < betapurge_ ){
946  beta = min( beta/coolingFactor_, betapurge_);
947  niter = 0;
948  while ((update(beta, tks, y, false, rho0) > 1.e-8) && (niter++ < maxIterations_)) {}
949  }
950 
951 
952  // eliminate insigificant vertices, this is more restrictive at higher T
953  while (purge(y, tks, rho0, beta)) {
954  niter = 0;
955  while (( update(beta, tks, y, true, rho0) > 1.e-6 ) && (niter++ < maxIterations_)) {
956  zorder(y);
957  }
958  }
959 
960 #ifdef VI_DEBUG
961  if (verbose_) {
962  update(beta, tks,y, true, rho0);
963  std::cout << " after purging " << std:: endl;
964  dump(beta, y, tks, 2);
965  }
966 #endif
967 
968  // optionally cool some more without doing anything, to make the assignment harder
969  while( beta < betastop_ ){
970  beta = min( beta/coolingFactor_, betastop_);
971  niter =0;
972  while ((update(beta, tks, y, true, rho0) > 1.e-8) && (niter++ < maxIterations_)) {}
973  zorder(y);
974  }
975 
976 #ifdef VI_DEBUG
977  if (verbose_) {
978  std::cout << "Final result, rho0=" << std::scientific << rho0 << endl;
979  dump(beta, y, tks, 2);
980  }
981 #endif
982 
983 
984  // new, merge here and not in "clusterize"
985  // final merging step
986  double betadummy = 1;
987  while( merge(y, betadummy) );
988 
989  // select significant tracks and use a TransientVertex as a container
990  GlobalError dummyError(0.01, 0, 0.01, 0., 0., 0.01);
991 
992  // ensure correct normalization of probabilities, should makes double assignment reasonably impossible
993  const unsigned int nv = y.getSize();
994  for (unsigned int k = 0; k < nv; k++)
995  if ( edm::isNotFinite(y.pk_[k]) || edm::isNotFinite(y.z_[k]) ) { y.pk_[k]=0; y.z_[k]=0;}
996 
997  for (unsigned int i = 0; i < nt; i++) // initialize
998  tks.Z_sum_[i] = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
999 
1000  // improve vectorization (does not require reduction ....)
1001  for (unsigned int k = 0; k < nv; k++) {
1002  for (unsigned int i = 0; i < nt; i++)
1003  tks.Z_sum_[i] += y.pk_[k] * local_exp(-beta * Eik(tks.z_[i], y.z_[k],tks.dz2_[i], tks.t_[i], y.t_[k],tks.dt2_[i]));
1004  }
1005 
1006 
1007  for (unsigned int k = 0; k < nv; k++) {
1008  GlobalPoint pos(0, 0, y.z_[k]);
1009 
1010  vector<reco::TransientTrack> vertexTracks;
1011  for (unsigned int i = 0; i < nt; i++) {
1012  if (tks.Z_sum_[i] > 1e-100) {
1013 
1014  double p = y.pk_[k] * local_exp(-beta * Eik(tks.z_[i], y.z_[k], tks.dz2_[i],
1015  tks.t_[i], y.t_[k], tks.dt2_[i] )) / tks.Z_sum_[i];
1016  if ((tks.pi_[i] > 0) && (p > mintrkweight_)) {
1017  vertexTracks.push_back(*(tks.tt[i]));
1018  tks.Z_sum_[i] = 0; // setting Z=0 excludes double assignment
1019  }
1020  }
1021  }
1022  TransientVertex v(pos, y.t_[k], dummyError, vertexTracks, 0);
1023  clusters.push_back(v);
1024  }
1025 
1026  return clusters;
1027 
1028 }
1029 
1030 vector<vector<reco::TransientTrack> > DAClusterizerInZT_vect::clusterize(
1031  const vector<reco::TransientTrack> & tracks) const {
1032 
1033 #ifdef VI_DEBUG
1034  if (verbose_) {
1035  std::cout << "###################################################" << endl;
1036  std::cout << "# vectorized DAClusterizerInZT_vect::clusterize nt=" << tracks.size() << endl;
1037  std::cout << "###################################################" << endl;
1038  }
1039 #endif
1040 
1041  vector<vector<reco::TransientTrack> > clusters;
1042  vector<TransientVertex> && pv = vertices(tracks);
1043 
1044 #ifdef VI_DEBUG
1045  if (verbose_) {
1046  std::cout << "# DAClusterizerInZT::clusterize pv.size=" << pv.size() << endl;
1047  }
1048 #endif
1049 
1050  if (pv.empty()) {
1051  return clusters;
1052  }
1053 
1054  // fill into clusters, don't merge
1055  for (auto k = pv.begin(); k != pv.end(); k++) {
1056  vector<reco::TransientTrack> aCluster = k->originalTracks();
1057  if (aCluster.size()>1){
1058  clusters.push_back(aCluster);
1059  }
1060  }
1061 
1062  return clusters;
1063 
1064  /*
1065  //leaving the merging version here in case we want to turn it back on the future
1066  // needs optimization and validation
1067  // fill into clusters and merge
1068  vector<reco::TransientTrack> aCluster = pv.begin()->originalTracks();
1069 
1070  for (auto k = pv.begin() + 1; k != pv.end(); k++) {
1071  if ( (std::abs(k->position().z() - (k - 1)->position().z()) > (2 * vertexSize_))
1072  ||(std::abs(k->time() - (k - 1)->time() ) > (2 * vertexSizeTime_)) ){ // still not perfect
1073  // close a cluster
1074  if (aCluster.size()>1){
1075  clusters.push_back(aCluster);
1076  }else{
1077  #ifdef VI_DEBUG
1078  if(verbose_){
1079  std::cout << " one track cluster at " << k->position().z() << " suppressed" << std::endl;
1080  }
1081  #endif
1082  }
1083  aCluster.clear();
1084  }
1085  for (unsigned int i = 0; i < k->originalTracks().size(); i++) {
1086  aCluster.push_back(k->originalTracks()[i]);
1087  }
1088 
1089  }
1090  clusters.emplace_back(std::move(aCluster));
1091 
1092  return clusters;
1093  */
1094 }
1095 
1096 
1097 
1098 void DAClusterizerInZT_vect::dump(const double beta, const vertex_t & y,
1099  const track_t & tks, int verbosity) const {
1100 #ifdef VI_DEBUG
1101  const unsigned int nv = y.getSize();
1102  const unsigned int nt = tks.getSize();
1103 
1104  std::vector< unsigned int > iz;
1105  for(unsigned int j=0; j<nt; j++){ iz.push_back(j); }
1106  std::sort(iz.begin(), iz.end(), [tks](unsigned int a, unsigned int b){ return tks.z_[a]<tks.z_[b];} );
1107  std::cout << std::endl;
1108  std::cout << "-----DAClusterizerInZT::dump ----" << nv << " clusters " << std::endl;
1109  string h = " ";
1110  std::cout << h << " k= ";
1111  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
1112  if (std::fabs(y.z_[ivertex]-zdumpcenter_) < zdumpwidth_){
1113  std::cout << setw(8) << fixed << ivertex;
1114  }
1115  }
1116  std::cout << endl;
1117 
1118  std::cout << h << " z= ";
1119  std::cout << setprecision(4);
1120  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
1121  if (std::fabs(y.z_[ivertex]-zdumpcenter_) < zdumpwidth_){
1122  std::cout << setw(8) << fixed << y.z_[ivertex];
1123  }
1124  }
1125  std::cout << endl;
1126 
1127  std::cout << h << " t= ";
1128  std::cout << setprecision(4);
1129  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
1130  if (std::fabs(y.z_[ivertex]-zdumpcenter_) < zdumpwidth_){
1131  std::cout << setw(8) << fixed << y.t_[ivertex];
1132  }
1133  }
1134  std::cout << endl;
1135 
1136  std::cout << "T=" << setw(15) << 1. / beta
1137  << " Tmin =" << setw(10) << 1./betamax_
1138  << " Tc= ";
1139  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
1140  if (std::fabs(y.z_[ivertex]-zdumpcenter_) < zdumpwidth_){
1141  double Tc = get_Tc(y, ivertex);
1142  std::cout << setw(8) << fixed << setprecision(1) << Tc;
1143  }
1144  }
1145  std::cout << endl;
1146 
1147  std::cout << h << "pk= ";
1148  double sumpk = 0;
1149  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
1150  sumpk += y.pk_[ivertex];
1151  if (std::fabs(y.z_[ivertex] - zdumpcenter_) > zdumpwidth_) continue;
1152  std::cout << setw(8) << setprecision(4) << fixed << y.pk_[ivertex];
1153  }
1154  std::cout << endl;
1155 
1156  std::cout << h << "nt= ";
1157  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
1158  sumpk += y.pk_[ivertex];
1159  if (std::fabs(y.z_[ivertex] - zdumpcenter_) > zdumpwidth_) continue;
1160  std::cout << setw(8) << setprecision(1) << fixed << y.pk_[ivertex]*nt;
1161  }
1162  std::cout << endl;
1163 
1164  if (verbosity > 0) {
1165  double E = 0, F = 0;
1166  std::cout << endl;
1167  std::cout
1168  << "---- z +/- dz t +/- dt ip +/-dip pt phi eta weights ----"
1169  << endl;
1170  std::cout << setprecision(4);
1171  for (unsigned int i0 = 0; i0 < nt; i0++) {
1172  unsigned int i = iz[i0];
1173  if (tks.Z_sum_[i] > 0) {
1174  F -= std::log(tks.Z_sum_[i]) / beta;
1175  }
1176  double tz = tks.z_[i];
1177 
1178  if( std::fabs(tz - zdumpcenter_) > zdumpwidth_) continue;
1179  std::cout << setw(4) << i << ")" << setw(8) << fixed << setprecision(4)
1180  << tz << " +/-" << setw(6) << sqrt(1./tks.dz2_[i]);
1181 
1182  if( tks.dt2_[i] >0){
1183  std::cout << setw(8) << fixed << setprecision(4)
1184  << tks.t_[i] << " +/-" << setw(6) << sqrt(1./tks.dt2_[i]);
1185  }else{
1186  std::cout << " ";
1187  }
1188 
1189  if (tks.tt[i]->track().quality(reco::TrackBase::highPurity)) {
1190  std::cout << " *";
1191  } else {
1192  std::cout << " ";
1193  }
1194  if (tks.tt[i]->track().hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1)) {
1195  std::cout << "+";
1196  } else {
1197  std::cout << "-";
1198  }
1199  std::cout << setw(1)
1200  << tks.tt[i]->track().hitPattern().pixelBarrelLayersWithMeasurement(); // see DataFormats/TrackReco/interface/HitPattern.h
1201  std::cout << setw(1)
1202  << tks.tt[i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
1203  std::cout << setw(1) << hex
1204  << tks.tt[i]->track().hitPattern().trackerLayersWithMeasurement()
1205  - tks.tt[i]->track().hitPattern().pixelLayersWithMeasurement()
1206  << dec;
1207  std::cout << "=" << setw(1) << hex
1208  << tks.tt[i]->track().hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS)
1209  << dec;
1210 
1211  Measurement1D IP =
1212  tks.tt[i]->stateAtBeamLine().transverseImpactParameter();
1213  std::cout << setw(8) << IP.value() << "+/-" << setw(6) << IP.error();
1214  std::cout << " " << setw(6) << setprecision(2)
1215  << tks.tt[i]->track().pt() * tks.tt[i]->track().charge();
1216  std::cout << " " << setw(5) << setprecision(2)
1217  << tks.tt[i]->track().phi() << " " << setw(5)
1218  << setprecision(2) << tks.tt[i]->track().eta();
1219 
1220  double sump = 0.;
1221  for (unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
1222  if (std::fabs(y.z_[ivertex]-zdumpcenter_) > zdumpwidth_) continue;
1223 
1224  if ((tks.pi_[i] > 0) && (tks.Z_sum_[i] > 0)) {
1225  //double p=pik(beta,tks[i],*k);
1226  double p = y.pk_[ivertex] * exp(-beta * Eik(tks.z_[i], y.z_[ivertex], tks.dz2_[i],
1227  tks.t_[i], y.t_[ivertex], tks.dt2_[i] )) / tks.Z_sum_[i];
1228  if (p > 0.0001) {
1229  std::cout << setw(8) << setprecision(3) << p;
1230  } else {
1231  std::cout << " . ";
1232  }
1233  E += p * Eik(tks.z_[i], y.z_[ivertex], tks.dz2_[i],
1234  tks.t_[i], y.t_[ivertex], tks.dt2_[i] );
1235  sump += p;
1236  } else {
1237  std::cout << " ";
1238  }
1239  }
1240  std::cout << endl;
1241  }
1242  std::cout << endl << "T=" << 1 / beta << " E=" << E << " n=" << y.getSize()
1243  << " F= " << F << endl << "----------" << endl;
1244  }
1245 #endif
1246 }
dbl * delta
Definition: mlp_gen.cc:36
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
float dt
Definition: AMPTWrapper.h:126
unsigned int insertOrdered(double z, double t, double pk)
const double w
Definition: UKUtility.cc:23
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const override
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double error() const
Definition: Measurement1D.h:27
void splitAll(vertex_t &y) const
DAClusterizerInZT_vect(const edm::ParameterSet &conf)
bool merge(vertex_t &y, double &beta) const
void addItem(double new_z, double new_t, double new_dz2, double new_dt2, const reco::TransientTrack *new_tt, double new_pi)
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks, const int verbosity=0) const
A arg
Definition: Factorize.h:38
void zorder(vertex_t &y) const
U second(std::pair< T, U > const &p)
double get_Tc(const vertex_t &y, int k) const
bool purge(vertex_t &, track_t &, double &, const double) const
T sqrt(T t)
Definition: SSEVec.h:18
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 find_nearest(double z, double t, vertex_t &y, unsigned int &k_min, double dz, double dt) const
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
T min(T a, T b)
Definition: MathUtil.h:58
double p2[4]
Definition: TauolaWrapper.h:90
#define LogTrace(id)
int nt
Definition: AMPTWrapper.h:32
double update(double beta, track_t &gtracks, vertex_t &gvertices, bool useRho0, const double &rho0) const
int k[5][pyjets_maxn]
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:152
double q1[4]
Definition: TauolaWrapper.h:87
void addItem(double new_z, double new_t, double new_pk)
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:88
double b
Definition: hdecay.h:120
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
double value() const
Definition: Measurement1D.h:25
double p1[4]
Definition: TauolaWrapper.h:89
#define update(a, b)
double a
Definition: hdecay.h:121
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
std::vector< const reco::TransientTrack * > tt
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:281
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0) const
double split
Definition: MVATrainer.cc:139
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
#define constexpr
def merge(dictlist, TELL=False)
Definition: MatrixUtil.py:194