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]);
67 for (
auto it = tracks.begin(); it!= tracks.end(); it++){
69 double t_z = ((*it).stateAtBeamLine().trackStateAtPCA()).
position().z();
70 double phi=((*it).stateAtBeamLine().trackStateAtPCA()).momentum().phi();
71 double tantheta =
std::tan( ((*it).stateAtBeamLine().trackStateAtPCA()).momentum().theta());
80 (*it).stateAtBeamLine().transverseImpactParameter();
85 tks.
AddItem(t_z, t_dz2, &(*it), t_pi);
90 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Track count " << tks.
GetSize() << std::endl;
99 double Eik(
double t_z,
double k_z,
double t_dz2) {
100 return std::pow(t_z - k_z, 2) / t_dz2;
105 vertex_t & gvertices,
bool useRho0,
double & rho0)
const {
111 const unsigned int nt = gtracks.
GetSize();
112 const unsigned int nv = gvertices.
GetSize();
125 auto kernel_calc_exp_arg = [
beta, nv ] (
const unsigned int itrack,
128 const double track_z = tracks._z[itrack];
129 const double botrack_dz2 = -beta/tracks._dz2[itrack];
132 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
133 auto mult_res = track_z - vertices._z[ivertex];
134 vertices._ei_cache[ivertex] = botrack_dz2 * ( mult_res * mult_res );
138 auto kernel_add_Z = [ nv, Z_init ] (
vertex_t const& vertices) ->
double
140 double ZTemp = Z_init;
141 for (
unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
142 ZTemp += vertices._pk[ivertex] * vertices._ei[ivertex];
147 auto kernel_calc_normalization = [
beta, nv ] (
const unsigned int track_num,
150 auto tmp_trk_pi = tks_vec._pi[track_num];
151 auto o_trk_Z_sum = 1./tks_vec._Z_sum[track_num];
152 auto o_trk_dz2 = 1./tks_vec._dz2[track_num];
153 auto tmp_trk_z = tks_vec._z[track_num];
154 auto obeta = -1./
beta;
157 for (
unsigned int k = 0;
k < nv; ++
k) {
158 y_vec._se[
k] += y_vec._ei[
k] * (tmp_trk_pi* o_trk_Z_sum);
159 auto w = y_vec._pk[
k] * y_vec._ei[
k] * (tmp_trk_pi*o_trk_Z_sum *o_trk_dz2);
161 y_vec._swz[
k] +=
w * tmp_trk_z;
162 y_vec._swE[
k] +=
w * y_vec._ei_cache[
k]*obeta;
167 for (
auto ivertex = 0U; ivertex < nv; ++ivertex) {
168 gvertices.
_se[ivertex] = 0.0;
169 gvertices.
_sw[ivertex] = 0.0;
170 gvertices.
_swz[ivertex] = 0.0;
171 gvertices.
_swE[ivertex] = 0.0;
178 Z_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
182 for (
auto itrack = 0U; itrack <
nt; ++itrack) {
183 kernel_calc_exp_arg(itrack, gtracks, gvertices);
186 gtracks.
_Z_sum[itrack] = kernel_add_Z(gvertices);
190 sumpi += gtracks.
_pi[itrack];
192 if (gtracks.
_Z_sum[itrack] > 0) {
193 kernel_calc_normalization(itrack, gtracks, gvertices);
198 auto kernel_calc_z = [ sumpi, nv,
this, useRho0 ] (
vertex_t & vertices ) ->
double {
202 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex ) {
203 if (vertices._sw[ivertex] > 0) {
204 auto znew = vertices._swz[ ivertex ] / vertices._sw[ ivertex ];
206 delta +=
std::pow( vertices._z[ ivertex ] - znew, 2 );
207 vertices._z[ ivertex ] = znew;
211 edm::LogInfo(
"sumw") <<
"invalid sum of weights in fit: " << vertices._sw[ivertex] << endl;
212 if (this->verbose_) {
213 LogDebug(
"DAClusterizerinZ_vectorized") <<
" a cluster melted away ? pk=" << vertices._pk[ ivertex ] <<
" sumw="
214 << vertices._sw[ivertex] << endl;
221 auto osumpi = 1./sumpi;
222 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex )
223 vertices._pk[ ivertex ] = vertices._pk[ ivertex ] * vertices._se[ ivertex ] * osumpi;
228 delta += kernel_calc_z(gvertices);
240 const unsigned int nv = y.
GetSize();
246 for (
unsigned int k = 0; (
k + 1) < nv;
k++) {
247 if (fabs(y.
_z[
k + 1] - y.
_z[
k]) < 2.e-2) {
250 double Tc=2*swE/(y.
_sw[
k]+y.
_sw[
k+1]);
277 const unsigned int nv = y.
GetSize();
279 if (nv < 2)
return false;
281 for (
unsigned int k = 0; (
k + 1) < nv;
k++) {
291 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Merging vertices k = " <<
k << std::endl;
304 const unsigned int nv = y.
GetSize();
311 unsigned int k0 = nv;
313 for (
unsigned int k = 0;
k < nv;
k++) {
317 double pmax = y.
_pk[
k] / (y.
_pk[
k] + rho0 * local_exp(-beta * dzCutOff_* dzCutOff_));
319 for (
unsigned int i = 0;
i <
nt;
i++) {
325 if ((p > 0.9 * pmax) && (tks.
_pi[
i] > 0)) {
331 if ((nUnique < 2) && (sump < sumpmin)) {
339 LogDebug(
"DAClusterizerinZ_vectorized") <<
"eliminating prototype at " << y.
_z[
k0] <<
" with sump="
356 const unsigned int nv = y.
GetSize();
358 for (
unsigned int k = 0;
k < nv;
k++) {
363 for (
unsigned int i = 0;
i <
nt;
i++) {
365 sumwz += w * tks.
_z[
i];
368 y.
_z[
k] = sumwz / sumw;
372 for (
unsigned int i = 0;
i <
nt;
i++) {
373 double dx = tks.
_z[
i] - y.
_z[
k];
378 double Tc = 2. * a /
b;
379 if (Tc > T0) T0 = Tc;
382 if (T0 > 1. / betamax) {
386 return betamax / coolingFactor_;
402 std::vector<std::pair<double, unsigned int> >
critical;
403 for(
unsigned int k=0;
k<nv;
k++){
406 critical.push_back( make_pair(Tc,
k));
409 if (critical.size()==0)
return false;
410 std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int> >() );
416 for(
unsigned int ic=0; ic<critical.size(); ic++){
417 unsigned int k=critical[ic].second;
419 double p1=0, z1=0, w1=0;
420 double p2=0, z2=0,
w2=0;
421 for(
unsigned int i=0;
i<
nt;
i++){
425 if(tks.
_z[
i]<y.
_z[k]){
426 p1+=
p; z1+=w*tks.
_z[
i]; w1+=
w;
432 if(w1>0){ z1=z1/w1;}
else{z1=y.
_z[
k]-
epsilon;}
436 if( ( k > 0 ) && ( y.
_z[k-1]>=z1 ) ){ z1=0.5*(y.
_z[
k]+y.
_z[k-1]); }
437 if( ( k+1 < nv) && ( y.
_z[k+1]<=z2 ) ){ z2=0.5*(y.
_z[
k]+y.
_z[k+1]); }
442 double pk1=p1*y.
_pk[
k]/(p1+
p2);
443 double pk2=p2*y.
_pk[
k]/(p1+
p2);
450 for(
unsigned int jc=ic; jc<critical.size(); jc++){
451 if (critical[jc].
second>k) {critical[jc].second++;}
462 const unsigned int nv = y.
GetSize();
469 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Before Split "<< std::endl;
473 for (
unsigned int k = 0;
k < nv;
k++) {
475 ( (
k == 0) || ( y.
_z[
k - 1] < (y.
_z[
k] - zsep)) ) &&
476 ( ((
k + 1) == nv) || ( y.
_z[
k + 1] > (y.
_z[
k] + zsep)) ) )
482 double new_pk = 0.5 * y.
pk[
k];
488 else if ( (y1.
GetSize() == 0 ) ||
505 LogDebug(
"DAClusterizerinZ_vectorized") <<
"After split " << std::endl;
512 vector<TransientVertex>
517 unsigned int nt = tracks.size();
520 vector<TransientVertex> clusters;
521 if (tks.
GetSize() == 0)
return clusters;
532 double beta = beta0(betamax_, tks, y);
533 if ( verbose_)
LogDebug(
"DAClusterizerinZ_vectorized") <<
"Beta0 is " << beta << std::endl;
536 while ((
update(beta, tks, y,
false, rho0) > 1.
e-6) &&
537 (niter++ < maxIterations_)) {}
540 while (beta < betamax_) {
542 update(beta, tks,y,
false, rho0);
543 while(
merge(y,beta)){
update(beta, tks,y,
false, rho0);}
545 beta=beta/coolingFactor_;
547 beta=beta/coolingFactor_;
553 while ((
update(beta, tks, y,
false, rho0) > 1.
e-6) &&
554 (niter++ < maxIterations_)) {}
560 update(beta, tks,y,
false, rho0);
561 while(
merge(y,beta)){
update(beta, tks,y,
false, rho0);}
563 while(
split(beta, tks,y) && (ntry++<10) ){
565 while((
update(beta, tks,y,
false, rho0)>1.
e-6) && (niter++ < maxIterations_)){}
567 update(beta, tks,y,
false, rho0);
571 while(
merge(y,beta)){
update(beta, tks,y,
false, rho0);}
576 LogDebug(
"DAClusterizerinZ_vectorized") <<
"dump after 1st merging " << endl;
577 dump(beta, y, tks, 2);
586 while ((
update(beta, tks, y,
true, rho0) > 1.
e-8) && (niter++ < maxIterations_)) {}
588 LogDebug(
"DAClusterizerinZ_vectorized") <<
"rho0=" << rho0 <<
" niter=" << niter << endl;
589 dump(beta, y, tks, 2);
595 LogDebug(
"DAClusterizerinZ_vectorized") <<
"dump after 2nd merging " << endl;
596 dump(beta, y, tks, 2);
600 while (beta <= betastop_) {
601 while (purge(y, tks, rho0, beta)) {
603 while ((
update(beta, tks, y,
true, rho0) > 1.
e-6) && (niter++ < maxIterations_)) {}
605 beta /= coolingFactor_;
607 while ((
update(beta, tks, y,
true, rho0) > 1.
e-6) && (niter++< maxIterations_)) {}
611 LogDebug(
"DAClusterizerinZ_vectorized") <<
"Final result, rho0=" << rho0 << endl;
612 dump(beta, y, tks, 2);
619 const unsigned int nv = y.
GetSize();
621 for (
unsigned int i = 0;
i <
nt;
i++) {
622 tks.
_Z_sum[
i] = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
623 for (
unsigned int k = 0;
k < nv;
k++) {
628 for (
unsigned int k = 0;
k < nv;
k++) {
631 vector<reco::TransientTrack> vertexTracks;
632 for (
unsigned int i = 0;
i <
nt;
i++) {
635 double p = y.
_pk[
k] * local_exp(-beta * Eik(tks.
_z[
i], y.
_z[
k],
637 if ((tks.
_pi[
i] > 0) && (p > 0.5)) {
639 vertexTracks.push_back(*(tks.
tt[
i]));
645 clusters.push_back(v);
653 const vector<reco::TransientTrack> &
tracks)
const {
656 std::cout <<
"###################################################" << endl;
657 std::cout <<
"# vectorized DAClusterizerInZ_vect::clusterize nt=" << tracks.size() << endl;
658 std::cout <<
"###################################################" << endl;
661 vector<vector<reco::TransientTrack> > clusters;
662 vector<TransientVertex> && pv = vertices(tracks);
665 LogDebug(
"DAClusterizerinZ_vectorized") <<
"# DAClusterizerInZ::clusterize pv.size=" << pv.size()
668 if (pv.size() == 0) {
673 vector<reco::TransientTrack> aCluster = pv.begin()->originalTracks();
675 for (
auto k = pv.begin() + 1;
k != pv.end();
k++) {
678 clusters.push_back(aCluster);
681 for (
unsigned int i = 0;
i <
k->originalTracks().size();
i++) {
682 aCluster.push_back(
k->originalTracks()[
i]);
686 clusters.emplace_back(std::move(aCluster));
699 const unsigned int nv = y.
GetSize();
702 LogDebug(
"DAClusterizerinZ_vectorized") <<
"-----DAClusterizerInZ::dump ----" << endl;
703 LogDebug(
"DAClusterizerinZ_vectorized") <<
"beta=" << beta <<
" betamax= " << betamax_ << endl;
704 LogDebug(
"DAClusterizerinZ_vectorized") <<
" z= ";
705 LogDebug(
"DAClusterizerinZ_vectorized") << setprecision(4);
706 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
707 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << fixed << y.
_z[ivertex];
709 LogDebug(
"DAClusterizerinZ_vectorized") << endl <<
"T=" << setw(15) << 1. / beta
711 LogDebug(
"DAClusterizerinZ_vectorized") << endl
714 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
715 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << setprecision(3) << fixed << y.
_pk[ivertex];
716 sumpk += y.
_pk[ivertex];
718 LogDebug(
"DAClusterizerinZ_vectorized") << endl;
722 LogDebug(
"DAClusterizerinZ_vectorized") << endl;
723 LogDebug(
"DAClusterizerinZ_vectorized")
724 <<
"---- z +/- dz ip +/-dip pt phi eta weights ----"
726 LogDebug(
"DAClusterizerinZ_vectorized") << setprecision(4);
727 for (
unsigned int i = 0;
i <
nt;
i++) {
731 double tz = tks.
_z[
i];
732 LogDebug(
"DAClusterizerinZ_vectorized") << setw(3) <<
i <<
")" << setw(8) << fixed << setprecision(4)
733 << tz <<
" +/-" << setw(6) <<
sqrt(tks.
_dz2[
i]);
736 LogDebug(
"DAClusterizerinZ_vectorized") <<
" *";
738 LogDebug(
"DAClusterizerinZ_vectorized") <<
" ";
740 if (tks.
tt[
i]->track().hitPattern().hasValidHitInFirstPixelBarrel()) {
741 LogDebug(
"DAClusterizerinZ_vectorized") <<
"+";
743 LogDebug(
"DAClusterizerinZ_vectorized") <<
"-";
745 LogDebug(
"DAClusterizerinZ_vectorized") << setw(1)
746 << tks.
tt[
i]->track().hitPattern().pixelBarrelLayersWithMeasurement();
747 LogDebug(
"DAClusterizerinZ_vectorized") << setw(1)
748 << tks.
tt[
i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
749 LogDebug(
"DAClusterizerinZ_vectorized") << setw(1) << hex
750 << tks.
tt[
i]->track().hitPattern().trackerLayersWithMeasurement()
751 - tks.
tt[
i]->track().hitPattern().pixelLayersWithMeasurement()
753 LogDebug(
"DAClusterizerinZ_vectorized") <<
"=" << setw(1) << hex
758 tks.
tt[
i]->stateAtBeamLine().transverseImpactParameter();
759 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << IP.
value() <<
"+/-" << setw(6) << IP.
error();
760 LogDebug(
"DAClusterizerinZ_vectorized") <<
" " << setw(6) << setprecision(2)
761 << tks.
tt[
i]->track().pt() * tks.
tt[
i]->track().charge();
762 LogDebug(
"DAClusterizerinZ_vectorized") <<
" " << setw(5) << setprecision(2)
763 << tks.
tt[
i]->track().phi() <<
" " << setw(5)
764 << setprecision(2) << tks.
tt[
i]->track().eta();
767 for (
unsigned int ivertex = 0; ivertex < nv; ++ ivertex) {
772 LogDebug(
"DAClusterizerinZ_vectorized") << setw(8) << setprecision(3) <<
p;
774 LogDebug(
"DAClusterizerinZ_vectorized") <<
" . ";
776 E += p * Eik(tks.
_z[
i], y.
_z[ivertex], tks.
_dz2[
i]);
779 LogDebug(
"DAClusterizerinZ_vectorized") <<
" ";
782 LogDebug(
"DAClusterizerinZ_vectorized") << endl;
784 LogDebug(
"DAClusterizerinZ_vectorized") << endl <<
"T=" << 1 / beta <<
" E=" << E <<
" n=" << y.
GetSize()
785 <<
" 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)
Sin< T >::type sin(const T &t)
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
static int position[TOTALCHAMBERS][3]
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)
Cos< T >::type cos(const T &t)
void InsertItem(unsigned int i, double new_z, double new_pk)
double *__restrict__ _swz
Tan< T >::type tan(const T &t)
Abs< T >::type abs(const T &t)
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
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks, const int verbosity=0) const
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