11 #include "vdt/vdtMath.h" 22 maxIterations_ = 1000;
36 coolingFactor_ = conf.
getParameter<
double>(
"coolingFactor");
39 uniquetrkweight_ = conf.
getParameter<
double>(
"uniquetrkweight");
40 uniquetrkminp_ = conf.
getParameter<
double>(
"uniquetrkminp");
43 convergence_mode_ = conf.
getParameter<
int>(
"convergence_mode");
47 block_size_ = conf.
getParameter<
unsigned int>(
"block_size");
48 overlap_frac_ = conf.
getParameter<
double>(
"overlap_frac");
51 std::cout <<
"DAClusterizerinZ_vect: mintrkweight = " << mintrkweight << std::endl;
52 std::cout <<
"DAClusterizerinZ_vect: uniquetrkweight = " << uniquetrkweight_ << std::endl;
53 std::cout <<
"DAClusterizerInZ_vect: uniquetrkminp = " << uniquetrkminp_ << std::endl;
54 std::cout <<
"DAClusterizerinZ_vect: zmerge = " << zmerge_ << std::endl;
55 std::cout <<
"DAClusterizerinZ_vect: Tmin = " <<
Tmin << std::endl;
57 std::cout <<
"DAClusterizerinZ_vect: Tstop = " <<
Tstop << std::endl;
58 std::cout <<
"DAClusterizerinZ_vect: vertexSize = " << vertexSize_ << std::endl;
59 std::cout <<
"DAClusterizerinZ_vect: coolingFactor = " << coolingFactor_ << std::endl;
60 std::cout <<
"DAClusterizerinZ_vect: d0CutOff = " << d0CutOff_ << std::endl;
61 std::cout <<
"DAClusterizerinZ_vect: dzCutOff = " << dzCutOff_ << std::endl;
62 std::cout <<
"DAClusterizerInZ_vect: zrange = " << sel_zrange_ << std::endl;
63 std::cout <<
"DAClusterizerinZ_vect: convergence mode = " << convergence_mode_ << std::endl;
64 std::cout <<
"DAClusterizerinZ_vect: delta_highT = " << delta_highT_ << std::endl;
65 std::cout <<
"DAClusterizerinZ_vect: delta_lowT = " << delta_lowT_ << std::endl;
66 std::cout <<
"DAClusterizerinZ_vect: DEBUGLEVEL " << DEBUGLEVEL << std::endl;
69 if (convergence_mode_ > 1) {
71 <<
"DAClusterizerInZ_vect: invalid convergence_mode " << convergence_mode_ <<
" reset to default " << 0;
72 convergence_mode_ = 0;
78 <<
"DAClusterizerInZ_vect: invalid Tmin " <<
Tmin <<
" reset to default " << 1. / betamax_;
85 <<
"DAClusterizerInZ_vect: invalid Tpurge " <<
Tpurge <<
" set to " <<
Tmin;
92 <<
"DAClusterizerInZ_vect: invalid Tstop " <<
Tstop <<
" set to " <<
max(1.,
Tpurge);
95 betastop_ = 1. /
Tstop;
99 inline double local_exp(
double const& inp) {
return vdt::fast_exp(inp); }
101 inline void local_exp_list_range(
double const* __restrict__ arg_inp,
102 double* __restrict__ arg_out,
105 for (
auto i = kmin;
i != kmax; ++
i)
106 arg_out[
i] = vdt::fast_exp(arg_inp[
i]);
112 if (!(nv == 999999)) {
118 if (!(
nt == 999999)) {
124 assert(
v.zvtx_vec.size() == nv);
125 assert(
v.rho_vec.size() == nv);
126 assert(
v.swz_vec.size() == nv);
127 assert(
v.exp_arg_vec.size() == nv);
128 assert(
v.exp_vec.size() == nv);
130 assert(
v.swz_vec.size() == nv);
131 assert(
v.swE_vec.size() == nv);
133 assert(
v.zvtx == &
v.zvtx_vec.front());
134 assert(
v.rho == &
v.rho_vec.front());
135 assert(
v.exp_arg == &
v.exp_arg_vec.front());
137 assert(
v.swz == &
v.swz_vec.front());
139 assert(
v.swE == &
v.swE_vec.front());
141 for (
unsigned int k = 0;
k < nv - 1;
k++) {
142 if (
v.zvtx_vec[
k] <=
v.zvtx_vec[
k + 1])
144 cout <<
" Z, cluster z-ordering assertion failure z[" <<
k <<
"] =" <<
v.zvtx_vec[
k] <<
" z[" <<
k + 1
145 <<
"] =" <<
v.zvtx_vec[
k + 1] << endl;
161 for (
unsigned int i = 0;
i <
nt;
i++) {
164 cout <<
"track vertex range assertion failure" <<
i <<
"/" <<
nt <<
" kmin,kmax=" << tks.
kmin[
i] <<
", " 165 << tks.
kmax[
i] <<
" nv=" << nv << endl;
168 for (
unsigned int i = 0;
i <
nt;
i++) {
177 for (
auto it =
tracks.begin(); it !=
tracks.end(); it++) {
178 if (!(*it).isValid())
181 double t_z = ((*it).stateAtBeamLine().trackStateAtPCA()).
position().z();
182 if (std::fabs(
t_z) > 1000.)
184 auto const& t_mom = (*it).stateAtBeamLine().trackStateAtPCA().momentum();
187 double t_dz2 =
std::pow((*it).track().dzError(), 2)
199 edm::LogWarning(
"DAClusterizerinZ_vect") <<
"rejected track t_tkwt " << t_tkwt;
208 tks.
osumtkwt = sumtkwt > 0 ? 1. / sumtkwt : 0.;
211 if (DEBUGLEVEL > 0) {
220 inline double Eik(
double t_z,
double k_z,
double t_dz2) {
return std::pow(
t_z - k_z, 2) * t_dz2; }
224 const unsigned int nv = gvertices.
getSize();
225 const unsigned int nt = gtracks.
getSize();
228 edm::LogWarning(
"DAClusterizerinZ_vect") <<
"empty cluster list in set_vtx_range";
232 for (
auto itrack = 0
U; itrack <
nt; ++itrack) {
236 unsigned int kmin =
min(nv - 1, gtracks.
kmin[itrack]);
239 while ((kmin > 0) && (gvertices.
zvtx[kmin - 1] >
zmin)) {
243 while ((kmin < (nv - 1)) && (gvertices.
zvtx[kmin] <
zmin)) {
249 unsigned int kmax =
min(nv - 1, gtracks.
kmax[itrack] - 1);
253 while ((kmax < (nv - 1)) && (gvertices.
zvtx[kmax + 1] <
zmax)) {
257 while ((kmax > 0) && (gvertices.
zvtx[kmax] >
zmax)) {
263 gtracks.
kmin[itrack] = kmin;
264 gtracks.
kmax[itrack] = kmax + 1;
267 gtracks.
kmax[itrack] =
min(nv,
max(kmin, kmax) + 1);
273 const unsigned int nt = gtracks.
getSize();
274 const unsigned int nv = gvertices.
getSize();
275 for (
auto itrack = 0
U; itrack <
nt; ++itrack) {
276 gtracks.
kmin[itrack] = 0;
277 gtracks.
kmax[itrack] = nv;
282 double beta,
track_t& gtracks,
vertex_t& gvertices,
const double rho0,
const bool updateTc)
const {
287 const unsigned int nt = gtracks.
getSize();
288 const unsigned int nv = gvertices.
getSize();
293 Z_init = rho0 * local_exp(-
beta * dzCutOff_ * dzCutOff_);
297 auto kernel_calc_exp_arg_range = [
beta](
const unsigned int itrack,
300 const unsigned int kmin,
301 const unsigned int kmax) {
302 const double track_z =
tracks.zpca[itrack];
303 const double botrack_dz2 = -
beta *
tracks.dz2[itrack];
306 for (
unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
307 auto mult_res = track_z -
vertices.zvtx[ivertex];
308 vertices.exp_arg[ivertex] = botrack_dz2 * (mult_res * mult_res);
312 auto kernel_add_Z_range = [Z_init](
313 vertex_t const&
vertices,
const unsigned int kmin,
const unsigned int kmax) ->
double {
314 double ZTemp = Z_init;
315 for (
unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
321 auto kernel_calc_normalization_range = [updateTc](
const unsigned int track_num,
324 const unsigned int kmin,
325 const unsigned int kmax) {
326 auto o_trk_sum_Z =
tracks.tkwt[track_num] /
tracks.sum_Z[track_num];
327 auto o_trk_dz2 =
tracks.dz2[track_num];
328 auto tmp_trk_z =
tracks.zpca[track_num];
333 for (
unsigned int k = kmin;
k < kmax; ++
k) {
343 for (
unsigned int k = kmin;
k < kmax; ++
k) {
353 for (
auto ivertex = 0
U; ivertex < nv; ++ivertex) {
354 gvertices.
se[ivertex] = 0.0;
355 gvertices.
sw[ivertex] = 0.0;
356 gvertices.
swz[ivertex] = 0.0;
357 gvertices.
swE[ivertex] = 0.0;
360 for (
auto ivertex = 0
U; ivertex < nv; ++ivertex) {
361 gvertices.
se[ivertex] = 0.0;
362 gvertices.
sw[ivertex] = 0.0;
363 gvertices.
swz[ivertex] = 0.0;
368 for (
auto itrack = 0
U; itrack <
nt; ++itrack) {
369 const unsigned int kmin = gtracks.
kmin[itrack];
370 const unsigned int kmax = gtracks.
kmax[itrack];
372 kernel_calc_exp_arg_range(itrack, gtracks, gvertices, kmin, kmax);
373 local_exp_list_range(gvertices.
exp_arg, gvertices.
exp, kmin, kmax);
374 gtracks.
sum_Z[itrack] = kernel_add_Z_range(gvertices, kmin, kmax);
377 gtracks.
sum_Z[itrack] = 0.0;
379 if (gtracks.
sum_Z[itrack] > 1.e-100) {
380 kernel_calc_normalization_range(itrack, gtracks, gvertices, kmin, kmax);
386 auto obeta = -1. /
beta;
387 for (
auto ivertex = 0
U; ivertex < nv; ++ivertex) {
388 gvertices.
swE[ivertex] *= obeta;
396 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
405 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex)
411 double delta = kernel_calc_z(gvertices);
419 unsigned int niter = 0;
421 double delta_max = delta_lowT_;
423 if (convergence_mode_ == 0) {
424 delta_max = delta_max0;
425 }
else if (convergence_mode_ == 1) {
429 set_vtx_range(
beta, tks,
v);
430 double delta_sum_range = 0;
431 std::vector<double> z0 =
v.zvtx_vec;
433 while (niter++ < maxIterations_) {
435 delta_sum_range +=
delta;
437 if (delta_sum_range > zrange_min_) {
438 for (
unsigned int k = 0;
k <
v.getSize();
k++) {
439 if (
std::abs(
v.zvtx_vec[
k] - z0[
k]) > zrange_min_) {
440 set_vtx_range(
beta, tks,
v);
448 if (
delta < delta_max) {
454 if (DEBUGLEVEL > 0) {
455 std::cout <<
"DAClusterizerInZ_vect.thermalize niter = " << niter <<
" at T = " << 1 /
beta 456 <<
" nv = " <<
v.getSize() << std::endl;
469 const unsigned int nv = y.getSize();
475 std::vector<std::pair<double, unsigned int>>
critical;
476 for (
unsigned int k = 0; (
k + 1) < nv;
k++) {
477 if (std::fabs(y.zvtx[
k + 1] - y.zvtx[
k]) < zmerge_) {
478 critical.push_back(make_pair(std::fabs(y.zvtx[
k + 1] - y.zvtx[
k]),
k));
484 std::stable_sort(
critical.begin(),
critical.end(), std::less<std::pair<double, unsigned int>>());
486 for (
unsigned int ik = 0; ik <
critical.size(); ik++) {
488 double rho = y.rho[
k] + y.rho[
k + 1];
492 if (DEBUGLEVEL > 1) {
493 std::cout <<
"merging " <<
fixed << setprecision(4) << y.zvtx[
k + 1] <<
" and " << y.zvtx[
k]
494 <<
" sw = " << y.sw[
k] + y.sw[
k + 1] << std::endl;
499 y.zvtx[
k] = (y.rho[
k] * y.zvtx[
k] + y.rho[
k + 1] * y.zvtx[
k + 1]) / rho;
501 y.zvtx[
k] = 0.5 * (y.zvtx[
k] + y.zvtx[
k + 1]);
504 y.sw[
k] += y.sw[
k + 1];
505 y.removeItem(
k + 1, tks);
506 set_vtx_range(
beta, tks, y);
515 constexpr
double eps = 1.e-100;
517 const unsigned int nv = y.getSize();
523 std::vector<double> sump_v(nv), arg_cache_v(nv), exp_cache_v(nv), pcut_cache_v(nv);
524 std::vector<int> nUnique_v(nv);
525 double* __restrict__ parg_cache;
526 double* __restrict__ pexp_cache;
527 double* __restrict__ ppcut_cache;
528 double* __restrict__ psump;
529 int* __restrict__ pnUnique;
530 int constexpr nunique_min_ = 2;
532 set_vtx_range(
beta, tks, y);
534 parg_cache = arg_cache_v.data();
535 pexp_cache = exp_cache_v.data();
536 ppcut_cache = pcut_cache_v.data();
537 psump = sump_v.data();
538 pnUnique = nUnique_v.data();
540 const auto rhoconst = rho0 * local_exp(-
beta * dzCutOff_ * dzCutOff_);
541 for (
unsigned int k = 0;
k < nv;
k++) {
542 const double pmax = y.rho[
k] / (y.rho[
k] + rhoconst);
543 ppcut_cache[
k] = uniquetrkweight_ * pmax;
546 for (
unsigned int i = 0;
i <
nt;
i++) {
547 const auto invZ = ((tks.
sum_Z[
i] >
eps) && (tks.
tkwt[
i] > uniquetrkminp_)) ? 1. / tks.
sum_Z[
i] : 0.;
548 const auto track_z = tks.
zpca[
i];
549 const auto botrack_dz2 = -
beta * tks.
dz2[
i];
550 const auto kmin = tks.
kmin[
i];
551 const auto kmax = tks.
kmax[
i];
553 for (
unsigned int k = kmin;
k < kmax;
k++) {
554 const auto mult_resz = track_z - y.zvtx[
k];
555 parg_cache[
k] = botrack_dz2 * (mult_resz * mult_resz);
558 local_exp_list_range(parg_cache, pexp_cache, kmin, kmax);
560 for (
unsigned int k = kmin;
k < kmax;
k++) {
561 const double p = y.rho[
k] * pexp_cache[
k] * invZ;
563 pnUnique[
k] += (
p > ppcut_cache[
k]) ? 1 : 0;
568 unsigned int k0 = nv;
569 for (
unsigned k = 0;
k < nv;
k++) {
570 if ((pnUnique[
k] < nunique_min_) && (psump[
k] < sumpmin)) {
579 if (DEBUGLEVEL > 1) {
580 std::cout <<
"eliminating prototype at " << std::setw(10) << std::setprecision(4) << y.zvtx[
k0]
581 <<
" with sump=" << sumpmin <<
" rho*nt =" << y.rho[
k0] *
nt << endl;
585 y.removeItem(
k0, tks);
586 set_vtx_range(
beta, tks, y);
597 const unsigned int nv = y.getSize();
599 for (
unsigned int k = 0;
k < nv;
k++) {
603 for (
unsigned int i = 0;
i <
nt;
i++) {
609 y.zvtx[
k] = sumwz / sumw;
613 for (
unsigned int i = 0;
i <
nt;
i++) {
614 double dx = tks.
zpca[
i] - y.zvtx[
k];
619 double Tc = 2. *
a /
b;
627 if (DEBUGLEVEL > 0) {
628 std::cout <<
"DAClusterizerInZ_vect.beta0: Tc = " << T0 << std::endl;
630 std::cout <<
"DAClusterizerInZ_vect.beta0: nstep = " << coolingsteps << std::endl;
634 if (T0 > 1. / betamax) {
637 return betamax *
std::pow(coolingFactor_, coolingsteps);
640 return betamax * coolingFactor_;
652 unsigned int nv = y.getSize();
656 std::vector<std::pair<double, unsigned int>>
critical;
657 for (
unsigned int k = 0;
k < nv;
k++) {
658 double Tc = 2 * y.swE[
k] / y.sw[
k];
666 std::stable_sort(
critical.begin(),
critical.end(), std::greater<std::pair<double, unsigned int>>());
671 for (
unsigned int ic = 0; ic <
critical.size(); ic++) {
675 double p1 = 0, z1 = 0,
w1 = 0;
676 double p2 = 0,
z2 = 0,
w2 = 0;
677 for (
unsigned int i = 0;
i <
nt;
i++) {
678 if (tks.
sum_Z[
i] > 1.e-100) {
680 double tl = tks.
zpca[
i] < y.zvtx[
k] ? 1. : 0.;
685 if (std::fabs(
arg) < 20) {
686 double t = local_exp(-
arg);
692 double w =
p * tks.
dz2[
i];
694 z1 +=
w * tl * tks.
zpca[
i];
714 if ((
k > 0) && (z1 < (0.6 * y.zvtx[
k] + 0.4 * y.zvtx[
k - 1]))) {
715 z1 = 0.6 * y.zvtx[
k] + 0.4 * y.zvtx[
k - 1];
717 if ((
k + 1 < nv) && (
z2 > (0.6 * y.zvtx[
k] + 0.4 * y.zvtx[
k + 1]))) {
718 z2 = 0.6 * y.zvtx[
k] + 0.4 * y.zvtx[
k + 1];
723 if (DEBUGLEVEL > 1) {
724 if (std::fabs(y.zvtx[
k] - zdumpcenter_) < zdumpwidth_) {
726 <<
std::fixed << std::setprecision(4) << y.zvtx[
k] <<
" --> " << z1 <<
"," <<
z2 <<
" [" <<
p1 740 double pk1 =
p1 * y.rho[
k] / (
p1 +
p2);
741 double pk2 =
p2 * y.rho[
k] / (
p1 +
p2);
744 y.insertItem(
k, z1, pk1, tks);
751 for (
unsigned int jc = ic; jc <
critical.size(); jc++) {
758 std::cout <<
"warning ! split rejected, too small." << endl;
781 clear_vtx_range(tks, y);
784 double beta = beta0(betamax_, tks, y);
790 thermalize(
beta, tks, y, delta_highT_);
794 double betafreeze = betamax_ *
sqrt(coolingFactor_);
796 while (
beta < betafreeze) {
803 thermalize(
beta, tks, y, delta_highT_);
809 if (DEBUGLEVEL > 0) {
810 std::cout <<
"DAClusterizerInZ_vect::vertices :" 811 <<
"last round of splitting" << std::endl;
815 set_vtx_range(
beta, tks, y);
819 set_vtx_range(
beta, tks, y);
823 unsigned int ntry = 0;
826 thermalize(
beta, tks, y, delta_highT_, rho0);
837 if (DEBUGLEVEL > 0) {
838 std::cout <<
"DAClusterizerInZ_vect::vertices :" 839 <<
"turning on outlier rejection at T=" << 1 /
beta << std::endl;
845 rho0 = y.getSize() > 1 ? 1. / y.getSize() : 1.;
846 for (
unsigned int a = 0;
a < 5;
a++) {
851 thermalize(
beta, tks, y, delta_lowT_, rho0);
855 if (DEBUGLEVEL > 0) {
856 std::cout <<
"DAClusterizerInZ_vect::vertices :" 857 <<
"merging with outlier rejection at T=" << 1 /
beta << std::endl;
865 set_vtx_range(
beta, tks, y);
871 if (DEBUGLEVEL > 0) {
872 std::cout <<
"DAClusterizerInZ_vect::vertices :" 873 <<
"after merging with outlier rejection at T=" << 1 /
beta << std::endl;
880 while (
beta < betapurge_) {
882 thermalize(
beta, tks, y, delta_lowT_, rho0);
887 if (DEBUGLEVEL > 0) {
888 std::cout <<
"DAClusterizerInZ_vect::vertices :" 889 <<
"purging at T=" << 1 /
beta << std::endl;
894 while (purge(y, tks, rho0,
beta)) {
895 thermalize(
beta, tks, y, delta_lowT_, rho0);
900 if (DEBUGLEVEL > 0) {
901 std::cout <<
"DAClusterizerInZ_vect::vertices :" 902 <<
"last cooling T=" << 1 /
beta << std::endl;
907 while (
beta < betastop_) {
909 thermalize(
beta, tks, y, delta_lowT_, rho0);
914 if (DEBUGLEVEL > 0) {
915 std::cout <<
"DAClusterizerInZ_vect::vertices :" 916 <<
"stop cooling at T=" << 1 /
beta << std::endl;
924 set_vtx_range(
beta, tks, y);
925 const unsigned int nv = y.getSize();
926 for (
unsigned int k = 0;
k < nv;
k++) {
933 const auto z_sum_init = rho0 * local_exp(-
beta * dzCutOff_ * dzCutOff_);
934 std::vector<std::vector<unsigned int>> vtx_track_indices(nv);
935 for (
unsigned int i = 0;
i <
nt;
i++) {
936 const auto kmin = tks.
kmin[
i];
937 const auto kmax = tks.
kmax[
i];
938 for (
auto k = kmin;
k < kmax;
k++) {
942 local_exp_list_range(y.exp_arg, y.exp, kmin, kmax);
944 tks.
sum_Z[
i] = z_sum_init;
945 for (
auto k = kmin;
k < kmax;
k++) {
946 tks.
sum_Z[
i] += y.rho[
k] * y.exp[
k];
948 const double invZ = tks.
sum_Z[
i] > 1
e-100 ? 1. / tks.
sum_Z[
i] : 0.0;
950 for (
auto k = kmin;
k < kmax;
k++) {
951 double p = y.rho[
k] * y.exp[
k] * invZ;
952 if (
p > mintrkweight_) {
954 vtx_track_indices[
k].push_back(
i);
961 GlobalError dummyError(0.01, 0, 0.01, 0., 0., 0.01);
962 for (
unsigned int k = 0;
k < nv;
k++) {
963 if (!vtx_track_indices[
k].
empty()) {
965 vector<reco::TransientTrack> vertexTracks;
966 for (
auto i : vtx_track_indices[
k]) {
967 vertexTracks.push_back(*(tks.
tt[
i]));
978 vector<reco::TransientTrack> sorted_tracks;
979 vector<pair<float, float>> vertices_tot;
980 for (
unsigned int i = 0;
i <
tracks.size();
i++) {
981 sorted_tracks.push_back(
tracks[
i]);
987 return (
a.stateAtBeamLine().trackStateAtPCA()).
position().z() <
988 (
b.stateAtBeamLine().trackStateAtPCA()).
position().z();
991 unsigned int nBlocks = (
unsigned int)std::floor(sorted_tracks.size() / (block_size_ * (1 - overlap_frac_)));
995 <<
"Warning nBlocks was 0 with ntracks = " << sorted_tracks.size() <<
" block_size = " << block_size_
996 <<
" and overlap fraction = " << overlap_frac_ <<
". Setting nBlocks = 1";
999 vector<reco::TransientTrack> block_tracks;
1000 unsigned int begin = (
unsigned int)(
block * block_size_ * (1 - overlap_frac_));
1001 unsigned int end = (
unsigned int)
std::min(begin + block_size_, (
unsigned int)sorted_tracks.size());
1002 for (
unsigned int i = begin;
i < end;
i++) {
1003 block_tracks.push_back(sorted_tracks[
i]);
1005 if (block_tracks.empty()) {
1010 std::cout <<
"Running vertices_in_blocks on" << std::endl;
1011 std::cout <<
"- block no." <<
block <<
" on " << nBlocks <<
" blocks " << std::endl;
1012 std::cout <<
"- block track size: " << sorted_tracks.size() <<
" - block size: " << block_size_ << std::endl;
1023 clear_vtx_range(tks, y);
1026 beta = beta0(betamax_, tks, y);
1032 thermalize(
beta, tks, y, delta_highT_);
1036 double betafreeze = betamax_ *
sqrt(coolingFactor_);
1037 while (
beta < betafreeze) {
1044 thermalize(
beta, tks, y, delta_highT_);
1050 if (DEBUGLEVEL > 0) {
1051 std::cout <<
"DAClusterizerInZSubCluster_vect::vertices :" 1052 <<
"last round of splitting" << std::endl;
1056 set_vtx_range(
beta, tks, y);
1060 set_vtx_range(
beta, tks, y);
1064 unsigned int ntry = 0;
1067 thermalize(
beta, tks, y, delta_highT_, rho0);
1078 if (DEBUGLEVEL > 0) {
1079 std::cout <<
"DAClusterizerInZSubCluster_vect::vertices :" 1080 <<
"turning on outlier rejection at T=" << 1 /
beta << std::endl;
1085 if (dzCutOff_ > 0) {
1086 rho0 = y.getSize() > 1 ? 1. / y.getSize() : 1.;
1087 for (
unsigned int a = 0;
a < 5;
a++) {
1092 thermalize(
beta, tks, y, delta_lowT_, rho0);
1096 if (DEBUGLEVEL > 0) {
1097 std::cout <<
"DAClusterizerInZSubCluster_vect::vertices :" 1098 <<
"merging with outlier rejection at T=" << 1 /
beta << std::endl;
1106 set_vtx_range(
beta, tks, y);
1112 if (DEBUGLEVEL > 0) {
1113 std::cout <<
"DAClusterizerInZSubCluster_vect::vertices :" 1114 <<
"after merging with outlier rejection at T=" << 1 /
beta << std::endl;
1121 while (
beta < betapurge_) {
1123 thermalize(
beta, tks, y, delta_lowT_, rho0);
1128 if (DEBUGLEVEL > 0) {
1129 std::cout <<
"DAClusterizerInZSubCluster_vect::vertices :" 1130 <<
"purging at T=" << 1 /
beta << std::endl;
1135 while (purge(y, tks, rho0,
beta)) {
1136 thermalize(
beta, tks, y, delta_lowT_, rho0);
1141 if (DEBUGLEVEL > 0) {
1142 std::cout <<
"DAClusterizerInZSubCluster_vect::vertices :" 1143 <<
"last cooling T=" << 1 /
beta << std::endl;
1148 while (
beta < betastop_) {
1150 thermalize(
beta, tks, y, delta_lowT_, rho0);
1155 if (DEBUGLEVEL > 0) {
1156 std::cout <<
"DAClusterizerInZSubCluster_vect::vertices :" 1157 <<
"stop cooling at T=" << 1 /
beta << std::endl;
1163 for (
unsigned int ivertex = 0; ivertex < y.getSize(); ivertex++) {
1164 if (y.zvtx_vec[ivertex] != 0 && y.rho_vec[ivertex] != 0) {
1165 vertices_tot.push_back(pair(y.zvtx_vec[ivertex], y.rho_vec[ivertex]));
1167 std::cout <<
"Found new vertex " << y.zvtx_vec[ivertex] <<
" , " << y.rho_vec[ivertex] << std::endl;
1175 [](
const pair<float, float>&
a,
const pair<float, float>&
b) ->
bool {
return a.first <
b.first; });
1179 const unsigned int nv = vertices_tot.size();
1180 const unsigned int nt = tracks_tot.getSize();
1182 for (
auto itrack = 0
U; itrack <
nt; ++itrack) {
1183 double zrange =
max(sel_zrange_ /
sqrt(
beta * tracks_tot.dz2[itrack]), zrange_min_);
1185 double zmin = tracks_tot.zpca[itrack] -
zrange;
1186 unsigned int kmin =
min(nv - 1, tracks_tot.kmin[itrack]);
1189 while ((kmin > 0) && (vertices_tot[kmin - 1].
first >
zmin)) {
1193 while ((kmin < (nv - 1)) && (vertices_tot[kmin].
first <
zmin)) {
1198 double zmax = tracks_tot.zpca[itrack] +
zrange;
1199 unsigned int kmax =
min(nv - 1, tracks_tot.kmax[itrack] - 1);
1203 while ((kmax < (nv - 1)) && (vertices_tot[kmax + 1].
first <
zmax)) {
1207 while ((kmax > 0) && (vertices_tot[kmax].
first >
zmax)) {
1213 tracks_tot.kmin[itrack] = kmin;
1214 tracks_tot.kmax[itrack] = kmax + 1;
1216 tracks_tot.kmin[itrack] =
max(0
U,
min(kmin, kmax));
1217 tracks_tot.kmax[itrack] =
min(nv,
max(kmin, kmax) + 1);
1221 rho0 = nv > 1 ? 1. / nv : 1.;
1222 const auto z_sum_init = rho0 * local_exp(-
beta * dzCutOff_ * dzCutOff_);
1224 std::vector<std::vector<unsigned int>> vtx_track_indices(nv);
1225 for (
unsigned int i = 0;
i <
nt;
i++) {
1226 const auto kmin = tracks_tot.kmin[
i];
1227 const auto kmax = tracks_tot.kmax[
i];
1229 unsigned int iMax = 10000;
1230 float sum_Z = z_sum_init;
1231 for (
auto k = kmin;
k < kmax;
k++) {
1232 float v_exp = local_exp(-
beta * Eik(tracks_tot.zpca[
i], vertices_tot[
k].first, tracks_tot.dz2[
i]));
1233 sum_Z += vertices_tot[
k].second * v_exp;
1235 double invZ = sum_Z > 1
e-100 ? 1. / sum_Z : 0.0;
1236 for (
auto k = kmin;
k < kmax && invZ != 0.0;
k++) {
1237 float v_exp = local_exp(-
beta * Eik(tracks_tot.zpca[
i], vertices_tot[
k].first, tracks_tot.dz2[
i]));
1238 double p = vertices_tot[
k].second * v_exp * invZ;
1239 if (
p > p_max &&
p > mintrkweight_) {
1244 if (iMax < vtx_track_indices.size())
1248 for (
auto itrack = 0
U; itrack <
nt; ++itrack) {
1249 std::cout <<
"itrack " << itrack <<
" , " << tracks_tot.kmin[itrack] <<
" , " << tracks_tot.kmax[itrack]
1256 GlobalError dummyError(0.01, 0, 0.01, 0., 0., 0.01);
1257 for (
unsigned int k = 0;
k < nv;
k++) {
1258 if (!vtx_track_indices[
k].
empty()) {
1260 vector<reco::TransientTrack> vertexTracks;
1261 for (
auto i : vtx_track_indices[
k]) {
1262 vertexTracks.push_back(*(tracks_tot.tt[
i]));
1264 std::cout << y.zvtx[
k] <<
"," << (*tks.tt[
i]).stateAtBeamLine().trackStateAtPCA().position().z() << std::endl;
1276 const vector<reco::TransientTrack>&
tracks)
const {
1277 vector<vector<reco::TransientTrack>>
clusters;
1279 vector<TransientVertex>
pv;
1280 if (runInBlocks_ and (block_size_ <
tracks.size()))
1286 if (DEBUGLEVEL > 0) {
1287 std::cout <<
"###################################################" << endl;
1288 std::cout <<
"# vectorized DAClusterizerInZ_vect::clusterize nt=" <<
tracks.size() << endl;
1289 std::cout <<
"# DAClusterizerInZ_vect::clusterize pv.size=" <<
pv.size() << endl;
1290 std::cout <<
"###################################################" << endl;
1299 vector<reco::TransientTrack> aCluster =
pv.begin()->originalTracks();
1301 for (
auto k =
pv.begin() + 1;
k !=
pv.end();
k++) {
1304 if (aCluster.size() > 1) {
1309 std::cout <<
" one track cluster at " <<
k->position().z() <<
" suppressed" << std::endl;
1314 for (
unsigned int i = 0;
i <
k->originalTracks().size();
i++) {
1315 aCluster.push_back(
k->originalTracks()[
i]);
1326 const unsigned int nv = y.getSize();
1331 update(
beta, tks_local, y_local, rho0,
true);
1333 std::vector<unsigned int> iz;
1334 for (
unsigned int j = 0;
j <
nt;
j++) {
1337 std::sort(iz.begin(), iz.end(), [tks](
unsigned int a,
unsigned int b) {
return tks.
zpca[
a] < tks.
zpca[
b]; });
1339 std::cout <<
"-----DAClusterizerInZ::dump ----" << nv <<
" clusters " << std::endl;
1341 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1342 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1343 std::cout <<
" " << setw(3) << ivertex <<
" ";
1349 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1350 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1355 <<
"T=" << setw(15) << 1. /
beta <<
" Tmin =" << setw(10) << 1. / betamax_
1357 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1358 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1359 double Tc = 2 * y_local.
swE[ivertex] / y_local.
sw[ivertex];
1367 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1368 sumpk += y.rho[ivertex];
1369 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1371 std::cout << setw(8) << setprecision(4) <<
fixed << y.rho[ivertex];
1376 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1377 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1384 double E = 0,
F = 0;
1386 std::cout <<
"---- z +/- dz ip +/-dip pt phi eta weights ----" << endl;
1388 for (
unsigned int i0 = 0; i0 <
nt; i0++) {
1389 unsigned int i = iz[i0];
1393 double tz = tks.
zpca[
i];
1395 if (std::fabs(tz - zdumpcenter_) > zdumpwidth_)
1397 std::cout << setw(4) <<
i <<
")" << setw(8) <<
fixed << setprecision(4) << tz <<
" +/-" << setw(6)
1399 if ((tks.
tt[
i] ==
nullptr)) {
1416 .pixelBarrelLayersWithMeasurement();
1417 std::cout << setw(1) << tks.
tt[
i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
1419 << tks.
tt[
i]->track().hitPattern().trackerLayersWithMeasurement() -
1420 tks.
tt[
i]->track().hitPattern().pixelLayersWithMeasurement()
1426 std::cout << setw(8) <<
IP.value() <<
"+/-" << setw(6) <<
IP.error();
1427 std::cout <<
" " << setw(6) << setprecision(2) << tks.
tt[
i]->track().pt() * tks.
tt[
i]->track().charge();
1428 std::cout <<
" " << setw(5) << setprecision(2) << tks.
tt[
i]->track().phi() <<
" " << setw(5) << setprecision(2)
1429 << tks.
tt[
i]->track().eta();
1433 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1434 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1439 double p = y.rho[ivertex] * local_exp(-
beta * Eik(tks.
zpca[
i], y.zvtx[ivertex], tks.
dz2[
i])) / tks.
sum_Z[
i];
1441 std::cout << setw(8) << setprecision(3) <<
p;
1445 E +=
p * Eik(tks.
zpca[
i], y.zvtx[ivertex], tks.
dz2[
i]);
1451 std::cout <<
" ( " << std::setw(3) << tks.
kmin[
i] <<
"," << std::setw(3) << tks.
kmax[
i] - 1 <<
" ) ";
1455 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1456 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1457 std::cout <<
" " << setw(3) << ivertex <<
" ";
1463 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1464 if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1470 <<
"T=" << 1 /
beta <<
" E=" << E <<
" n=" << y.getSize() <<
" F= " <<
F << endl
1471 <<
"----------" << endl;
1477 desc.addUntracked<
double>(
"zdumpcenter", 0.);
1478 desc.addUntracked<
double>(
"zdumpwidth", 20.);
1479 desc.add<
double>(
"d0CutOff", 3.0);
1480 desc.add<
double>(
"Tmin", 2.0);
1481 desc.add<
double>(
"delta_lowT", 0.001);
1482 desc.add<
double>(
"zmerge", 0.01);
1483 desc.add<
double>(
"dzCutOff", 3.0);
1484 desc.add<
double>(
"Tpurge", 2.0);
1485 desc.add<
int>(
"convergence_mode", 0);
1486 desc.add<
double>(
"delta_highT", 0.01);
1487 desc.add<
double>(
"Tstop", 0.5);
1488 desc.add<
double>(
"coolingFactor", 0.6);
1489 desc.add<
double>(
"vertexSize", 0.006);
1490 desc.add<
double>(
"uniquetrkweight", 0.8);
1491 desc.add<
double>(
"uniquetrkminp", 0.0);
1492 desc.add<
double>(
"zrange", 4.0);
1493 desc.add<
bool>(
"runInBlocks",
false);
1494 desc.add<
unsigned int>(
"block_size", 512);
1495 desc.add<
double>(
"overlap_frac", 0.5);
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
void addItemSorted(double new_zpca, double new_dz2, const reco::TransientTrack *new_tt, double new_tkwt)
T getParameter(std::string const &) const
bool merge(vertex_t &y, track_t &tks, double &beta) const
common ppss p3p6s2 common epss epspn46 common const1 w2
void clear_vtx_range(track_t >racks, vertex_t &gvertices) const
constexpr bool isNotFinite(T x)
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
bool purge(vertex_t &, track_t &, double &, const double) const
double *__restrict__ sum_Z
static void fillPSetDescription(edm::ParameterSetDescription &desc)
double *__restrict__ tkwt
std::vector< const reco::TransientTrack * > tt
std::vector< double > dz2_vec
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const override
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const
std::vector< double > tkwt_vec
T getUntrackedParameter(std::string const &, T const &) const
U second(std::pair< T, U > const &p)
std::vector< unsigned int > kmin
double *__restrict__ zpca
double update(double beta, track_t >racks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
std::vector< unsigned int > kmax
double *__restrict__ exp_arg
int merge(int argc, char *argv[])
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0, const double rho0=0.) const
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
double *__restrict__ zvtx
Abs< T >::type abs(const T &t)
unsigned int getSize() const
unsigned int getSize() const
std::vector< TransientVertex > vertices_in_blocks(const std::vector< reco::TransientTrack > &tracks) const
def split(sequence, size)
std::vector< double > sum_Z_vec
DAClusterizerInZ_vect(const edm::ParameterSet &conf)
weight_default_t w1[2000]
unsigned int thermalize(double beta, track_t >racks, vertex_t &gvertices, const double delta_max, const double rho0=0.) const
void set_vtx_range(double beta, track_t >racks, vertex_t &gvertices) const
deadvectors [0] push_back({0.0175431, 0.538005, 6.80997, 13.29})
static int position[264][3]
std::vector< double > zpca_vec
Log< level::Warning, false > LogWarning
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
void verify(const vertex_t &v, const track_t &tks, unsigned int nv=999999, unsigned int nt=999999) const
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const
Power< A, B >::type pow(const A &a, const B &b)