11 #include "vdt/vdtMath.h"
32 coolingFactor_ = conf.
getParameter<
double> (
"coolingFactor");
35 coolingFactor_=-coolingFactor_;
42 LogDebug(
"DAClusterizerinZ_vectorized") <<
"DAClusterizerInZ: invalid Tmin" << Tmin
43 <<
" reset do default " << 1. / betamax_ << endl;
51 inline double local_exp(
double const& inp) {
52 return vdt::fast_exp( inp );
55 inline void local_exp_list(
double const * __restrict__ arg_inp,
double * __restrict__ arg_out,
const int arg_arr_size) {
56 for(
int i=0;
i!=arg_arr_size; ++
i ) arg_out[
i]=vdt::fast_exp(arg_inp[
i]);
66 LogDebug(
"DAClusterizerinZ_vectorized") <<
"input values";
68 for (
auto it = tracks.begin(); it!= tracks.end(); it++){
69 if (!(*it).isValid())
continue;
71 double t_z = ((*it).stateAtBeamLine().trackStateAtPCA()).
position().z();
73 auto const & t_mom = (*it).stateAtBeamLine().trackStateAtPCA().momentum();
84 (*it).stateAtBeamLine().transverseImpactParameter();
88 LogTrace(
"DAClusterizerinZ_vectorized") << t_z <<
' '<< t_dz2 <<
' '<< t_pi;
89 tks.
AddItem(t_z, t_dz2, &(*it), t_pi);
94 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Track count " << tks.
GetSize() << std::endl;
103 double Eik(
double t_z,
double k_z,
double t_dz2) {
104 return std::pow(t_z - k_z, 2) * t_dz2;
109 vertex_t & gvertices,
bool useRho0,
double & rho0)
const {
115 const unsigned int nt = gtracks.
GetSize();
116 const unsigned int nv = gvertices.
GetSize();
129 auto kernel_calc_exp_arg = [
beta, nv ] (
const unsigned int itrack,
132 const double track_z = tracks._z[itrack];
133 const double botrack_dz2 = -beta*tracks._dz2[itrack];
136 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
137 auto mult_res = track_z -
vertices._z[ivertex];
138 vertices._ei_cache[ivertex] = botrack_dz2 * ( mult_res * mult_res );
144 double ZTemp = Z_init;
145 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
151 auto kernel_calc_normalization = [
beta, nv ] (
const unsigned int track_num,
154 auto tmp_trk_pi = tks_vec._pi[track_num];
155 auto o_trk_Z_sum = 1./tks_vec._Z_sum[track_num];
156 auto o_trk_dz2 = tks_vec._dz2[track_num];
157 auto tmp_trk_z = tks_vec._z[track_num];
158 auto obeta = -1./
beta;
161 for (
unsigned int k = 0;
k < nv; ++
k) {
162 y_vec._se[
k] += y_vec._ei[
k] * (tmp_trk_pi* o_trk_Z_sum);
163 auto w = y_vec._pk[
k] * y_vec._ei[
k] * (tmp_trk_pi*o_trk_Z_sum *o_trk_dz2);
165 y_vec._swz[
k] +=
w * tmp_trk_z;
166 y_vec._swE[
k] +=
w * y_vec._ei_cache[
k]*obeta;
171 for (
auto ivertex = 0U; ivertex < nv; ++ivertex) {
172 gvertices.
_se[ivertex] = 0.0;
173 gvertices.
_sw[ivertex] = 0.0;
174 gvertices.
_swz[ivertex] = 0.0;
175 gvertices.
_swE[ivertex] = 0.0;
182 Z_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
186 for (
auto itrack = 0U; itrack <
nt; ++itrack) {
187 kernel_calc_exp_arg(itrack, gtracks, gvertices);
190 gtracks.
_Z_sum[itrack] = kernel_add_Z(gvertices);
194 sumpi += gtracks.
_pi[itrack];
196 if (gtracks.
_Z_sum[itrack] > 0) {
197 kernel_calc_normalization(itrack, gtracks, gvertices);
202 auto kernel_calc_z = [ sumpi, nv,
this, useRho0 ] (
vertex_t &
vertices ) ->
double {
206 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex ) {
216 if (this->verbose_) {
217 LogDebug(
"DAClusterizerinZ_vectorized") <<
" a cluster melted away ? pk=" <<
vertices._pk[ ivertex ] <<
" sumw="
225 auto osumpi = 1./sumpi;
226 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex )
232 delta += kernel_calc_z(gvertices);
244 const unsigned int nv = y.
GetSize();
250 for (
unsigned int k = 0; (
k + 1) < nv;
k++) {
251 if (fabs(y.
_z[
k + 1] - y.
_z[
k]) < 2.e-2) {
254 double Tc=2*swE/(y.
_sw[
k]+y.
_sw[
k+1]);
281 const unsigned int nv = y.
GetSize();
283 if (nv < 2)
return false;
285 for (
unsigned int k = 0; (
k + 1) < nv;
k++) {
295 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Merging vertices k = " <<
k << std::endl;
308 const unsigned int nv = y.
GetSize();
315 unsigned int k0 = nv;
317 for (
unsigned int k = 0;
k < nv;
k++) {
321 double pmax = y.
_pk[
k] / (y.
_pk[
k] + rho0 * local_exp(-beta * dzCutOff_* dzCutOff_));
323 for (
unsigned int i = 0;
i <
nt;
i++) {
329 if ((p > 0.9 * pmax) && (tks.
_pi[
i] > 0)) {
335 if ((nUnique < 2) && (sump < sumpmin)) {
343 LogDebug(
"DAClusterizerinZ_vectorized") <<
"eliminating prototype at " << y.
_z[
k0] <<
" with sump="
360 const unsigned int nv = y.
GetSize();
362 for (
unsigned int k = 0;
k < nv;
k++) {
367 for (
unsigned int i = 0;
i <
nt;
i++) {
369 sumwz += w * tks.
_z[
i];
372 y.
_z[
k] = sumwz / sumw;
376 for (
unsigned int i = 0;
i <
nt;
i++) {
377 double dx = tks.
_z[
i] - y.
_z[
k];
382 double Tc = 2. * a /
b;
383 if (Tc > T0) T0 = Tc;
386 if (T0 > 1. / betamax) {
390 return betamax / coolingFactor_;
406 std::vector<std::pair<double, unsigned int> >
critical;
407 for(
unsigned int k=0;
k<nv;
k++){
410 critical.push_back( make_pair(Tc,
k));
413 if (critical.size()==0)
return false;
414 std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int> >() );
420 for(
unsigned int ic=0; ic<critical.size(); ic++){
421 unsigned int k=critical[ic].second;
423 double p1=0, z1=0, w1=0;
424 double p2=0, z2=0,
w2=0;
425 for(
unsigned int i=0;
i<
nt;
i++){
429 if(tks.
_z[
i]<y.
_z[k]){
430 p1+=
p; z1+=w*tks.
_z[
i]; w1+=
w;
436 if(w1>0){ z1=z1/w1;}
else{z1=y.
_z[
k]-
epsilon;}
440 if( ( k > 0 ) && ( y.
_z[k-1]>=z1 ) ){ z1=0.5*(y.
_z[
k]+y.
_z[k-1]); }
441 if( ( k+1 < nv) && ( y.
_z[k+1]<=z2 ) ){ z2=0.5*(y.
_z[
k]+y.
_z[k+1]); }
446 double pk1=p1*y.
_pk[
k]/(p1+
p2);
447 double pk2=p2*y.
_pk[
k]/(p1+
p2);
454 for(
unsigned int jc=ic; jc<critical.size(); jc++){
455 if (critical[jc].
second>k) {critical[jc].second++;}
466 const unsigned int nv = y.
GetSize();
473 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Before Split "<< std::endl;
477 for (
unsigned int k = 0;
k < nv;
k++) {
479 ( (
k == 0) || ( y.
_z[
k - 1] < (y.
_z[
k] - zsep)) ) &&
480 ( ((
k + 1) == nv) || ( y.
_z[
k + 1] > (y.
_z[
k] + zsep)) ) )
486 double new_pk = 0.5 * y.
pk[
k];
492 else if ( (y1.
GetSize() == 0 ) ||
509 LogDebug(
"DAClusterizerinZ_vectorized") <<
"After split " << std::endl;
516 vector<TransientVertex>
525 if (tks.
GetSize() == 0)
return clusters;
536 double beta = beta0(betamax_, tks, y);
537 if ( verbose_)
LogDebug(
"DAClusterizerinZ_vectorized") <<
"Beta0 is " << beta << std::endl;
540 while ((
update(beta, tks, y,
false, rho0) > 1.
e-6) &&
541 (niter++ < maxIterations_)) {}
544 while (beta < betamax_) {
546 update(beta, tks,y,
false, rho0);
547 while(
merge(y,beta)){
update(beta, tks,y,
false, rho0);}
549 beta=beta/coolingFactor_;
551 beta=beta/coolingFactor_;
557 while ((
update(beta, tks, y,
false, rho0) > 1.
e-6) &&
558 (niter++ < maxIterations_)) {}
564 update(beta, tks,y,
false, rho0);
565 while(
merge(y,beta)){
update(beta, tks,y,
false, rho0);}
567 while(
split(beta, tks,y) && (ntry++<10) ){
569 while((
update(beta, tks,y,
false, rho0)>1.
e-6) && (niter++ < maxIterations_)){}
571 update(beta, tks,y,
false, rho0);
575 while(
merge(y,beta)){
update(beta, tks,y,
false, rho0);}
580 LogDebug(
"DAClusterizerinZ_vectorized") <<
"dump after 1st merging " << endl;
581 dump(beta, y, tks, 2);
590 while ((
update(beta, tks, y,
true, rho0) > 1.
e-8) && (niter++ < maxIterations_)) {}
592 LogDebug(
"DAClusterizerinZ_vectorized") <<
"rho0=" << rho0 <<
" niter=" << niter << endl;
593 dump(beta, y, tks, 2);
599 LogDebug(
"DAClusterizerinZ_vectorized") <<
"dump after 2nd merging " << endl;
600 dump(beta, y, tks, 2);
604 while (beta <= betastop_) {
605 while (purge(y, tks, rho0, beta)) {
607 while ((
update(beta, tks, y,
true, rho0) > 1.
e-6) && (niter++ < maxIterations_)) {}
609 beta /= coolingFactor_;
611 while ((
update(beta, tks, y,
true, rho0) > 1.
e-6) && (niter++< maxIterations_)) {}
615 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Final result, rho0=" << rho0 << endl;
616 dump(beta, y, tks, 2);
623 const unsigned int nv = y.
GetSize();
624 for (
unsigned int k = 0;
k < nv;
k++)
627 for (
unsigned int i = 0;
i <
nt;
i++)
628 tks.
_Z_sum[
i] = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
631 for (
unsigned int k = 0;
k < nv;
k++) {
632 for (
unsigned int i = 0;
i <
nt;
i++)
637 for (
unsigned int k = 0;
k < nv;
k++) {
640 vector<reco::TransientTrack> vertexTracks;
641 for (
unsigned int i = 0;
i <
nt;
i++) {
644 double p = y.
_pk[
k] * local_exp(-beta * Eik(tks.
_z[
i], y.
_z[
k],
646 if ((tks.
_pi[
i] > 0) && (p > 0.5)) {
647 vertexTracks.push_back(*(tks.
tt[
i]));
653 clusters.push_back(v);
661 const vector<reco::TransientTrack> &
tracks)
const {
664 std::cout <<
"###################################################" << endl;
665 std::cout <<
"# vectorized DAClusterizerInZ_vect::clusterize nt=" << tracks.size() << endl;
666 std::cout <<
"###################################################" << endl;
669 vector<vector<reco::TransientTrack> >
clusters;
670 vector<TransientVertex> &&
pv =
vertices(tracks);
673 LogDebug(
"DAClusterizerinZ_vectorized") <<
"# DAClusterizerInZ::clusterize pv.size=" << pv.size()
676 if (pv.size() == 0) {
681 vector<reco::TransientTrack> aCluster = pv.begin()->originalTracks();
683 for (
auto k = pv.begin() + 1;
k != pv.end();
k++) {
686 clusters.push_back(aCluster);
689 for (
unsigned int i = 0;
i <
k->originalTracks().size();
i++) {
690 aCluster.push_back(
k->originalTracks()[
i]);
694 clusters.emplace_back(
std::move(aCluster));
707 const unsigned int nv = y.
GetSize();
710 LogDebug(
"DAClusterizerinZ_vectorized") <<
"-----DAClusterizerInZ::dump ----" << endl;
711 LogDebug(
"DAClusterizerinZ_vectorized") <<
"beta=" << beta <<
" betamax= " << betamax_ << endl;
712 LogDebug(
"DAClusterizerinZ_vectorized") <<
" z= ";
713 LogDebug(
"DAClusterizerinZ_vectorized") << setprecision(4);
714 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
715 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << fixed << y.
_z[ivertex];
717 LogDebug(
"DAClusterizerinZ_vectorized") << endl <<
"T=" << setw(15) << 1. / beta
719 LogDebug(
"DAClusterizerinZ_vectorized") << endl
722 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
723 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << setprecision(3) << fixed << y.
_pk[ivertex];
724 sumpk += y.
_pk[ivertex];
726 LogDebug(
"DAClusterizerinZ_vectorized") << endl;
730 LogDebug(
"DAClusterizerinZ_vectorized") << endl;
731 LogDebug(
"DAClusterizerinZ_vectorized")
732 <<
"---- z +/- dz ip +/-dip pt phi eta weights ----"
734 LogDebug(
"DAClusterizerinZ_vectorized") << setprecision(4);
735 for (
unsigned int i = 0;
i <
nt;
i++) {
739 double tz = tks.
_z[
i];
740 LogDebug(
"DAClusterizerinZ_vectorized") << setw(3) <<
i <<
")" << setw(8) << fixed << setprecision(4)
741 << tz <<
" +/-" << setw(6) <<
sqrt(1./tks.
_dz2[
i]);
744 LogDebug(
"DAClusterizerinZ_vectorized") <<
" *";
746 LogDebug(
"DAClusterizerinZ_vectorized") <<
" ";
748 if (tks.
tt[
i]->track().hitPattern().hasValidHitInFirstPixelBarrel()) {
749 LogDebug(
"DAClusterizerinZ_vectorized") <<
"+";
751 LogDebug(
"DAClusterizerinZ_vectorized") <<
"-";
753 LogDebug(
"DAClusterizerinZ_vectorized") << setw(1)
754 << tks.
tt[
i]->track().hitPattern().pixelBarrelLayersWithMeasurement();
755 LogDebug(
"DAClusterizerinZ_vectorized") << setw(1)
756 << tks.
tt[
i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
757 LogDebug(
"DAClusterizerinZ_vectorized") << setw(1) << hex
758 << tks.
tt[
i]->track().hitPattern().trackerLayersWithMeasurement()
759 - tks.
tt[
i]->track().hitPattern().pixelLayersWithMeasurement()
761 LogDebug(
"DAClusterizerinZ_vectorized") <<
"=" << setw(1) << hex
766 tks.
tt[
i]->stateAtBeamLine().transverseImpactParameter();
767 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << IP.
value() <<
"+/-" << setw(6) << IP.
error();
768 LogDebug(
"DAClusterizerinZ_vectorized") <<
" " << setw(6) << setprecision(2)
769 << tks.
tt[
i]->track().pt() * tks.
tt[
i]->track().charge();
770 LogDebug(
"DAClusterizerinZ_vectorized") <<
" " << setw(5) << setprecision(2)
771 << tks.
tt[
i]->track().phi() <<
" " << setw(5)
772 << setprecision(2) << tks.
tt[
i]->track().eta();
775 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
780 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << setprecision(3) <<
p;
782 LogDebug(
"DAClusterizerinZ_vectorized") <<
" . ";
784 E += p * Eik(tks.
_z[
i], y.
_z[ivertex], tks.
_dz2[
i]);
787 LogDebug(
"DAClusterizerinZ_vectorized") <<
" ";
790 LogDebug(
"DAClusterizerinZ_vectorized") << endl;
792 LogDebug(
"DAClusterizerinZ_vectorized") << endl <<
"T=" << 1 / beta <<
" E=" << E <<
" n=" << y.
GetSize()
793 <<
" F= " <<
F << endl <<
"----------" << endl;
double *__restrict__ _ei_cache
void AddItem(double new_z, double new_dz2, const reco::TransientTrack *new_tt, double new_pi)
double *__restrict__ _swE
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
common ppss p3p6s2 common epss epspn46 common const1 w2
void splitAll(vertex_t &y) const
void AddItem(double new_z, double new_pk)
bool split(const double beta, track_t &t, vertex_t &y) const
track_t fill(const std::vector< reco::TransientTrack > &tracks) const
unsigned int GetSize() const
std::vector< const reco::TransientTrack * > tt
bool merge(vertex_t &) const
double *__restrict__ _dz2
std::vector< std::vector< reco::TransientTrack > > clusterize(const std::vector< reco::TransientTrack > &tracks) const
U second(std::pair< T, U > const &p)
void InsertItem(unsigned int i, double new_z, double new_pk)
double *__restrict__ _swz
Abs< T >::type abs(const T &t)
double BeamWidthX() const
beam width X
double update(double beta, track_t >racks, vertex_t &gvertices, bool useRho0, double &rho0) const
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
void RemoveItem(unsigned int i)
DAClusterizerInZ_vect(const edm::ParameterSet &conf)
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0) const
double BeamWidthY() const
beam width Y
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks, const int verbosity=0) const
static int position[264][3]
unsigned int GetSize() const
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
bool purge(vertex_t &, track_t &, double &, const double) const
Power< A, B >::type pow(const A &a, const B &b)
double *__restrict__ _Z_sum