CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Static Public Member Functions | Private Attributes
DAClusterizerInZ_vect Class Referencefinal

#include <DAClusterizerInZ_vect.h>

Inheritance diagram for DAClusterizerInZ_vect:
TrackClusterizerInZ

Classes

struct  track_t
 
struct  vertex_t
 

Public Member Functions

double beta0 (const double betamax, track_t const &tks, vertex_t const &y) const
 
void clear_vtx_range (track_t &gtracks, vertex_t &gvertices) const
 
std::vector< std::vector< reco::TransientTrack > > clusterize (const std::vector< reco::TransientTrack > &tracks) const override
 
 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 rho0=0.) const
 
track_t fill (const std::vector< reco::TransientTrack > &tracks) const
 
std::vector< TransientVertexfill_vertices (double beta, double rho0, track_t &tracks, vertex_t &vertices) const
 
bool merge (vertex_t &y, track_t &tks, double &beta) const
 
bool purge (vertex_t &, track_t &, double &, const double) const
 
void set_vtx_range (double beta, track_t &gtracks, vertex_t &gvertices) const
 
bool split (const double beta, track_t &t, vertex_t &y, double threshold=1.) const
 
unsigned int thermalize (double beta, track_t &gtracks, vertex_t &gvertices, const double delta_max, const double rho0=0.) const
 
double update (double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
 
void verify (const vertex_t &v, const track_t &tks, unsigned int nv=999999, unsigned int nt=999999) const
 
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &tracks) const override
 
std::vector< TransientVertexvertices_in_blocks (const std::vector< reco::TransientTrack > &tracks) const
 
std::vector< TransientVertexvertices_no_blocks (const std::vector< reco::TransientTrack > &tracks) const
 
- Public Member Functions inherited from TrackClusterizerInZ
 TrackClusterizerInZ ()=default
 
 TrackClusterizerInZ (const edm::ParameterSet &conf)
 
virtual ~TrackClusterizerInZ ()=default
 

Static Public Member Functions

static void fillPSetDescription (edm::ParameterSetDescription &desc)
 

Private Attributes

double betamax_
 
double betapurge_
 
double betastop_
 
unsigned int block_size_
 
unsigned int convergence_mode_
 
double coolingFactor_
 
double d0CutOff_
 
double delta_highT_
 
double delta_lowT_
 
double dzCutOff_
 
unsigned int maxIterations_
 
double mintrkweight_
 
double overlap_frac_
 
bool runInBlocks_
 
double sel_zrange_
 
double uniquetrkminp_
 
double uniquetrkweight_
 
double vertexSize_
 
double zdumpcenter_
 
double zdumpwidth_
 
double zmerge_
 
const double zrange_min_ = 0.1
 

Detailed Description

Description: separates event tracks into clusters along the beam line

Version which auto-vectorizes with gcc 4.6 or newer

Definition at line 21 of file DAClusterizerInZ_vect.h.

Constructor & Destructor Documentation

◆ DAClusterizerInZ_vect()

DAClusterizerInZ_vect::DAClusterizerInZ_vect ( const edm::ParameterSet conf)

Definition at line 20 of file DAClusterizerInZ_vect.cc.

References gather_cfg::cout, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), SiStripPI::max, HLT_2024v12_cff::Tmin, HLT_2024v12_cff::Tpurge, and HLT_2024v12_cff::Tstop.

20  {
21  // hardcoded parameters
22  maxIterations_ = 1000;
23  mintrkweight_ = 0.5;
24 
25  // configurable debug output
26 #ifdef DEBUG
27  zdumpcenter_ = conf.getUntrackedParameter<double>("zdumpcenter", 0.);
28  zdumpwidth_ = conf.getUntrackedParameter<double>("zdumpwidth", 20.);
29 #endif
30 
31  // configurable parameters
32  double Tmin = conf.getParameter<double>("Tmin");
33  double Tpurge = conf.getParameter<double>("Tpurge");
34  double Tstop = conf.getParameter<double>("Tstop");
35  vertexSize_ = conf.getParameter<double>("vertexSize");
36  coolingFactor_ = conf.getParameter<double>("coolingFactor");
37  d0CutOff_ = conf.getParameter<double>("d0CutOff");
38  dzCutOff_ = conf.getParameter<double>("dzCutOff");
39  uniquetrkweight_ = conf.getParameter<double>("uniquetrkweight");
40  uniquetrkminp_ = conf.getParameter<double>("uniquetrkminp");
41  zmerge_ = conf.getParameter<double>("zmerge");
42  sel_zrange_ = conf.getParameter<double>("zrange");
43  convergence_mode_ = conf.getParameter<int>("convergence_mode");
44  delta_lowT_ = conf.getParameter<double>("delta_lowT");
45  delta_highT_ = conf.getParameter<double>("delta_highT");
46  runInBlocks_ = conf.getParameter<bool>("runInBlocks");
47  block_size_ = conf.getParameter<unsigned int>("block_size");
48  overlap_frac_ = conf.getParameter<double>("overlap_frac");
49 
50 #ifdef DEBUG
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;
56  std::cout << "DAClusterizerinZ_vect: Tpurge = " << Tpurge << 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 
67  std::cout << "DAClusterizerinZ_vect: run in blocks = " << runInBlocks_ << std::endl;
68  std::cout << "DAClusterizerinZ_vect: block_size = " << block_size_ << std::endl;
69  std::cout << "DAClusterizerinZ_vect: overlap_fraction = " << overlap_frac_ << std::endl;
70  std::cout << "DAClusterizerinZ_vect: DEBUGLEVEL " << DEBUGLEVEL << std::endl;
71 #endif
72 
73  if (convergence_mode_ > 1) {
74  edm::LogWarning("DAClusterizerinZ_vect")
75  << "DAClusterizerInZ_vect: invalid convergence_mode " << convergence_mode_ << " reset to default " << 0;
77  }
78 
79  if (Tmin == 0) {
80  betamax_ = 1.0;
81  edm::LogWarning("DAClusterizerinZ_vect")
82  << "DAClusterizerInZ_vect: invalid Tmin " << Tmin << " reset to default " << 1. / betamax_;
83  } else {
84  betamax_ = 1. / Tmin;
85  }
86 
87  if ((Tpurge > Tmin) || (Tpurge == 0)) {
88  edm::LogWarning("DAClusterizerinZ_vect")
89  << "DAClusterizerInZ_vect: invalid Tpurge " << Tpurge << " set to " << Tmin;
90  Tpurge = Tmin;
91  }
92  betapurge_ = 1. / Tpurge;
93 
94  if ((Tstop > Tpurge) || (Tstop == 0)) {
95  edm::LogWarning("DAClusterizerinZ_vect")
96  << "DAClusterizerInZ_vect: invalid Tstop " << Tstop << " set to " << max(1., Tpurge);
97  Tstop = max(1., Tpurge);
98  }
99  betastop_ = 1. / Tstop;
100 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
T getUntrackedParameter(std::string const &, T const &) const
Log< level::Warning, false > LogWarning

Member Function Documentation

◆ beta0()

double DAClusterizerInZ_vect::beta0 ( const double  betamax,
track_t const &  tks,
vertex_t const &  y 
) const

Definition at line 596 of file DAClusterizerInZ_vect.cc.

References a, b, gather_cfg::cout, PVValHelper::dx, DAClusterizerInZ_vect::track_t::dz2, DAClusterizerInZ_vect::track_t::getSize(), mps_fire::i, createfilelist::int, dqmdumpme::k, dqm-mbProfile::log, nt, funct::pow(), DAClusterizerInZ_vect::track_t::tkwt, w(), and DAClusterizerInZ_vect::track_t::zpca.

596  {
597  double T0 = 0; // max Tc for beta=0
598  // estimate critical temperature from beta=0 (T=inf)
599  const unsigned int nt = tks.getSize();
600  const unsigned int nv = y.getSize();
601 
602  for (unsigned int k = 0; k < nv; k++) {
603  // vertex fit at T=inf
604  double sumwz = 0;
605  double sumw = 0;
606  for (unsigned int i = 0; i < nt; i++) {
607  double w = tks.tkwt[i] * tks.dz2[i];
608  sumwz += w * tks.zpca[i];
609  sumw += w;
610  }
611 
612  y.zvtx[k] = sumwz / sumw;
613 
614  // estimate Tcrit
615  double a = 0, b = 0;
616  for (unsigned int i = 0; i < nt; i++) {
617  double dx = tks.zpca[i] - y.zvtx[k];
618  double w = tks.tkwt[i] * tks.dz2[i];
619  a += w * std::pow(dx, 2) * tks.dz2[i];
620  b += w;
621  }
622  double Tc = 2. * a / b; // the critical temperature of this vertex
623 
624  if (Tc > T0)
625  T0 = Tc;
626 
627  } // vertex loop (normally there should be only one vertex at beta=0)
628 
629 #ifdef DEBUG
630  if (DEBUGLEVEL > 0) {
631  std::cout << "DAClusterizerInZ_vect.beta0: Tc = " << T0 << std::endl;
632  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
633  std::cout << "DAClusterizerInZ_vect.beta0: nstep = " << coolingsteps << std::endl;
634  }
635 #endif
636 
637  if (T0 > 1. / betamax) {
638  int coolingsteps = 1 - int(std::log(T0 * betamax) / std::log(coolingFactor_));
639 
640  return betamax * std::pow(coolingFactor_, coolingsteps);
641  } else {
642  // ensure at least one annealing step
643  return betamax * coolingFactor_;
644  }
645 }
T w() const
int nt
Definition: AMPTWrapper.h:42
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ clear_vtx_range()

void DAClusterizerInZ_vect::clear_vtx_range ( track_t gtracks,
vertex_t gvertices 
) const

Definition at line 275 of file DAClusterizerInZ_vect.cc.

References DAClusterizerInZ_vect::track_t::getSize(), DAClusterizerInZ_vect::vertex_t::getSize(), DAClusterizerInZ_vect::track_t::kmax, DAClusterizerInZ_vect::track_t::kmin, nt, and mitigatedMETSequence_cff::U.

275  {
276  const unsigned int nt = gtracks.getSize();
277  const unsigned int nv = gvertices.getSize();
278  for (auto itrack = 0U; itrack < nt; ++itrack) {
279  gtracks.kmin[itrack] = 0;
280  gtracks.kmax[itrack] = nv;
281  }
282 }
int nt
Definition: AMPTWrapper.h:42

◆ clusterize()

vector< vector< reco::TransientTrack > > DAClusterizerInZ_vect::clusterize ( const std::vector< reco::TransientTrack > &  tracks) const
overridevirtual

Implements TrackClusterizerInZ.

Definition at line 1338 of file DAClusterizerInZ_vect.cc.

References funct::abs(), bsc_activity_cfg::clusters, gather_cfg::cout, mps_fire::i, dqmdumpme::k, eostools::move(), position, DiMuonV_cfg::tracks, and AlignmentTracksFromVertexSelector_cfi::vertices.

1339  {
1340  vector<vector<reco::TransientTrack>> clusters;
1341  vector<TransientVertex>&& pv = vertices(tracks);
1342 
1343 #ifdef DEBUG
1344  if (DEBUGLEVEL > 0) {
1345  std::cout << "###################################################" << endl;
1346  std::cout << "# vectorized DAClusterizerInZ_vect::clusterize nt=" << tracks.size() << endl;
1347  std::cout << "# DAClusterizerInZ_vect::clusterize pv.size=" << pv.size() << endl;
1348  std::cout << "###################################################" << endl;
1349  }
1350 #endif
1351 
1352  if (pv.empty()) {
1353  return clusters;
1354  }
1355 
1356  // fill into clusters and merge
1357  vector<reco::TransientTrack> aCluster = pv.begin()->originalTracks();
1358 
1359  for (auto k = pv.begin() + 1; k != pv.end(); k++) {
1360  if (std::abs(k->position().z() - (k - 1)->position().z()) > (2 * vertexSize_)) {
1361  // close a cluster
1362  if (aCluster.size() > 1) {
1363  clusters.push_back(aCluster);
1364  }
1365 #ifdef DEBUG
1366  else {
1367  std::cout << " one track cluster at " << k->position().z() << " suppressed" << std::endl;
1368  }
1369 #endif
1370  aCluster.clear();
1371  }
1372  for (unsigned int i = 0; i < k->originalTracks().size(); i++) {
1373  aCluster.push_back(k->originalTracks()[i]);
1374  }
1375  }
1376  clusters.emplace_back(std::move(aCluster));
1377 
1378  return clusters;
1379 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static int position[264][3]
Definition: ReadPGInfo.cc:289
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override
def move(src, dest)
Definition: eostools.py:511

◆ dump()

void DAClusterizerInZ_vect::dump ( const double  beta,
const vertex_t y,
const track_t tks,
const int  verbosity = 0,
const double  rho0 = 0. 
) const

Definition at line 1381 of file DAClusterizerInZ_vect.cc.

References a, b, HLT_2024v12_cff::beta, gather_cfg::cout, TauDecayModes::dec, DAClusterizerInZ_vect::track_t::dz2, F(), alignBH_cfg::fixed, DAClusterizerInZ_vect::track_t::getSize(), reco::TrackBase::highPurity, mps_fire::i, listHistos::IP, dqmiolumiharvest::j, DAClusterizerInZ_vect::track_t::kmax, DAClusterizerInZ_vect::track_t::kmin, dqm-mbProfile::log, reco::HitPattern::MISSING_OUTER_HITS, nt, AlCaHLTBitMon_ParallelJobs::p, GeomDetEnumerators::PixelBarrel, jetUpdater_cfi::sort, mathSSE::sqrt(), DAClusterizerInZ_vect::track_t::sum_Z, DAClusterizerInZ_vect::vertex_t::sw, DAClusterizerInZ_vect::vertex_t::swE, DAClusterizerInZ_vect::track_t::tkwt, DAClusterizerInZ_vect::track_t::tt, update, verbosity, and DAClusterizerInZ_vect::track_t::zpca.

1382  {
1383 #ifdef DEBUG
1384  const unsigned int nv = y.getSize();
1385  const unsigned int nt = tks.getSize();
1386  // make a local copies of clusters and tracks to update Tc [ = y_local.swE / y_local.sw ]
1387  vertex_t y_local = y;
1388  track_t tks_local = tks;
1389  update(beta, tks_local, y_local, rho0, true);
1390 
1391  std::vector<unsigned int> iz;
1392  for (unsigned int j = 0; j < nt; j++) {
1393  iz.push_back(j);
1394  }
1395  std::sort(iz.begin(), iz.end(), [tks](unsigned int a, unsigned int b) { return tks.zpca[a] < tks.zpca[b]; });
1396  std::cout << std::endl;
1397  std::cout << "-----DAClusterizerInZ::dump ----" << nv << " clusters " << std::endl;
1398  std::cout << " ";
1399  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1400  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1401  std::cout << " " << setw(3) << ivertex << " ";
1402  }
1403  }
1404  std::cout << endl;
1405  std::cout << " z= ";
1406  std::cout << setprecision(4);
1407  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1408  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1409  std::cout << setw(8) << fixed << y.zvtx[ivertex];
1410  }
1411  }
1412  std::cout << endl
1413  << "T=" << setw(15) << 1. / beta << " Tmin =" << setw(10) << 1. / betamax_
1414  << " Tc= ";
1415  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1416  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1417  double Tc = 2 * y_local.swE[ivertex] / y_local.sw[ivertex];
1418  std::cout << setw(8) << fixed << setprecision(1) << Tc;
1419  }
1420  }
1421  std::cout << endl;
1422 
1423  std::cout << " pk= ";
1424  double sumpk = 0;
1425  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1426  sumpk += y.rho[ivertex];
1427  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1428  continue;
1429  std::cout << setw(8) << setprecision(4) << fixed << y.rho[ivertex];
1430  }
1431  std::cout << endl;
1432 
1433  std::cout << " nt= ";
1434  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1435  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1436  continue;
1437  std::cout << setw(8) << setprecision(1) << fixed << y.rho[ivertex] * nt;
1438  }
1439  std::cout << endl;
1440 
1441  if (verbosity > 0) {
1442  double E = 0, F = 0;
1443  std::cout << endl;
1444  std::cout << "---- z +/- dz ip +/-dip pt phi eta weights ----" << endl;
1445  std::cout << setprecision(4);
1446  for (unsigned int i0 = 0; i0 < nt; i0++) {
1447  unsigned int i = iz[i0];
1448  if (tks.sum_Z[i] > 0) {
1449  F -= std::log(tks.sum_Z[i]) / beta;
1450  }
1451  double tz = tks.zpca[i];
1452 
1453  if (std::fabs(tz - zdumpcenter_) > zdumpwidth_)
1454  continue;
1455  std::cout << setw(4) << i << ")" << setw(8) << fixed << setprecision(4) << tz << " +/-" << setw(6)
1456  << sqrt(1. / tks.dz2[i]);
1457  if ((tks.tt[i] == nullptr)) {
1458  std::cout << " effective track ";
1459  } else {
1460  if (tks.tt[i]->track().quality(reco::TrackBase::highPurity)) {
1461  std::cout << " *";
1462  } else {
1463  std::cout << " ";
1464  }
1465  if (tks.tt[i]->track().hitPattern().hasValidHitInPixelLayer(PixelSubdetector::SubDetector::PixelBarrel, 1)) {
1466  std::cout << "+";
1467  } else {
1468  std::cout << "-";
1469  }
1470  std::cout << setw(1)
1471  << tks.tt[i]
1472  ->track()
1473  .hitPattern()
1474  .pixelBarrelLayersWithMeasurement(); // see DataFormats/TrackReco/interface/HitPattern.h
1475  std::cout << setw(1) << tks.tt[i]->track().hitPattern().pixelEndcapLayersWithMeasurement();
1476  std::cout << setw(1) << hex
1477  << tks.tt[i]->track().hitPattern().trackerLayersWithMeasurement() -
1478  tks.tt[i]->track().hitPattern().pixelLayersWithMeasurement()
1479  << dec;
1480  std::cout << "=" << setw(1) << hex
1481  << tks.tt[i]->track().hitPattern().numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS) << dec;
1482 
1483  Measurement1D IP = tks.tt[i]->stateAtBeamLine().transverseImpactParameter();
1484  std::cout << setw(8) << IP.value() << "+/-" << setw(6) << IP.error();
1485  std::cout << " " << setw(6) << setprecision(2) << tks.tt[i]->track().pt() * tks.tt[i]->track().charge();
1486  std::cout << " " << setw(5) << setprecision(2) << tks.tt[i]->track().phi() << " " << setw(5) << setprecision(2)
1487  << tks.tt[i]->track().eta();
1488  } // not a dummy track
1489 
1490  double sump = 0.;
1491  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1492  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) > zdumpwidth_)
1493  continue;
1494 
1495  if ((tks.tkwt[i] > 0) && (tks.sum_Z[i] > 0)) {
1496  //double p=pik(beta,tks[i],*k);
1497  double p = y.rho[ivertex] * local_exp(-beta * Eik(tks.zpca[i], y.zvtx[ivertex], tks.dz2[i])) / tks.sum_Z[i];
1498  if (p > 0.0001) {
1499  std::cout << setw(8) << setprecision(3) << p;
1500  } else {
1501  std::cout << " . ";
1502  }
1503  E += p * Eik(tks.zpca[i], y.zvtx[ivertex], tks.dz2[i]);
1504  sump += p;
1505  } else {
1506  std::cout << " ";
1507  }
1508  }
1509  std::cout << " ( " << std::setw(3) << tks.kmin[i] << "," << std::setw(3) << tks.kmax[i] - 1 << " ) ";
1510  std::cout << endl;
1511  }
1512  std::cout << " ";
1513  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1514  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1515  std::cout << " " << setw(3) << ivertex << " ";
1516  }
1517  }
1518  std::cout << endl;
1519  std::cout << " z= ";
1520  std::cout << setprecision(4);
1521  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
1522  if (std::fabs(y.zvtx[ivertex] - zdumpcenter_) < zdumpwidth_) {
1523  std::cout << setw(8) << fixed << y.zvtx[ivertex];
1524  }
1525  }
1526  std::cout << endl;
1527  std::cout << endl
1528  << "T=" << 1 / beta << " E=" << E << " n=" << y.getSize() << " F= " << F << endl
1529  << "----------" << endl;
1530  }
1531 #endif
1532 }
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
T sqrt(T t)
Definition: SSEVec.h:19
int nt
Definition: AMPTWrapper.h:42
const int verbosity
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
static uInt32 F(BLOWFISH_CTX *ctx, uInt32 x)
Definition: blowfish.cc:163

◆ fill()

DAClusterizerInZ_vect::track_t DAClusterizerInZ_vect::fill ( const std::vector< reco::TransientTrack > &  tracks) const

Definition at line 177 of file DAClusterizerInZ_vect.cc.

References DAClusterizerInZ_vect::track_t::addItemSorted(), HLT_2024v12_cff::atIP, pwdgSkimBPark_cfi::beamSpot, gather_cfg::cout, geometryDiff::epsilon, DAClusterizerInZ_vect::track_t::extractRaw(), DAClusterizerInZ_vect::track_t::getSize(), edm::isNotFinite(), ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, SiStripPI::min, DAClusterizerInZ_vect::track_t::osumtkwt, position, funct::pow(), AlignmentPI::t_z, and DiMuonV_cfg::tracks.

177  {
178  // prepare track data for clustering
179  track_t tks;
180  double sumtkwt = 0.;
181  for (auto it = tracks.begin(); it != tracks.end(); it++) {
182  if (!(*it).isValid())
183  continue;
184  double t_tkwt = 1.;
185  double t_z = ((*it).stateAtBeamLine().trackStateAtPCA()).position().z();
186  if (std::fabs(t_z) > 1000.)
187  continue;
188  auto const& t_mom = (*it).stateAtBeamLine().trackStateAtPCA().momentum();
189  // get the beam-spot
190  reco::BeamSpot beamspot = (it->stateAtBeamLine()).beamSpot();
191  double t_dz2 = std::pow((*it).track().dzError(), 2) // track errror
192  + (std::pow(beamspot.BeamWidthX() * t_mom.x(), 2) + std::pow(beamspot.BeamWidthY() * t_mom.y(), 2)) *
193  std::pow(t_mom.z(), 2) / std::pow(t_mom.perp2(), 2) // beam spot width
194  + std::pow(vertexSize_, 2); // intrinsic vertex size, safer for outliers and short lived decays
195  t_dz2 = 1. / t_dz2;
196  if (edm::isNotFinite(t_dz2) || t_dz2 < std::numeric_limits<double>::min())
197  continue;
198  if (d0CutOff_ > 0) {
199  Measurement1D atIP = (*it).stateAtBeamLine().transverseImpactParameter(); // error contains beamspot
200  t_tkwt = 1. / (1. + local_exp(std::pow(atIP.value() / atIP.error(), 2) -
201  std::pow(d0CutOff_, 2))); // reduce weight for high ip tracks
203  edm::LogWarning("DAClusterizerinZ_vect") << "rejected track t_tkwt " << t_tkwt;
204  continue; // usually is > 0.99
205  }
206  }
207  tks.addItemSorted(t_z, t_dz2, &(*it), t_tkwt);
208  sumtkwt += t_tkwt;
209  }
210 
211  tks.extractRaw();
212  tks.osumtkwt = sumtkwt > 0 ? 1. / sumtkwt : 0.;
213 
214 #ifdef DEBUG
215  if (DEBUGLEVEL > 0) {
216  std::cout << "Track count (Z) " << tks.getSize() << std::endl;
217  }
218 #endif
219  return tks;
220 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
static int position[264][3]
Definition: ReadPGInfo.cc:289
Log< level::Warning, false > LogWarning
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ fill_vertices()

vector< TransientVertex > DAClusterizerInZ_vect::fill_vertices ( double  beta,
double  rho0,
track_t tracks,
vertex_t vertices 
) const

Definition at line 1240 of file DAClusterizerInZ_vect.cc.

References HLT_2024v12_cff::beta, cms::cuda::bs, bsc_activity_cfg::clusters, DAClusterizerInZ_vect::track_t::dz2, MillePedeFileConverter_cfg::e, relativeConstraints::empty, submitPVResolutionJobs::err, DAClusterizerInZ_vect::track_t::getSize(), mps_fire::i, edm::isNotFinite(), dqmiolumiharvest::j, dqmdumpme::k, DAClusterizerInZ_vect::track_t::kmax, DAClusterizerInZ_vect::track_t::kmin, nt, AlCaHLTBitMon_ParallelJobs::p, funct::pow(), DAClusterizerInZ_vect::track_t::sum_Z, DAClusterizerInZ_vect::track_t::tt, findQualityFiles::v, w(), and DAClusterizerInZ_vect::track_t::zpca.

1240  {
1241  // select significant tracks and use a TransientVertex as a container
1242 
1243  set_vtx_range(beta, tks, y);
1244  const unsigned int nv = y.getSize();
1245  for (unsigned int k = 0; k < nv; k++) {
1246  if (edm::isNotFinite(y.rho[k]) || edm::isNotFinite(y.zvtx[k])) {
1247  y.rho[k] = 0;
1248  y.zvtx[k] = 0;
1249  }
1250  }
1251 
1252  // ensure consistent assignment probabillities and make a hard assignment
1253  const unsigned int nt = tks.getSize();
1254  const auto z_sum_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
1255  std::vector<std::vector<unsigned int>> vtx_track_indices(nv);
1256  std::vector<std::vector<float>> vtx_track_weights(nv);
1257  for (unsigned int i = 0; i < nt; i++) {
1258  const auto kmin = tks.kmin[i];
1259  const auto kmax = tks.kmax[i];
1260  for (auto k = kmin; k < kmax; k++) {
1261  y.exp_arg[k] = -beta * Eik(tks.zpca[i], y.zvtx[k], tks.dz2[i]);
1262  }
1263 
1264  local_exp_list_range(y.exp_arg, y.exp, kmin, kmax);
1265 
1266  tks.sum_Z[i] = z_sum_init;
1267  for (auto k = kmin; k < kmax; k++) {
1268  tks.sum_Z[i] += y.rho[k] * y.exp[k];
1269  }
1270  const double invZ = tks.sum_Z[i] > 1e-100 ? 1. / tks.sum_Z[i] : 0.0;
1271 
1272  double pmax = -1;
1273  unsigned int k_pmax = 0;
1274  for (auto k = kmin; k < kmax; k++) {
1275  double p = y.rho[k] * y.exp[k] * invZ;
1276  if (p > pmax) {
1277  pmax = p;
1278  k_pmax = k;
1279  }
1280  }
1281 
1282  if (pmax > mintrkweight_) {
1283  // assign to the cluster with the highest assignment weight, if it is at least mintrkweight_
1284  vtx_track_indices[k_pmax].push_back(i);
1285  vtx_track_weights[k_pmax].push_back(pmax);
1286  }
1287  }
1288 
1289  // fill transient vertices
1290  // the position is normally not used, probably not optimal when Tstop <> 2, anyway
1291  vector<TransientVertex> clusters;
1292  for (unsigned int k = 0; k < nv; k++) {
1293  double sump = 0;
1294  double sumw = 0;
1295  double sumwp = 0, sumwz = 0;
1296  if (!vtx_track_indices[k].empty()) {
1297  vector<reco::TransientTrack> vertexTracks;
1299  unsigned int j = 0;
1300  for (auto i : vtx_track_indices[k]) {
1301  auto p = vtx_track_weights[k][j];
1302  vertexTracks.push_back(*(tks.tt[i]));
1303  trkWeightMap[vertexTracks[j]] = p;
1304  auto w = p * tks.dz2[i];
1305  sump += p;
1306  sumw += w;
1307  sumwp += w * p;
1308  sumwz += w * tks.zpca[i];
1309  j++;
1310  }
1311  float zerror_squared = 1.; //
1312  if ((sumw > 0) && (sumwp > 0)) {
1313  zerror_squared = sumwp / (sumw * sumw);
1314  y.zvtx[k] = sumwz / sumw;
1315  }
1316 
1317  reco::BeamSpot bs = vertexTracks[0].stateAtBeamLine().beamSpot();
1318  GlobalPoint pos(bs.x(y.zvtx[k]), bs.y(y.zvtx[k]), y.zvtx[k]);
1319  const float xerror_squared = pow(bs.BeamWidthX(), 2);
1320  const float yerror_squared = pow(bs.BeamWidthY(), 2);
1321  GlobalError err(xerror_squared, 0, yerror_squared, 0., 0., zerror_squared);
1322  TransientVertex v(pos, err, vertexTracks, 0, 2 * sump - 3.);
1323  v.weightMap(trkWeightMap);
1324  clusters.push_back(v);
1325  }
1326  }
1327 
1328  return clusters;
1329 }
T w() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
std::map< reco::TransientTrack, float > TransientTrackToFloatMap
int nt
Definition: AMPTWrapper.h:42
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ fillPSetDescription()

void DAClusterizerInZ_vect::fillPSetDescription ( edm::ParameterSetDescription desc)
static

Definition at line 1534 of file DAClusterizerInZ_vect.cc.

References submitPVResolutionJobs::desc.

Referenced by PrimaryVertexProducer::fillDescriptions(), PrimaryVertexValidation::fillDescriptions(), and DAClusterizerInZT_vect::fillPSetDescription().

1534  {
1535  desc.addUntracked<double>("zdumpcenter", 0.);
1536  desc.addUntracked<double>("zdumpwidth", 20.);
1537  desc.add<double>("d0CutOff", 3.0);
1538  desc.add<double>("Tmin", 2.0);
1539  desc.add<double>("delta_lowT", 0.001);
1540  desc.add<double>("zmerge", 0.01);
1541  desc.add<double>("dzCutOff", 3.0);
1542  desc.add<double>("Tpurge", 2.0);
1543  desc.add<int>("convergence_mode", 0);
1544  desc.add<double>("delta_highT", 0.01);
1545  desc.add<double>("Tstop", 0.5);
1546  desc.add<double>("coolingFactor", 0.6);
1547  desc.add<double>("vertexSize", 0.006);
1548  desc.add<double>("uniquetrkweight", 0.8);
1549  desc.add<double>("uniquetrkminp", 0.0);
1550  desc.add<double>("zrange", 4.0);
1551  desc.add<bool>("runInBlocks", false);
1552  desc.add<unsigned int>("block_size", 10000);
1553  desc.add<double>("overlap_frac", 0.0);
1554 }

◆ merge()

bool DAClusterizerInZ_vect::merge ( vertex_t y,
track_t tks,
double &  beta 
) const

Definition at line 468 of file DAClusterizerInZ_vect.cc.

References cms::cuda::assert(), HLT_2024v12_cff::beta, gather_cfg::cout, hlt_jetmet_dqm_QT_fromfile_cfg::critical, alignBH_cfg::fixed, and dqmdumpme::k.

468  {
469  // merge clusters that collapsed or never separated,
470  // only merge if the estimated critical temperature of the merged vertex is below the current temperature
471  // return true if vertices were merged, false otherwise
472  const unsigned int nv = y.getSize();
473 
474  if (nv < 2)
475  return false;
476 
477  // merge the smallest distance clusters first
478  std::vector<std::pair<double, unsigned int>> critical;
479  for (unsigned int k = 0; (k + 1) < nv; k++) {
480  if (std::fabs(y.zvtx[k + 1] - y.zvtx[k]) < zmerge_) {
481  critical.push_back(make_pair(std::fabs(y.zvtx[k + 1] - y.zvtx[k]), k));
482  }
483  }
484  if (critical.empty())
485  return false;
486 
487  std::stable_sort(critical.begin(), critical.end(), std::less<std::pair<double, unsigned int>>());
488 
489  for (unsigned int ik = 0; ik < critical.size(); ik++) {
490  unsigned int k = critical[ik].second;
491  double rho = y.rho[k] + y.rho[k + 1];
492 
493 #ifdef DEBUG
494  assert((k + 1) < nv);
495  if (DEBUGLEVEL > 1) {
496  std::cout << "merging " << fixed << setprecision(4) << y.zvtx[k + 1] << " and " << y.zvtx[k]
497  << " sw = " << y.sw[k] + y.sw[k + 1] << std::endl;
498  }
499 #endif
500 
501  if (rho > 0) {
502  y.zvtx[k] = (y.rho[k] * y.zvtx[k] + y.rho[k + 1] * y.zvtx[k + 1]) / rho;
503  } else {
504  y.zvtx[k] = 0.5 * (y.zvtx[k] + y.zvtx[k + 1]);
505  }
506  y.rho[k] = rho;
507  y.sw[k] += y.sw[k + 1];
508  y.removeItem(k + 1, tks);
509  set_vtx_range(beta, tks, y);
510  y.extractRaw();
511  return true;
512  }
513 
514  return false;
515 }
assert(be >=bs)
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const

◆ purge()

bool DAClusterizerInZ_vect::purge ( vertex_t y,
track_t tks,
double &  rho0,
const double  beta 
) const

Definition at line 517 of file DAClusterizerInZ_vect.cc.

References cms::cuda::assert(), HLT_2024v12_cff::beta, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), gather_cfg::cout, DAClusterizerInZ_vect::track_t::dz2, HLT_2024v12_cff::eps, cms::cuda::for(), DAClusterizerInZ_vect::track_t::getSize(), mps_fire::i, dqmdumpme::k, reco::ParticleMasses::k0, DAClusterizerInZ_vect::track_t::kmax, DAClusterizerInZ_vect::track_t::kmin, nt, AlCaHLTBitMon_ParallelJobs::p, DAClusterizerInZ_vect::track_t::sum_Z, DAClusterizerInZ_vect::track_t::tkwt, and DAClusterizerInZ_vect::track_t::zpca.

517  {
518  constexpr double eps = 1.e-100;
519  // eliminate clusters with only one significant/unique track
520  const unsigned int nv = y.getSize();
521  const unsigned int nt = tks.getSize();
522 
523  if (nv < 2)
524  return false;
525 
526  std::vector<double> sump_v(nv), arg_cache_v(nv), exp_cache_v(nv), pcut_cache_v(nv);
527  std::vector<int> nUnique_v(nv);
528  double* __restrict__ parg_cache;
529  double* __restrict__ pexp_cache;
530  double* __restrict__ ppcut_cache;
531  double* __restrict__ psump;
532  int* __restrict__ pnUnique;
533  int constexpr nunique_min_ = 2;
534 
535  set_vtx_range(beta, tks, y);
536 
537  parg_cache = arg_cache_v.data();
538  pexp_cache = exp_cache_v.data();
539  ppcut_cache = pcut_cache_v.data();
540  psump = sump_v.data();
541  pnUnique = nUnique_v.data();
542 
543  const auto rhoconst = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
544  for (unsigned int k = 0; k < nv; k++) {
545  const double pmax = y.rho[k] / (y.rho[k] + rhoconst);
546  ppcut_cache[k] = uniquetrkweight_ * pmax;
547  }
548 
549  for (unsigned int i = 0; i < nt; i++) {
550  const auto invZ = ((tks.sum_Z[i] > eps) && (tks.tkwt[i] > uniquetrkminp_)) ? 1. / tks.sum_Z[i] : 0.;
551  const auto track_z = tks.zpca[i];
552  const auto botrack_dz2 = -beta * tks.dz2[i];
553  const auto kmin = tks.kmin[i];
554  const auto kmax = tks.kmax[i];
555 
556  for (unsigned int k = kmin; k < kmax; k++) {
557  const auto mult_resz = track_z - y.zvtx[k];
558  parg_cache[k] = botrack_dz2 * (mult_resz * mult_resz);
559  }
560 
561  local_exp_list_range(parg_cache, pexp_cache, kmin, kmax);
562 
563  for (unsigned int k = kmin; k < kmax; k++) {
564  const double p = y.rho[k] * pexp_cache[k] * invZ;
565  psump[k] += p;
566  pnUnique[k] += (p > ppcut_cache[k]) ? 1 : 0;
567  }
568  }
569 
570  double sumpmin = nt;
571  unsigned int k0 = nv;
572  for (unsigned k = 0; k < nv; k++) {
573  if ((pnUnique[k] < nunique_min_) && (psump[k] < sumpmin)) {
574  sumpmin = psump[k];
575  k0 = k;
576  }
577  }
578 
579  if (k0 != nv) {
580 #ifdef DEBUG
581  assert(k0 < y.getSize());
582  if (DEBUGLEVEL > 1) {
583  std::cout << "eliminating prototype at " << std::setw(10) << std::setprecision(4) << y.zvtx[k0]
584  << " with sump=" << sumpmin << " rho*nt =" << y.rho[k0] * nt << " pnUnique=" << pnUnique[k0] << endl;
585  }
586 #endif
587 
588  y.removeItem(k0, tks);
589  set_vtx_range(beta, tks, y);
590  return true;
591  } else {
592  return false;
593  }
594 }
for(int i=first, nt=offsets[nh];i< nt;i+=gridDim.x *blockDim.x)
assert(be >=bs)
int nt
Definition: AMPTWrapper.h:42
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const

◆ set_vtx_range()

void DAClusterizerInZ_vect::set_vtx_range ( double  beta,
track_t gtracks,
vertex_t gvertices 
) const

Definition at line 226 of file DAClusterizerInZ_vect.cc.

References HLT_2024v12_cff::beta, DAClusterizerInZ_vect::track_t::dz2, DAClusterizerInZ_vect::track_t::getSize(), DAClusterizerInZ_vect::vertex_t::getSize(), DAClusterizerInZ_vect::track_t::kmax, DAClusterizerInZ_vect::track_t::kmin, SiStripPI::max, SiStripPI::min, nt, mathSSE::sqrt(), mitigatedMETSequence_cff::U, SiStripMonitorCluster_cfi::zmax, SiStripMonitorCluster_cfi::zmin, DAClusterizerInZ_vect::track_t::zpca, OfflinePixel3DPrimaryVertices_cfi::zrange, and DAClusterizerInZ_vect::vertex_t::zvtx.

226  {
227  const unsigned int nv = gvertices.getSize();
228  const unsigned int nt = gtracks.getSize();
229 
230  if (nv == 0) {
231  edm::LogWarning("DAClusterizerinZ_vect") << "empty cluster list in set_vtx_range";
232  return;
233  }
234 
235  for (auto itrack = 0U; itrack < nt; ++itrack) {
236  double zrange = max(sel_zrange_ / sqrt(beta * gtracks.dz2[itrack]), zrange_min_);
237 
238  double zmin = gtracks.zpca[itrack] - zrange;
239  unsigned int kmin = min(nv - 1, gtracks.kmin[itrack]);
240  // find the smallest vertex_z that is larger than zmin
241  if (gvertices.zvtx[kmin] > zmin) {
242  while ((kmin > 0) && (gvertices.zvtx[kmin - 1] > zmin)) {
243  kmin--;
244  }
245  } else {
246  while ((kmin < (nv - 1)) && (gvertices.zvtx[kmin] < zmin)) {
247  kmin++;
248  }
249  }
250 
251  double zmax = gtracks.zpca[itrack] + zrange;
252  unsigned int kmax = min(nv - 1, gtracks.kmax[itrack] - 1);
253  // note: kmax points to the last vertex in the range, while gtracks.kmax points to the entry BEHIND the last vertex
254  // find the largest vertex_z that is smaller than zmax
255  if (gvertices.zvtx[kmax] < zmax) {
256  while ((kmax < (nv - 1)) && (gvertices.zvtx[kmax + 1] < zmax)) {
257  kmax++;
258  }
259  } else {
260  while ((kmax > 0) && (gvertices.zvtx[kmax] > zmax)) {
261  kmax--;
262  }
263  }
264 
265  if (kmin <= kmax) {
266  gtracks.kmin[itrack] = kmin;
267  gtracks.kmax[itrack] = kmax + 1;
268  } else {
269  gtracks.kmin[itrack] = max(0U, min(kmin, kmax));
270  gtracks.kmax[itrack] = min(nv, max(kmin, kmax) + 1);
271  }
272  }
273 }
T sqrt(T t)
Definition: SSEVec.h:19
int nt
Definition: AMPTWrapper.h:42
Log< level::Warning, false > LogWarning

◆ split()

bool DAClusterizerInZ_vect::split ( const double  beta,
track_t t,
vertex_t y,
double  threshold = 1. 
) const

Definition at line 647 of file DAClusterizerInZ_vect.cc.

References cms::cuda::assert(), HLT_2024v12_cff::beta, ALPAKA_ACCELERATOR_NAMESPACE::brokenline::constexpr(), gather_cfg::cout, hlt_jetmet_dqm_QT_fromfile_cfg::critical, DAClusterizerInZ_vect::track_t::dz2, MillePedeFileConverter_cfg::e, geometryDiff::epsilon, alignBH_cfg::fixed, DAClusterizerInZ_vect::track_t::getSize(), mps_fire::i, dqmdumpme::k, nt, AlCaHLTBitMon_ParallelJobs::p, LaserDQM_cfg::p1, SiStripOfflineCRack_cfg::p2, edm::second(), submitPVValidationJobs::split(), mathSSE::sqrt(), DAClusterizerInZ_vect::track_t::sum_Z, submitPVValidationJobs::t, DiMuonV_cfg::threshold, DAClusterizerInZ_vect::track_t::tkwt, update, w(), w2, testProducerWithPsetDescEmpty_cfi::z2, and DAClusterizerInZ_vect::track_t::zpca.

647  {
648  // split only critical vertices (Tc >~ T=1/beta <==> beta*Tc>~1)
649  // returns true if at least one cluster was split
650 
651  // update the sums needed for Tc, rho0 is never non-zero while splitting is still active
652  update(beta, tks, y, 0., true); // the "true" flag enables Tc evaluation
653 
654  constexpr double epsilon = 1e-3; // minimum split size
655  unsigned int nv = y.getSize();
656 
657  // avoid left-right biases by splitting highest Tc first
658 
659  std::vector<std::pair<double, unsigned int>> critical;
660  for (unsigned int k = 0; k < nv; k++) {
661  double Tc = 2 * y.swE[k] / y.sw[k];
662  if (beta * Tc > threshold) {
663  critical.push_back(make_pair(Tc, k));
664  }
665  }
666  if (critical.empty())
667  return false;
668 
669  std::stable_sort(critical.begin(), critical.end(), std::greater<std::pair<double, unsigned int>>());
670 
671  bool split = false;
672  const unsigned int nt = tks.getSize();
673 
674  for (unsigned int ic = 0; ic < critical.size(); ic++) {
675  unsigned int k = critical[ic].second;
676 
677  // estimate subcluster positions and weight
678  double p1 = 0, z1 = 0, w1 = 0;
679  double p2 = 0, z2 = 0, w2 = 0;
680  for (unsigned int i = 0; i < nt; i++) {
681  if (tks.sum_Z[i] > 1.e-100) {
682  // winner-takes-all, usually overestimates splitting
683  double tl = tks.zpca[i] < y.zvtx[k] ? 1. : 0.;
684  double tr = 1. - tl;
685 
686  // soften it, especially at low T
687  double arg = (tks.zpca[i] - y.zvtx[k]) * sqrt(beta * tks.dz2[i]);
688  if (std::fabs(arg) < 20) {
689  double t = local_exp(-arg);
690  tl = t / (t + 1.);
691  tr = 1 / (t + 1.);
692  }
693 
694  double p = y.rho[k] * tks.tkwt[i] * local_exp(-beta * Eik(tks.zpca[i], y.zvtx[k], tks.dz2[i])) / tks.sum_Z[i];
695  double w = p * tks.dz2[i];
696  p1 += p * tl;
697  z1 += w * tl * tks.zpca[i];
698  w1 += w * tl;
699  p2 += p * tr;
700  z2 += w * tr * tks.zpca[i];
701  w2 += w * tr;
702  }
703  }
704 
705  if (w1 > 0) {
706  z1 = z1 / w1;
707  } else {
708  z1 = y.zvtx[k] - epsilon;
709  }
710  if (w2 > 0) {
711  z2 = z2 / w2;
712  } else {
713  z2 = y.zvtx[k] + epsilon;
714  }
715 
716  // reduce split size if there is not enough room
717  if ((k > 0) && (z1 < (0.6 * y.zvtx[k] + 0.4 * y.zvtx[k - 1]))) {
718  z1 = 0.6 * y.zvtx[k] + 0.4 * y.zvtx[k - 1];
719  }
720  if ((k + 1 < nv) && (z2 > (0.6 * y.zvtx[k] + 0.4 * y.zvtx[k + 1]))) {
721  z2 = 0.6 * y.zvtx[k] + 0.4 * y.zvtx[k + 1];
722  }
723 
724 #ifdef DEBUG
725  assert(k < nv);
726  if (DEBUGLEVEL > 1) {
727  if (std::fabs(y.zvtx[k] - zdumpcenter_) < zdumpwidth_) {
728  std::cout << " T= " << std::setw(8) << 1. / beta << " Tc= " << critical[ic].first << " splitting "
729  << std::fixed << std::setprecision(4) << y.zvtx[k] << " --> " << z1 << "," << z2 << " [" << p1
730  << "," << p2 << "]";
731  if (std::fabs(z2 - z1) > epsilon) {
732  std::cout << std::endl;
733  } else {
734  std::cout << " rejected " << std::endl;
735  }
736  }
737  }
738 #endif
739 
740  // split if the new subclusters are significantly separated
741  if ((z2 - z1) > epsilon) {
742  split = true;
743  double pk1 = p1 * y.rho[k] / (p1 + p2);
744  double pk2 = p2 * y.rho[k] / (p1 + p2);
745  y.zvtx[k] = z2;
746  y.rho[k] = pk2;
747  y.insertItem(k, z1, pk1, tks);
748  if (k == 0)
749  y.extractRaw();
750 
751  nv++;
752 
753  // adjust remaining pointers
754  for (unsigned int jc = ic; jc < critical.size(); jc++) {
755  if (critical[jc].second >= k) {
756  critical[jc].second++;
757  }
758  }
759  } else {
760 #ifdef DEBUG
761  std::cout << "warning ! split rejected, too small." << endl;
762 #endif
763  }
764  }
765 
766  return split;
767 }
common ppss p3p6s2 common epss epspn46 common const1 w2
Definition: inclppp.h:1
T w() const
assert(be >=bs)
A arg
Definition: Factorize.h:31
U second(std::pair< T, U > const &p)
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
T sqrt(T t)
Definition: SSEVec.h:19
int nt
Definition: AMPTWrapper.h:42
bool split(const double beta, track_t &t, vertex_t &y, double threshold=1.) const

◆ thermalize()

unsigned int DAClusterizerInZ_vect::thermalize ( double  beta,
track_t gtracks,
vertex_t gvertices,
const double  delta_max,
const double  rho0 = 0. 
) const

Definition at line 420 of file DAClusterizerInZ_vect.cc.

References funct::abs(), HLT_2024v12_cff::beta, gather_cfg::cout, dumpMFGeometry_cfg::delta, GCP_Ntuples_cfg::dump, dqmdumpme::k, SiStripPI::max, mathSSE::sqrt(), update, and findQualityFiles::v.

421  {
422  unsigned int niter = 0;
423  double delta = 0;
424  double delta_max = delta_lowT_;
425 
426  if (convergence_mode_ == 0) {
427  delta_max = delta_max0;
428  } else if (convergence_mode_ == 1) {
429  delta_max = delta_lowT_ / sqrt(std::max(beta, 1.0));
430  }
431 
432  set_vtx_range(beta, tks, v);
433  double delta_sum_range = 0; // accumulate max(|delta-z|) as a lower bound
434  std::vector<double> z0 = v.zvtx_vec;
435 
436  while (niter++ < maxIterations_) {
437  delta = update(beta, tks, v, rho0);
438  delta_sum_range += delta;
439 
440  if (delta_sum_range > zrange_min_) {
441  for (unsigned int k = 0; k < v.getSize(); k++) {
442  if (std::abs(v.zvtx_vec[k] - z0[k]) > zrange_min_) {
443  set_vtx_range(beta, tks, v);
444  delta_sum_range = 0;
445  z0 = v.zvtx_vec;
446  break;
447  }
448  }
449  }
450 
451  if (delta < delta_max) {
452  break;
453  }
454  }
455 
456 #ifdef DEBUG
457  if (DEBUGLEVEL > 0) {
458  std::cout << "DAClusterizerInZ_vect.thermalize niter = " << niter << " at T = " << 1 / beta
459  << " nv = " << v.getSize() << std::endl;
460  if (DEBUGLEVEL > 2)
461  dump(beta, v, tks, 0, rho0);
462  }
463 #endif
464 
465  return niter;
466 }
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
T sqrt(T t)
Definition: SSEVec.h:19
void dump(const double beta, const vertex_t &y, const track_t &tks, const int verbosity=0, const double rho0=0.) const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const

◆ update()

double DAClusterizerInZ_vect::update ( double  beta,
track_t gtracks,
vertex_t gvertices,
const double  rho0 = 0,
const bool  updateTc = false 
) const

Definition at line 284 of file DAClusterizerInZ_vect.cc.

References funct::abs(), HLT_2024v12_cff::beta, dumpMFGeometry_cfg::delta, DAClusterizerInZ_vect::vertex_t::exp, DAClusterizerInZ_vect::vertex_t::exp_arg, DAClusterizerInZ_vect::track_t::getSize(), DAClusterizerInZ_vect::vertex_t::getSize(), edm::isNotFinite(), dqmdumpme::k, DAClusterizerInZ_vect::track_t::kmax, DAClusterizerInZ_vect::track_t::kmin, SiStripPI::max, nt, DAClusterizerInZ_vect::track_t::osumtkwt, DAClusterizerInZ_vect::vertex_t::se, DAClusterizerInZ_vect::track_t::sum_Z, DAClusterizerInZ_vect::vertex_t::sw, DAClusterizerInZ_vect::vertex_t::swE, DAClusterizerInZ_vect::vertex_t::swz, DiMuonV_cfg::tracks, mitigatedMETSequence_cff::U, AlignmentTracksFromVertexSelector_cfi::vertices, and w().

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

285  {
286  // update weights and vertex positions
287  // returns the maximum of changes of vertex positions
288  // sums needed for Tc are only updated if updateTC == true
289 
290  const unsigned int nt = gtracks.getSize();
291  const unsigned int nv = gvertices.getSize();
292  auto osumtkwt = gtracks.osumtkwt;
293 
294  double Z_init = 0;
295  if (rho0 > 0) {
296  Z_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_); // cut-off
297  }
298 
299  // define kernels
300  auto kernel_calc_exp_arg_range = [beta](const unsigned int itrack,
301  track_t const& tracks,
302  vertex_t const& vertices,
303  const unsigned int kmin,
304  const unsigned int kmax) {
305  const double track_z = tracks.zpca[itrack];
306  const double botrack_dz2 = -beta * tracks.dz2[itrack];
307 
308  // auto-vectorized
309  for (unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
310  auto mult_res = track_z - vertices.zvtx[ivertex];
311  vertices.exp_arg[ivertex] = botrack_dz2 * (mult_res * mult_res);
312  }
313  };
314 
315  auto kernel_add_Z_range = [Z_init](
316  vertex_t const& vertices, const unsigned int kmin, const unsigned int kmax) -> double {
317  double ZTemp = Z_init;
318  for (unsigned int ivertex = kmin; ivertex < kmax; ++ivertex) {
319  ZTemp += vertices.rho[ivertex] * vertices.exp[ivertex];
320  }
321  return ZTemp;
322  };
323 
324  auto kernel_calc_normalization_range = [updateTc](const unsigned int track_num,
325  track_t& tracks,
326  vertex_t& vertices,
327  const unsigned int kmin,
328  const unsigned int kmax) {
329  auto o_trk_sum_Z = tracks.tkwt[track_num] / tracks.sum_Z[track_num];
330  auto o_trk_dz2 = tracks.dz2[track_num];
331  auto tmp_trk_z = tracks.zpca[track_num];
332 
333  // auto-vectorized
334  if (updateTc) {
335 #pragma GCC ivdep
336  for (unsigned int k = kmin; k < kmax; ++k) {
337  vertices.se[k] += vertices.exp[k] * o_trk_sum_Z;
338  auto w = vertices.rho[k] * vertices.exp[k] * (o_trk_sum_Z * o_trk_dz2);
339  vertices.sw[k] += w;
340  vertices.swz[k] += w * tmp_trk_z;
341  vertices.swE[k] += w * vertices.exp_arg[k];
342  }
343  } else {
344  // same loop but without updating sWE
345 #pragma GCC ivdep
346  for (unsigned int k = kmin; k < kmax; ++k) {
347  vertices.se[k] += vertices.exp[k] * o_trk_sum_Z;
348  auto w = vertices.rho[k] * vertices.exp[k] * (o_trk_sum_Z * o_trk_dz2);
349  vertices.sw[k] += w;
350  vertices.swz[k] += w * tmp_trk_z;
351  }
352  }
353  };
354 
355  if (updateTc) {
356  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
357  gvertices.se[ivertex] = 0.0;
358  gvertices.sw[ivertex] = 0.0;
359  gvertices.swz[ivertex] = 0.0;
360  gvertices.swE[ivertex] = 0.0;
361  }
362  } else {
363  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
364  gvertices.se[ivertex] = 0.0;
365  gvertices.sw[ivertex] = 0.0;
366  gvertices.swz[ivertex] = 0.0;
367  }
368  }
369 
370  // loop over tracks
371  for (auto itrack = 0U; itrack < nt; ++itrack) {
372  const unsigned int kmin = gtracks.kmin[itrack];
373  const unsigned int kmax = gtracks.kmax[itrack];
374 
375  kernel_calc_exp_arg_range(itrack, gtracks, gvertices, kmin, kmax);
376  local_exp_list_range(gvertices.exp_arg, gvertices.exp, kmin, kmax);
377  gtracks.sum_Z[itrack] = kernel_add_Z_range(gvertices, kmin, kmax);
378 
379  if (edm::isNotFinite(gtracks.sum_Z[itrack]))
380  gtracks.sum_Z[itrack] = 0.0;
381 
382  if (gtracks.sum_Z[itrack] > 1.e-100) {
383  kernel_calc_normalization_range(itrack, gtracks, gvertices, kmin, kmax);
384  }
385  }
386 
387  // (un-)apply the factor -beta which is needed in exp_arg, but not in swE
388  if (updateTc) {
389  auto obeta = -1. / beta;
390  for (auto ivertex = 0U; ivertex < nv; ++ivertex) {
391  gvertices.swE[ivertex] *= obeta;
392  }
393  }
394 
395  // now update z and rho
396  auto kernel_calc_z = [osumtkwt, nv](vertex_t& vertices) -> double {
397  double delta = 0;
398  // does not vectorize
399  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex) {
400  if (vertices.sw[ivertex] > 0) {
401  auto znew = vertices.swz[ivertex] / vertices.sw[ivertex];
402  // prevents from vectorizing if
403  delta = max(std::abs(vertices.zvtx[ivertex] - znew), delta);
404  vertices.zvtx[ivertex] = znew;
405  }
406  }
407 
408  for (unsigned int ivertex = 0; ivertex < nv; ++ivertex)
409  vertices.rho[ivertex] = vertices.rho[ivertex] * vertices.se[ivertex] * osumtkwt;
410 
411  return delta;
412  };
413 
414  double delta = kernel_calc_z(gvertices);
415 
416  // return how much the prototypes moved
417  return delta;
418 }
T w() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int nt
Definition: AMPTWrapper.h:42
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &tracks) const override

◆ verify()

void DAClusterizerInZ_vect::verify ( const vertex_t v,
const track_t tks,
unsigned int  nv = 999999,
unsigned int  nt = 999999 
) const

Definition at line 115 of file DAClusterizerInZ_vect.cc.

References cms::cuda::assert(), gather_cfg::cout, DAClusterizerInZ_vect::track_t::dz2, DAClusterizerInZ_vect::track_t::dz2_vec, DAClusterizerInZ_vect::track_t::getSize(), mps_fire::i, dqmdumpme::k, DAClusterizerInZ_vect::track_t::kmax, DAClusterizerInZ_vect::track_t::kmin, nt, DAClusterizerInZ_vect::track_t::sum_Z, DAClusterizerInZ_vect::track_t::sum_Z_vec, DAClusterizerInZ_vect::track_t::tkwt, DAClusterizerInZ_vect::track_t::tkwt_vec, DAClusterizerInZ_vect::track_t::tt, findQualityFiles::v, DAClusterizerInZ_vect::track_t::zpca, and DAClusterizerInZ_vect::track_t::zpca_vec.

115  {
116  if (!(nv == 999999)) {
117  assert(nv == v.getSize());
118  } else {
119  nv = v.getSize();
120  }
121 
122  if (!(nt == 999999)) {
123  assert(nt == tks.getSize());
124  } else {
125  nt = tks.getSize();
126  }
127 
128  assert(v.zvtx_vec.size() == nv);
129  assert(v.rho_vec.size() == nv);
130  assert(v.swz_vec.size() == nv);
131  assert(v.exp_arg_vec.size() == nv);
132  assert(v.exp_vec.size() == nv);
133  assert(v.se_vec.size() == nv);
134  assert(v.swz_vec.size() == nv);
135  assert(v.swE_vec.size() == nv);
136 
137  assert(v.zvtx == &v.zvtx_vec.front());
138  assert(v.rho == &v.rho_vec.front());
139  assert(v.exp_arg == &v.exp_arg_vec.front());
140  assert(v.sw == &v.sw_vec.front());
141  assert(v.swz == &v.swz_vec.front());
142  assert(v.se == &v.se_vec.front());
143  assert(v.swE == &v.swE_vec.front());
144 
145  for (unsigned int k = 0; k < nv - 1; k++) {
146  if (v.zvtx_vec[k] <= v.zvtx_vec[k + 1])
147  continue;
148  cout << " Z, cluster z-ordering assertion failure z[" << k << "] =" << v.zvtx_vec[k] << " z[" << k + 1
149  << "] =" << v.zvtx_vec[k + 1] << endl;
150  }
151 
152  assert(nt == tks.zpca_vec.size());
153  assert(nt == tks.dz2_vec.size());
154  assert(nt == tks.tt.size());
155  assert(nt == tks.tkwt_vec.size());
156  assert(nt == tks.sum_Z_vec.size());
157  assert(nt == tks.kmin.size());
158  assert(nt == tks.kmax.size());
159 
160  assert(tks.zpca == &tks.zpca_vec.front());
161  assert(tks.dz2 == &tks.dz2_vec.front());
162  assert(tks.tkwt == &tks.tkwt_vec.front());
163  assert(tks.sum_Z == &tks.sum_Z_vec.front());
164 
165  for (unsigned int i = 0; i < nt; i++) {
166  if ((tks.kmin[i] < tks.kmax[i]) && (tks.kmax[i] <= nv))
167  continue;
168  cout << "track vertex range assertion failure" << i << "/" << nt << " kmin,kmax=" << tks.kmin[i] << ", "
169  << tks.kmax[i] << " nv=" << nv << endl;
170  }
171 
172  for (unsigned int i = 0; i < nt; i++) {
173  assert((tks.kmin[i] < tks.kmax[i]) && (tks.kmax[i] <= nv));
174  }
175 }
assert(be >=bs)
int nt
Definition: AMPTWrapper.h:42

◆ vertices()

vector< TransientVertex > DAClusterizerInZ_vect::vertices ( const std::vector< reco::TransientTrack > &  tracks) const
overridevirtual

Implements TrackClusterizerInZ.

Definition at line 1331 of file DAClusterizerInZ_vect.cc.

References DiMuonV_cfg::tracks.

1331  {
1332  if (runInBlocks_ and (block_size_ < tracks.size())) //doesn't bother if low number of tracks
1333  return vertices_in_blocks(tracks);
1334  else
1335  return vertices_no_blocks(tracks);
1336 }
std::vector< TransientVertex > vertices_no_blocks(const std::vector< reco::TransientTrack > &tracks) const
std::vector< TransientVertex > vertices_in_blocks(const std::vector< reco::TransientTrack > &tracks) const

◆ vertices_in_blocks()

vector< TransientVertex > DAClusterizerInZ_vect::vertices_in_blocks ( const std::vector< reco::TransientTrack > &  tracks) const

Definition at line 928 of file DAClusterizerInZ_vect.cc.

References a, funct::abs(), b, HLT_2024v12_cff::beta, groupFilesInBlocks::block, bsc_activity_cfg::clusters, gather_cfg::cout, GCP_Ntuples_cfg::dump, MillePedeFileConverter_cfg::e, relativeConstraints::empty, DAClusterizerInZ_vect::track_t::extractRaw(), ntuplemaker::fill, dqmdumpme::first, mps_fire::i, createfilelist::int, dqmdumpme::k, SiStripPI::max, merge(), SiStripPI::min, nt, AlCaHLTBitMon_ParallelJobs::p, position, jetUpdater_cfi::sort, submitPVValidationJobs::split(), mathSSE::sqrt(), DiMuonV_cfg::threshold, DiMuonV_cfg::tracks, mitigatedMETSequence_cff::U, update, findQualityFiles::v, SiStripMonitorCluster_cfi::zmax, SiStripMonitorCluster_cfi::zmin, and OfflinePixel3DPrimaryVertices_cfi::zrange.

928  {
929  vector<reco::TransientTrack> sorted_tracks;
930  vector<pair<float, float>> vertices_tot; // z, rho for each vertex
931  for (unsigned int i = 0; i < tracks.size(); i++) {
932  sorted_tracks.push_back(tracks[i]);
933  }
934  double rho0, beta;
935  std::sort(sorted_tracks.begin(),
936  sorted_tracks.end(),
937  [](const reco::TransientTrack& a, const reco::TransientTrack& b) -> bool {
938  return (a.stateAtBeamLine().trackStateAtPCA()).position().z() <
939  (b.stateAtBeamLine().trackStateAtPCA()).position().z();
940  });
941 
942  unsigned int nBlocks = (unsigned int)std::floor(sorted_tracks.size() / (block_size_ * (1 - overlap_frac_)));
943  if (nBlocks < 1) {
944  nBlocks = 1;
945  edm::LogWarning("DAClusterizerinZ_vect")
946  << "Warning nBlocks was 0 with ntracks = " << sorted_tracks.size() << " block_size = " << block_size_
947  << " and overlap fraction = " << overlap_frac_ << ". Setting nBlocks = 1";
948  }
949  for (unsigned int block = 0; block < nBlocks; block++) {
950  vector<reco::TransientTrack> block_tracks;
951  unsigned int begin = (unsigned int)(block * block_size_ * (1 - overlap_frac_));
952  unsigned int end = (unsigned int)std::min(begin + block_size_, (unsigned int)sorted_tracks.size());
953  for (unsigned int i = begin; i < end; i++) {
954  block_tracks.push_back(sorted_tracks[i]);
955  }
956  if (block_tracks.empty()) {
957  continue;
958  }
959 
960 #ifdef DEBUG
961  std::cout << "Running vertices_in_blocks on" << std::endl;
962  std::cout << "- block no." << block << " on " << nBlocks << " blocks " << std::endl;
963  std::cout << "- block track size: " << sorted_tracks.size() << " - block size: " << block_size_ << std::endl;
964 #endif
965  track_t&& tks = fill(block_tracks);
966  tks.extractRaw();
967 
968  rho0 = 0.0; // start with no outlier rejection
969 
970  vertex_t y; // the vertex prototypes
971 
972  // initialize:single vertex at infinite temperature
973  y.addItem(0, 1.0);
974  clear_vtx_range(tks, y);
975 
976  // estimate first critical temperature
977  beta = beta0(betamax_, tks, y);
978 #ifdef DEBUG
979  if (DEBUGLEVEL > 0)
980  std::cout << "Beta0 is " << beta << std::endl;
981 #endif
982 
983  thermalize(beta, tks, y, delta_highT_);
984 
985  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
986 
987  double betafreeze = betamax_ * sqrt(coolingFactor_);
988  while (beta < betafreeze) {
989  while (merge(y, tks, beta)) {
990  update(beta, tks, y, rho0, false);
991  }
992  split(beta, tks, y);
993 
995  thermalize(beta, tks, y, delta_highT_);
996  }
997 
998 #ifdef DEBUG
999  verify(y, tks);
1000 
1001  if (DEBUGLEVEL > 0) {
1002  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1003  << "last round of splitting" << std::endl;
1004  }
1005 #endif
1006 
1007  set_vtx_range(beta, tks, y);
1008  update(beta, tks, y, rho0, false);
1009 
1010  while (merge(y, tks, beta)) {
1011  set_vtx_range(beta, tks, y);
1012  update(beta, tks, y, rho0, false);
1013  }
1014 
1015  unsigned int ntry = 0;
1016  double threshold = 1.0;
1017  while (split(beta, tks, y, threshold) && (ntry++ < 10)) {
1018  thermalize(beta, tks, y, delta_highT_, rho0); // rho0 = 0. here
1019  while (merge(y, tks, beta)) {
1020  update(beta, tks, y, rho0, false);
1021  }
1022 
1023  // relax splitting a bit to reduce multiple split-merge cycles of the same cluster
1024  threshold *= 1.1;
1025  }
1026 
1027 #ifdef DEBUG
1028  verify(y, tks);
1029  if (DEBUGLEVEL > 0) {
1030  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1031  << "turning on outlier rejection at T=" << 1 / beta << std::endl;
1032  }
1033 #endif
1034 
1035  // switch on outlier rejection at T=Tmin
1036  if (dzCutOff_ > 0) {
1037  rho0 = y.getSize() > 1 ? 1. / y.getSize() : 1.;
1038  for (unsigned int a = 0; a < 5; a++) {
1039  update(beta, tks, y, a * rho0 / 5.); // adiabatic turn-on
1040  }
1041  }
1042 
1043  thermalize(beta, tks, y, delta_lowT_, rho0);
1044 
1045 #ifdef DEBUG
1046  verify(y, tks);
1047  if (DEBUGLEVEL > 0) {
1048  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1049  << "merging with outlier rejection at T=" << 1 / beta << std::endl;
1050  }
1051  if (DEBUGLEVEL > 2)
1052  dump(beta, y, tks, 2, rho0);
1053 #endif
1054 
1055  // merge again (some cluster split by outliers collapse here)
1056  while (merge(y, tks, beta)) {
1057  set_vtx_range(beta, tks, y);
1058  update(beta, tks, y, rho0, false);
1059  }
1060 
1061 #ifdef DEBUG
1062  verify(y, tks);
1063  if (DEBUGLEVEL > 0) {
1064  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1065  << "after merging with outlier rejection at T=" << 1 / beta << std::endl;
1066  }
1067  if (DEBUGLEVEL > 2)
1068  dump(beta, y, tks, 2, rho0);
1069 #endif
1070 
1071  // go down to the purging temperature (if it is lower than tmin)
1072  while (beta < betapurge_) {
1074  thermalize(beta, tks, y, delta_lowT_, rho0);
1075  }
1076 
1077 #ifdef DEBUG
1078  verify(y, tks);
1079  if (DEBUGLEVEL > 0) {
1080  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1081  << "purging at T=" << 1 / beta << std::endl;
1082  }
1083 #endif
1084 
1085  // eliminate insignificant vertices, this is more restrictive at higher T
1086  while (purge(y, tks, rho0, beta)) {
1087  thermalize(beta, tks, y, delta_lowT_, rho0);
1088  }
1089 
1090 #ifdef DEBUG
1091  verify(y, tks);
1092  if (DEBUGLEVEL > 0) {
1093  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1094  << "last cooling T=" << 1 / beta << std::endl;
1095  }
1096 #endif
1097 
1098  // optionally cool some more without doing anything, to make the assignment harder
1099  while (beta < betastop_) {
1101  thermalize(beta, tks, y, delta_lowT_, rho0);
1102  }
1103 
1104 #ifdef DEBUG
1105  verify(y, tks);
1106  if (DEBUGLEVEL > 0) {
1107  std::cout << "DAClusterizerInZSubCluster_vect::vertices :"
1108  << "stop cooling at T=" << 1 / beta << std::endl;
1109  }
1110  if (DEBUGLEVEL > 2)
1111  dump(beta, y, tks, 2, rho0);
1112 #endif
1113 
1114  for (unsigned int ivertex = 0; ivertex < y.getSize(); ivertex++) {
1115  if (y.zvtx_vec[ivertex] != 0 && y.rho_vec[ivertex] != 0) {
1116  vertices_tot.push_back(pair(y.zvtx_vec[ivertex], y.rho_vec[ivertex]));
1117 #ifdef DEBUG
1118  std::cout << "Found new vertex " << y.zvtx_vec[ivertex] << " , " << y.rho_vec[ivertex] << std::endl;
1119 #endif
1120  }
1121  }
1122  }
1123 
1124  std::sort(vertices_tot.begin(),
1125  vertices_tot.end(),
1126  [](const pair<float, float>& a, const pair<float, float>& b) -> bool { return a.first < b.first; });
1127 
1128  // reassign tracks to vertices
1129  track_t&& tracks_tot = fill(tracks);
1130  const unsigned int nv = vertices_tot.size();
1131  const unsigned int nt = tracks_tot.getSize();
1132 
1133  for (auto itrack = 0U; itrack < nt; ++itrack) {
1134  double zrange = max(sel_zrange_ / sqrt(beta * tracks_tot.dz2[itrack]), zrange_min_);
1135 
1136  double zmin = tracks_tot.zpca[itrack] - zrange;
1137  unsigned int kmin = min(nv - 1, tracks_tot.kmin[itrack]);
1138  // find the smallest vertex_z that is larger than zmin
1139  if (vertices_tot[kmin].first > zmin) {
1140  while ((kmin > 0) && (vertices_tot[kmin - 1].first > zmin)) {
1141  kmin--;
1142  }
1143  } else {
1144  while ((kmin < (nv - 1)) && (vertices_tot[kmin].first < zmin)) {
1145  kmin++;
1146  }
1147  }
1148 
1149  double zmax = tracks_tot.zpca[itrack] + zrange;
1150  unsigned int kmax = min(nv - 1, tracks_tot.kmax[itrack] - 1);
1151  // note: kmax points to the last vertex in the range, while gtracks.kmax points to the entry BEHIND the last vertex
1152  // find the largest vertex_z that is smaller than zmax
1153  if (vertices_tot[kmax].first < zmax) {
1154  while ((kmax < (nv - 1)) && (vertices_tot[kmax + 1].first < zmax)) {
1155  kmax++;
1156  }
1157  } else {
1158  while ((kmax > 0) && (vertices_tot[kmax].first > zmax)) {
1159  kmax--;
1160  }
1161  }
1162 
1163  if (kmin <= kmax) {
1164  tracks_tot.kmin[itrack] = kmin;
1165  tracks_tot.kmax[itrack] = kmax + 1;
1166  } else {
1167  tracks_tot.kmin[itrack] = max(0U, min(kmin, kmax));
1168  tracks_tot.kmax[itrack] = min(nv, max(kmin, kmax) + 1);
1169  }
1170  }
1171 
1172  rho0 = nv > 1 ? 1. / nv : 1.;
1173  const auto z_sum_init = rho0 * local_exp(-beta * dzCutOff_ * dzCutOff_);
1174 
1175  std::vector<std::vector<unsigned int>> vtx_track_indices(nv);
1176  for (unsigned int i = 0; i < nt; i++) {
1177  const auto kmin = tracks_tot.kmin[i];
1178  const auto kmax = tracks_tot.kmax[i];
1179  double p_max = -1;
1180  unsigned int iMax = 10000; //just a "big" number w.r.t. number of vertices
1181  float sum_Z = z_sum_init;
1182  for (auto k = kmin; k < kmax; k++) {
1183  float v_exp = local_exp(-beta * Eik(tracks_tot.zpca[i], vertices_tot[k].first, tracks_tot.dz2[i]));
1184  sum_Z += vertices_tot[k].second * v_exp;
1185  }
1186  double invZ = sum_Z > 1e-100 ? 1. / sum_Z : 0.0;
1187  for (auto k = kmin; k < kmax && invZ != 0.0; k++) {
1188  float v_exp = local_exp(-beta * Eik(tracks_tot.zpca[i], vertices_tot[k].first, tracks_tot.dz2[i]));
1189  double p = vertices_tot[k].second * v_exp * invZ;
1190  if (p > p_max && p > mintrkweight_) {
1191  p_max = p;
1192  iMax = k;
1193  }
1194  }
1195  if (iMax < vtx_track_indices.size()) {
1196  vtx_track_indices[iMax].push_back(i);
1197  }
1198  }
1199 #ifdef DEBUG
1200  for (auto itrack = 0U; itrack < nt; ++itrack) {
1201  std::cout << "itrack " << itrack << " , " << tracks_tot.kmin[itrack] << " , " << tracks_tot.kmax[itrack]
1202  << std::endl;
1203  }
1204 #endif
1205 
1206  vector<TransientVertex> clusters;
1207  if (nv == 0) {
1208  return clusters;
1209  }
1210 
1211  GlobalError dummyError(0.01, 0, 0.01, 0., 0., 0.01);
1212  vector<reco::TransientTrack> vertexTracks;
1213 
1214  for (unsigned int k = 0; k < nv; k++) {
1215  if (!vtx_track_indices[k].empty()) {
1216  for (auto i : vtx_track_indices[k]) {
1217  vertexTracks.push_back(*(tracks_tot.tt[i]));
1218 #ifdef DEBUG
1219  std::cout << vertices_tot[k].first << ","
1220  << (*(tracks_tot.tt[i])).stateAtBeamLine().trackStateAtPCA().position().z() << std::endl;
1221 #endif
1222  }
1223  }
1224 
1225  // implement what clusterize() did before : merge left-to-right if distance < 2 * vertexSize_
1226  if ((k + 1 == nv) || (abs(vertices_tot[k + 1].first - vertices_tot[k].first) > (2 * vertexSize_))) {
1227  // close a cluster
1228  if (vertexTracks.size() > 1) {
1229  GlobalPoint pos(0, 0, vertices_tot[k].first); // only usable with subsequent fit
1230  TransientVertex v(pos, dummyError, vertexTracks, 0);
1231  clusters.push_back(v);
1232  }
1233  vertexTracks.clear();
1234  }
1235  }
1236 
1237  return clusters;
1238 } // end of vertices_in_blocks
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
bool merge(vertex_t &y, track_t &tks, double &beta) const
void clear_vtx_range(track_t &gtracks, vertex_t &gvertices) const
bool purge(vertex_t &, track_t &, double &, const double) const
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
T sqrt(T t)
Definition: SSEVec.h:19
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int nt
Definition: AMPTWrapper.h:42
unsigned int thermalize(double beta, track_t &gtracks, vertex_t &gvertices, const double delta_max, const double rho0=0.) const
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
static int position[264][3]
Definition: ReadPGInfo.cc:289
Log< level::Warning, false > LogWarning
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

◆ vertices_no_blocks()

vector< TransientVertex > DAClusterizerInZ_vect::vertices_no_blocks ( const std::vector< reco::TransientTrack > &  tracks) const

Definition at line 769 of file DAClusterizerInZ_vect.cc.

References a, HLT_2024v12_cff::beta, bsc_activity_cfg::clusters, gather_cfg::cout, GCP_Ntuples_cfg::dump, DAClusterizerInZ_vect::track_t::extractRaw(), ntuplemaker::fill, DAClusterizerInZ_vect::track_t::getSize(), merge(), SiStripPI::min, submitPVValidationJobs::split(), mathSSE::sqrt(), DiMuonV_cfg::threshold, DiMuonV_cfg::tracks, and update.

769  {
770  track_t&& tks = fill(tracks);
771  tks.extractRaw();
772 
773  double rho0 = 0.0; // start with no outlier rejection
774 
775  vector<TransientVertex> clusters;
776  if (tks.getSize() == 0)
777  return clusters;
778 
779  vertex_t y; // the vertex prototypes
780 
781  // initialize:single vertex at infinite temperature
782  y.addItem(0, 1.0);
783  clear_vtx_range(tks, y);
784 
785  // estimate first critical temperature
786  double beta = beta0(betamax_, tks, y);
787 #ifdef DEBUG
788  if (DEBUGLEVEL > 0)
789  std::cout << "Beta0 is " << beta << std::endl;
790 #endif
791 
792  thermalize(beta, tks, y, delta_highT_);
793 
794  // annealing loop, stop when T<Tmin (i.e. beta>1/Tmin)
795 
796  double betafreeze = betamax_ * sqrt(coolingFactor_);
797 
798  while (beta < betafreeze) {
799  while (merge(y, tks, beta)) {
800  update(beta, tks, y, rho0, false);
801  }
802  split(beta, tks, y);
803 
805  thermalize(beta, tks, y, delta_highT_);
806  }
807 
808 #ifdef DEBUG
809  verify(y, tks);
810 
811  if (DEBUGLEVEL > 0) {
812  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
813  << "last round of splitting" << std::endl;
814  }
815 #endif
816 
817  set_vtx_range(beta, tks, y);
818  update(beta, tks, y, rho0, false);
819 
820  while (merge(y, tks, beta)) {
821  set_vtx_range(beta, tks, y);
822  update(beta, tks, y, rho0, false);
823  }
824 
825  unsigned int ntry = 0;
826  double threshold = 1.0;
827  while (split(beta, tks, y, threshold) && (ntry++ < 10)) {
828  thermalize(beta, tks, y, delta_highT_, rho0); // rho0 = 0. here
829  while (merge(y, tks, beta)) {
830  update(beta, tks, y, rho0, false);
831  }
832 
833  // relax splitting a bit to reduce multiple split-merge cycles of the same cluster
834  threshold *= 1.1;
835  }
836 
837 #ifdef DEBUG
838  verify(y, tks);
839  if (DEBUGLEVEL > 0) {
840  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
841  << "turning on outlier rejection at T=" << 1 / beta << std::endl;
842  }
843 #endif
844 
845  // switch on outlier rejection at T=Tmin
846  if (dzCutOff_ > 0) {
847  rho0 = y.getSize() > 1 ? 1. / y.getSize() : 1.;
848  for (unsigned int a = 0; a < 5; a++) {
849  update(beta, tks, y, a * rho0 / 5.); // adiabatic turn-on
850  }
851  }
852 
853  thermalize(beta, tks, y, delta_lowT_, rho0);
854 
855 #ifdef DEBUG
856  verify(y, tks);
857  if (DEBUGLEVEL > 0) {
858  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
859  << "merging with outlier rejection at T=" << 1 / beta << std::endl;
860  }
861  if (DEBUGLEVEL > 2)
862  dump(beta, y, tks, 2, rho0);
863 #endif
864 
865  // merge again (some cluster split by outliers collapse here)
866  while (merge(y, tks, beta)) {
867  set_vtx_range(beta, tks, y);
868  update(beta, tks, y, rho0, false);
869  }
870 
871 #ifdef DEBUG
872  verify(y, tks);
873  if (DEBUGLEVEL > 0) {
874  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
875  << "after merging with outlier rejection at T=" << 1 / beta << std::endl;
876  }
877  if (DEBUGLEVEL > 2)
878  dump(beta, y, tks, 2, rho0);
879 #endif
880 
881  // go down to the purging temperature (if it is lower than tmin)
882  while (beta < betapurge_) {
884  thermalize(beta, tks, y, delta_lowT_, rho0);
885  }
886 
887 #ifdef DEBUG
888  verify(y, tks);
889  if (DEBUGLEVEL > 0) {
890  std::cout << "DAClusterizerInZ_vect::vertices :"
891  << "purging at T=" << 1 / beta << std::endl;
892  }
893 #endif
894 
895  // eliminate insignificant vertices, this is more restrictive at higher T
896  while (purge(y, tks, rho0, beta)) {
897  thermalize(beta, tks, y, delta_lowT_, rho0);
898  }
899 
900 #ifdef DEBUG
901  verify(y, tks);
902  if (DEBUGLEVEL > 0) {
903  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
904  << "last cooling T=" << 1 / beta << std::endl;
905  }
906 #endif
907 
908  // optionally cool some more without doing anything, to make the assignment harder
909  while (beta < betastop_) {
911  thermalize(beta, tks, y, delta_lowT_, rho0);
912  }
913 
914 #ifdef DEBUG
915  verify(y, tks);
916  if (DEBUGLEVEL > 0) {
917  std::cout << "DAClusterizerInZ_vect::vertices_no_blocks :"
918  << "stop cooling at T=" << 1 / beta << std::endl;
919  }
920  if (DEBUGLEVEL > 2)
921  dump(beta, y, tks, 2, rho0);
922 #endif
923 
924  // assign tracks and fill into transient vertices
925  return fill_vertices(beta, rho0, tks, y);
926 }
double beta0(const double betamax, track_t const &tks, vertex_t const &y) const
bool merge(vertex_t &y, track_t &tks, double &beta) const
void clear_vtx_range(track_t &gtracks, vertex_t &gvertices) const
bool purge(vertex_t &, track_t &, double &, const double) const
std::vector< TransientVertex > fill_vertices(double beta, double rho0, track_t &tracks, vertex_t &vertices) const
double update(double beta, track_t &gtracks, vertex_t &gvertices, const double rho0=0, const bool updateTc=false) const
T sqrt(T t)
Definition: SSEVec.h:19
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
unsigned int thermalize(double beta, track_t &gtracks, vertex_t &gvertices, const double delta_max, const double rho0=0.) const
void set_vtx_range(double beta, track_t &gtracks, vertex_t &gvertices) const
double a
Definition: hdecay.h:121
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

Member Data Documentation

◆ betamax_

double DAClusterizerInZ_vect::betamax_
private

Definition at line 211 of file DAClusterizerInZ_vect.h.

◆ betapurge_

double DAClusterizerInZ_vect::betapurge_
private

Definition at line 220 of file DAClusterizerInZ_vect.h.

◆ betastop_

double DAClusterizerInZ_vect::betastop_
private

Definition at line 212 of file DAClusterizerInZ_vect.h.

◆ block_size_

unsigned int DAClusterizerInZ_vect::block_size_
private

Definition at line 230 of file DAClusterizerInZ_vect.h.

◆ convergence_mode_

unsigned int DAClusterizerInZ_vect::convergence_mode_
private

Definition at line 222 of file DAClusterizerInZ_vect.h.

◆ coolingFactor_

double DAClusterizerInZ_vect::coolingFactor_
private

Definition at line 210 of file DAClusterizerInZ_vect.h.

◆ d0CutOff_

double DAClusterizerInZ_vect::d0CutOff_
private

Definition at line 214 of file DAClusterizerInZ_vect.h.

◆ delta_highT_

double DAClusterizerInZ_vect::delta_highT_
private

Definition at line 223 of file DAClusterizerInZ_vect.h.

◆ delta_lowT_

double DAClusterizerInZ_vect::delta_lowT_
private

Definition at line 224 of file DAClusterizerInZ_vect.h.

◆ dzCutOff_

double DAClusterizerInZ_vect::dzCutOff_
private

Definition at line 213 of file DAClusterizerInZ_vect.h.

◆ maxIterations_

unsigned int DAClusterizerInZ_vect::maxIterations_
private

Definition at line 209 of file DAClusterizerInZ_vect.h.

◆ mintrkweight_

double DAClusterizerInZ_vect::mintrkweight_
private

Definition at line 216 of file DAClusterizerInZ_vect.h.

◆ overlap_frac_

double DAClusterizerInZ_vect::overlap_frac_
private

Definition at line 231 of file DAClusterizerInZ_vect.h.

◆ runInBlocks_

bool DAClusterizerInZ_vect::runInBlocks_
private

Definition at line 229 of file DAClusterizerInZ_vect.h.

◆ sel_zrange_

double DAClusterizerInZ_vect::sel_zrange_
private

Definition at line 226 of file DAClusterizerInZ_vect.h.

◆ uniquetrkminp_

double DAClusterizerInZ_vect::uniquetrkminp_
private

Definition at line 218 of file DAClusterizerInZ_vect.h.

◆ uniquetrkweight_

double DAClusterizerInZ_vect::uniquetrkweight_
private

Definition at line 217 of file DAClusterizerInZ_vect.h.

◆ vertexSize_

double DAClusterizerInZ_vect::vertexSize_
private

Definition at line 208 of file DAClusterizerInZ_vect.h.

◆ zdumpcenter_

double DAClusterizerInZ_vect::zdumpcenter_
private

Definition at line 205 of file DAClusterizerInZ_vect.h.

◆ zdumpwidth_

double DAClusterizerInZ_vect::zdumpwidth_
private

Definition at line 206 of file DAClusterizerInZ_vect.h.

◆ zmerge_

double DAClusterizerInZ_vect::zmerge_
private

Definition at line 219 of file DAClusterizerInZ_vect.h.

◆ zrange_min_

const double DAClusterizerInZ_vect::zrange_min_ = 0.1
private

Definition at line 227 of file DAClusterizerInZ_vect.h.