CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes
AdaptiveChisquarePrimaryVertexFitter Class Reference

#include <AdaptiveChisquarePrimaryVertexFitter.h>

Inheritance diagram for AdaptiveChisquarePrimaryVertexFitter:
PrimaryVertexFitterBase

Classes

struct  TrackInfo
 

Public Types

using Error3 = ROOT::Math::SMatrix< double, 3 >
 

Public Member Functions

 AdaptiveChisquarePrimaryVertexFitter (double chicutoff=2.5, double zcutoff=1.0, double mintrkweight=0.4, bool multivertexfit=false)
 
std::vector< TransientVertexfit (const std::vector< reco::TransientTrack > &, const std::vector< TransientVertex > &, const reco::BeamSpot &, const bool) override
 
 ~AdaptiveChisquarePrimaryVertexFitter () override=default
 
- Public Member Functions inherited from PrimaryVertexFitterBase
 PrimaryVertexFitterBase (const edm::ParameterSet &conf)
 
 PrimaryVertexFitterBase ()
 
virtual ~PrimaryVertexFitterBase ()=default
 

Protected Member Functions

bool clean ()
 
void fill_trackinfo (const std::vector< reco::TransientTrack > &, const reco::BeamSpot &)
 
void fill_weights (const reco::BeamSpot &, const double beta=1.)
 
Error3 get_inverse_beam_covariance (const reco::BeamSpot &)
 
TransientVertex get_TransientVertex (const unsigned int, const std::vector< std::pair< unsigned int, float >> &, const std::vector< reco::TransientTrack > &, const float, const reco::BeamSpot &)
 
void make_vtx_trk_map (const double)
 
TransientVertex refit (const TransientVertex &, const reco::BeamSpot &, const bool)
 
void remove_vertex (unsigned int)
 
double track_in_vertex_chsq (const TrackInfo &, const double, const double, const double)
 
double update (const reco::BeamSpot &, float beam_weight, const bool fill_covariances=false)
 
void verify ()
 
std::vector< TransientVertexvertices (const std::vector< reco::TransientTrack > &, const std::vector< TransientVertex > &, const reco::BeamSpot &, const bool)
 

Protected Attributes

double chi_cutoff_
 
std::vector< Error3covv_
 
double min_trackweight_
 
double multivertexfit_
 
std::vector< unsigned int > tkfirstv_
 
std::vector< unsigned int > tkmap_
 
std::vector< double > tkweight_
 
std::vector< TrackInfotrackinfo_
 
std::vector< double > xv_
 
std::vector< double > yv_
 
double z_cutoff_
 
std::vector< double > zv_
 

Detailed Description

Description: simultaneous chisquared fit of primary vertices

Definition at line 15 of file AdaptiveChisquarePrimaryVertexFitter.h.

Member Typedef Documentation

◆ Error3

using AdaptiveChisquarePrimaryVertexFitter::Error3 = ROOT::Math::SMatrix<double, 3>

Definition at line 28 of file AdaptiveChisquarePrimaryVertexFitter.h.

Constructor & Destructor Documentation

◆ AdaptiveChisquarePrimaryVertexFitter()

AdaptiveChisquarePrimaryVertexFitter::AdaptiveChisquarePrimaryVertexFitter ( double  chicutoff = 2.5,
double  zcutoff = 1.0,
double  mintrkweight = 0.4,
bool  multivertexfit = false 
)

Definition at line 8 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References chi_cutoff_, min_trackweight_, multivertexfit_, and z_cutoff_.

12  : chi_cutoff_(chicutoff), z_cutoff_(zcutoff), min_trackweight_(mintrkweight), multivertexfit_(multivertexfit) {
13 #ifdef PVTX_DEBUG
14  edm::LogInfo("AdaptiveChisquarePrimaryVertexFitter")
15  << "instantiating AdaptiveChisquarePrimaryVertexFitter with "
16  << " chi cutoff =" << chi_cutoff_ << " z cutoff =" << z_cutoff_ << " min_trackweight=" << min_trackweight_
17  << " multivertex mode " << multivertexfit_ << std::endl;
18 #endif
19 }
Log< level::Info, false > LogInfo

◆ ~AdaptiveChisquarePrimaryVertexFitter()

AdaptiveChisquarePrimaryVertexFitter::~AdaptiveChisquarePrimaryVertexFitter ( )
overridedefault

Member Function Documentation

◆ clean()

bool AdaptiveChisquarePrimaryVertexFitter::clean ( )
protected

Definition at line 177 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References funct::abs(), PVValHelper::dz, dqmiolumiharvest::j, dqmdumpme::k, SiStripPI::min, remove_vertex(), tkfirstv_, tkweight_, xv_, and zv_.

Referenced by vertices().

177  {
178  /* in multi-vertex fitting, nearby vertices can fall on top of each other,
179  even when the initial seeds don't, some kind of duplicate removal is required
180  the approach in this method is similar to the method applied in clustering:
181  at least two tracks with a weight above a threshold (trkweight_threshold) are required.
182  vertices that don't fulfill this are either insignficant or very close
183  to another vertex
184  */
185  const double trkweight_threshold = 0.7;
186  unsigned int nv = xv_.size();
187  if (nv < 2)
188  return false;
189 
190  // sum of weights per vertex
191  std::vector<double> wsumhi(nv, 0);
192  for (unsigned int k = 0; k < nv; k++) {
193  for (unsigned int j = tkfirstv_[k]; j < tkfirstv_[k + 1]; j++) {
194  if (tkweight_[j] > trkweight_threshold)
195  wsumhi[k] += tkweight_[j];
196  }
197  }
198 
199  double dzmin = 0;
200  unsigned int k_dzmin = 0;
201  for (unsigned int k = 0; k < nv - 1; k++) {
202  double dz = std::abs(zv_[k + 1] - zv_[k]);
203  if ((k == 0) || (dz < dzmin)) {
204  dzmin = dz;
205  k_dzmin = k;
206  }
207  }
208 
209  if ((std::abs(dzmin) < 0.0200) && (std::min(wsumhi[k_dzmin], wsumhi[k_dzmin + 1]) < 0.5)) {
210  if (wsumhi[k_dzmin] < wsumhi[k_dzmin + 1]) {
211  remove_vertex(k_dzmin);
212  } else {
213  remove_vertex(k_dzmin + 1);
214  }
215  }
216 
217  return true;
218 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22

◆ fill_trackinfo()

void AdaptiveChisquarePrimaryVertexFitter::fill_trackinfo ( const std::vector< reco::TransientTrack > &  tracks,
const reco::BeamSpot beamSpot 
)
protected

Definition at line 34 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References AdaptiveChisquarePrimaryVertexFitter::TrackInfo::b1, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::b2, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::C, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::dzError, PerigeeConversions::ftsToPerigeeError(), AdaptiveChisquarePrimaryVertexFitter::TrackInfo::g, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::H1, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::H2, dqmdumpme::k, MainPageGenerator::l, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::S11, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::S12, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::S22, trackinfo_, DiMuonV_cfg::tracks, and AdaptiveChisquarePrimaryVertexFitter::TrackInfo::zpca.

Referenced by refit(), and vertices().

35  {
36  /* fill track information used during fits into arrays, parallell to the list of input tracks */
37 
38  trackinfo_.clear();
39  trackinfo_.reserve(tracks.size());
40 
41  for (auto &trk : tracks) {
42  TrackInfo ti;
43  // F1,F2 are the perigee parameters (3,4)
44  const auto tspca = trk.stateAtBeamLine().trackStateAtPCA(); // freeTrajectoryState
45  const auto tspca_pe = PerigeeConversions::ftsToPerigeeError(tspca);
46  const auto momentum = tspca.momentum();
47  auto const cos_phi = momentum.x() / momentum.perp();
48  auto const sin_phi = momentum.y() / momentum.perp();
49  auto const tan_lambda = momentum.z() / momentum.perp();
50 
51  // covariance matrix of (F1,F2)
52  double cov11 = tspca_pe.covarianceMatrix()(3, 3);
53  double cov22 = tspca_pe.covarianceMatrix()(4, 4);
54  double cov12 = tspca_pe.covarianceMatrix()(3, 4);
55 
56  // S = cov^{-1}
57  double DetV = cov11 * cov22 - cov12 * cov12;
58  if (fabs(DetV) < 1.e-16) {
59  edm::LogWarning("AdaptiveChisquarePrimaryVertexFitter")
60  << "Warning, det(V) almost vanishes : " << DetV << " !! This should not happen!" << std::endl;
61  ti.S11 = 0;
62  ti.S22 = 0;
63  ti.S12 = 0;
64  } else {
65  ti.S11 = cov22 / DetV;
66  ti.S22 = cov11 / DetV;
67  ti.S12 = -cov12 / DetV;
68  }
69  ti.b1 = tspca.position().x() * sin_phi - tspca.position().y() * cos_phi;
70  ti.H1[0] = -sin_phi;
71  ti.H1[1] = cos_phi;
72  ti.H1[2] = 0;
73  ti.b2 = tspca.position().z() - (tspca.position().x() * cos_phi + tspca.position().y() * sin_phi) * tan_lambda;
74  ti.H2[0] = cos_phi * tan_lambda;
75  ti.H2[1] = sin_phi * tan_lambda;
76  ti.H2[2] = -1.;
77 
78  for (int k = 0; k < 3; k++) {
79  double SH1k = (ti.S11 * ti.H1[k] + ti.S12 * ti.H2[k]);
80  double SH2k = (ti.S12 * ti.H1[k] + ti.S22 * ti.H2[k]);
81  ti.g[k] = ti.b1 * SH1k + ti.b2 * SH2k;
82  for (int l = 0; l < 3; l++) {
83  ti.C(l, k) = ti.H1[l] * SH1k + ti.H2[l] * SH2k;
84  }
85  }
86 
87  ti.zpca = tspca.position().z();
88  ti.dzError = trk.track().dzError();
89  trackinfo_.push_back(ti);
90  }
91 }
PerigeeTrajectoryError ftsToPerigeeError(const FTS &originalFTS)
Log< level::Warning, false > LogWarning

◆ fill_weights()

void AdaptiveChisquarePrimaryVertexFitter::fill_weights ( const reco::BeamSpot beamspot,
const double  beta = 1. 
)
protected

Definition at line 141 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References cms::cuda::assert(), HLT_2024v11_cff::beta, chi_cutoff_, MillePedeFileConverter_cfg::e, mps_fire::i, dqmiolumiharvest::j, dqmdumpme::k, nt, tkfirstv_, tkmap_, tkweight_, track_in_vertex_chsq(), trackinfo_, xv_, yv_, and zv_.

Referenced by refit(), and vertices().

141  {
142  // multi-vertex version
143  unsigned const int nt = trackinfo_.size();
144  unsigned const int nv = xv_.size();
145  const double beta_over_2 = 0.5 * beta;
146  const double argmax = beta_over_2 * chi_cutoff_ * chi_cutoff_ * 5;
147  const double Z_cutoff = vdt::fast_exp(-beta_over_2 * chi_cutoff_ * chi_cutoff_);
148 
149  std::vector<double> Z_track(nt, Z_cutoff);
150 
151  // evaluate and cache track-vertex assignment chi**2 for all clusters and sum up Z
152  for (unsigned int k = 0; k < nv; k++) {
153  for (unsigned int j = tkfirstv_[k]; j < tkfirstv_[k + 1]; j++) {
154  const unsigned int i = tkmap_[j];
155  double arg = beta_over_2 * track_in_vertex_chsq(trackinfo_[i], xv_[k], yv_[k], zv_[k]);
156  if (arg < argmax) {
157  const double e = vdt::fast_exp(-arg);
158  tkweight_[j] = e; // must later be normalized by the proper Z_track[i]
159  Z_track[i] += e; // sum up exponentials to normalize
160  } else {
161  tkweight_[j] = 0.;
162  }
163  }
164  }
165 
166  // now we have the partition function, Z_i and can evaluate assignment probabilities (aka weights)
167  for (unsigned int j = 0; j < tkmap_.size(); j++) {
168  const unsigned int i = tkmap_[j];
169 #ifdef PVT_DEBUG
170  assert((i < nt) && "tkmap out of range");
171  assert((tkmap_.size() == tkweight_.size()) && "map and list not aliged");
172 #endif
173  tkweight_[j] /= Z_track[i];
174  }
175 }
assert(be >=bs)
A arg
Definition: Factorize.h:31
int nt
Definition: AMPTWrapper.h:42
double track_in_vertex_chsq(const TrackInfo &, const double, const double, const double)

◆ fit()

std::vector< TransientVertex > AdaptiveChisquarePrimaryVertexFitter::fit ( const std::vector< reco::TransientTrack > &  tracks,
const std::vector< TransientVertex > &  clusters,
const reco::BeamSpot beamspot,
const bool  useBeamConstraint 
)
overridevirtual

Implements PrimaryVertexFitterBase.

Definition at line 526 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References bsc_activity_cfg::clusters, multivertexfit_, FSQDQM_cfi::pvs, refit(), DiMuonV_cfg::tracks, HLT_2024v11_cff::useBeamConstraint, and vertices().

Referenced by trackingPlots.Iteration::modules().

529  {
530  if (multivertexfit_) {
532 
533  } else {
534  // fit the clusters one-by-one using the tracklist of the clusters (ignores the "tracks" argument)
535  std::vector<TransientVertex> pvs;
536  pvs.reserve(clusters.size());
537  for (auto &cluster : clusters) {
538  if (cluster.originalTracks().size() > (useBeamConstraint ? 0 : 1)) {
539  auto pv = refit(cluster, beamspot, useBeamConstraint);
540  if (pv.isValid()) {
541  pvs.emplace_back(pv);
542  }
543  }
544  }
545  return pvs;
546  }
547 }
std::vector< TransientVertex > vertices(const std::vector< reco::TransientTrack > &, const std::vector< TransientVertex > &, const reco::BeamSpot &, const bool)
TransientVertex refit(const TransientVertex &, const reco::BeamSpot &, const bool)

◆ get_inverse_beam_covariance()

AdaptiveChisquarePrimaryVertexFitter::Error3 AdaptiveChisquarePrimaryVertexFitter::get_inverse_beam_covariance ( const reco::BeamSpot beamspot)
protected

Definition at line 327 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References funct::pow().

Referenced by get_TransientVertex(), and update().

328  {
329  auto SBeam = beamspot.rotatedCovariance3D();
330  if (SBeam.Invert()) {
331  return SBeam;
332  } else {
333  edm::LogWarning("AdaptiveChisquarePrimaryVertexFitter")
334  << "Warning, beam-spot covariance matrix inversion failed " << std::endl;
335  Error3 S0;
336  S0(0, 0) = 1. / pow(beamspot.BeamWidthX(), 2);
337  S0(1, 1) = 1. / pow(beamspot.BeamWidthY(), 2);
338  S0(2, 2) = 1. / pow(beamspot.sigmaZ(), 2);
339  return S0;
340  }
341 }
Log< level::Warning, false > LogWarning
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ get_TransientVertex()

TransientVertex AdaptiveChisquarePrimaryVertexFitter::get_TransientVertex ( const unsigned int  k,
const std::vector< std::pair< unsigned int, float >> &  vertex_track_weights,
const std::vector< reco::TransientTrack > &  tracks,
const float  beam_weight,
const reco::BeamSpot beamspot 
)
protected

Definition at line 343 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References nano_mu_local_reco_cff::chi2, covv_, PVValHelper::dx, PVValHelper::dy, PVValHelper::dz, get_inverse_beam_covariance(), mps_fire::i, dqmdumpme::k, min_trackweight_, ntuplemaker::posError, track_in_vertex_chsq(), trackinfo_, DiMuonV_cfg::tracks, L1BJetProducer_cff::vtx, xv_, yv_, and zv_.

Referenced by refit(), and vertices().

348  {
349  const GlobalPoint pos(xv_[k], yv_[k], zv_[k]);
350  const GlobalError posError(
351  covv_[k](0, 0), covv_[k](1, 0), covv_[k](1, 1), covv_[k](2, 0), covv_[k](2, 1), covv_[k](2, 2));
352  float chi2 = 0.;
353  float vtx_ndof = -3.;
354  if (beam_weight > 0) {
355  // add beam-spot chi**2 and degrees of freedom
356  vtx_ndof = 3 * beam_weight;
357  const auto S = get_inverse_beam_covariance(beamspot);
358  const double dx = xv_[k] - beamspot.x0();
359  const double dy = yv_[k] - beamspot.y0();
360  const double dz = zv_[k] - beamspot.z0();
361  chi2 = beam_weight * (S(0, 0) * dx * dx + S(1, 1) * dy * dy + 2 * S(0, 1) * dx * dy + S(2, 2) * dz * dz +
362  2 * S(0, 2) * dx * dz + 2 * S(1, 2) * dy * dz);
363  }
364 
365  std::vector<reco::TransientTrack> vertex_tracks;
367  for (const auto &tk : vertex_track_weights) {
368  const unsigned int i = tk.first;
369  const float track_weight = tk.second;
370  if (track_weight >= min_trackweight_) {
371  vertex_tracks.emplace_back(tracks[i]);
372  trkWeightMap[tracks[i]] = track_weight;
373  vtx_ndof += 2 * track_weight;
374  chi2 += track_weight * track_in_vertex_chsq(trackinfo_[i], xv_[k], yv_[k], zv_[k]);
375  }
376  }
377 
378  auto vtx = TransientVertex(pos, posError, vertex_tracks, chi2, vtx_ndof);
379  vtx.weightMap(trkWeightMap);
380 
381  return vtx;
382 }
std::map< reco::TransientTrack, float > TransientTrackToFloatMap
double track_in_vertex_chsq(const TrackInfo &, const double, const double, const double)

◆ make_vtx_trk_map()

void AdaptiveChisquarePrimaryVertexFitter::make_vtx_trk_map ( const double  zrange_scale)
protected

Definition at line 93 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References funct::abs(), mps_fire::i, dqmdumpme::k, nt, tkfirstv_, tkmap_, tkweight_, trackinfo_, xv_, yv_, z_cutoff_, OfflinePixel3DPrimaryVertices_cfi::zrange, and zv_.

Referenced by vertices().

93  {
94  unsigned const int nv = xv_.size();
95  unsigned const int nt = trackinfo_.size();
96 
97 #ifdef PVTX_DEBUG
98  if (nv < 1) {
99  edm::LogWarning("AdaptiveChisquarePrimaryVertexFit") << " empty vertex list with " << nt << " tracks" << std::endl;
100  return;
101  }
102 #endif
103 
104  // parallel lists for track to vertex mapping, tracks are not sorted
105  tkmap_.clear(); // index in trackinfo_
106  tkweight_.clear(); // weight in vertex
107  tkfirstv_.clear(); // each vertex k owns a section of those list : tkfirstv_[k] .. tkfirstv_[k+1]-1
108 
109  if (nv == 1) {
110  // always accept all tracks for a single vertex fit
111  tkfirstv_.push_back(0);
112  tkfirstv_.push_back(nt);
113  tkweight_.assign(nt, 0.);
114  tkmap_.reserve(nt);
115  for (unsigned int i = 0; i < nt; i++) {
116  tkmap_.emplace_back(i);
117  }
118  return;
119  }
120 
121  // n > 1
122  tkmap_.reserve(nv * 100);
123  tkweight_.reserve(nv * 100);
124  for (unsigned int k = 0; k < nv; k++) {
125  tkfirstv_.emplace_back(tkmap_.size());
126  for (unsigned int i = 0; i < nt; i++) {
127  auto &ti = trackinfo_[i];
128  const double zrange = zrange_scale * ti.dzError;
129  if (std::abs(zv_[k] - ti.zpca) < z_cutoff_) {
130  const double dztrk = ti.b2 + xv_[k] * ti.H2[0] + yv_[k] * ti.H2[1] - zv_[k];
131  if (std::abs(dztrk) < zrange) {
132  tkmap_.emplace_back(i);
133  tkweight_.emplace_back(0.);
134  }
135  }
136  }
137  }
138  tkfirstv_.emplace_back(tkmap_.size()); // extra entry, simplifies loops, every vertex has a "successor" now
139 }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int nt
Definition: AMPTWrapper.h:42
Log< level::Warning, false > LogWarning

◆ refit()

TransientVertex AdaptiveChisquarePrimaryVertexFitter::refit ( const TransientVertex cluster,
const reco::BeamSpot beamspot,
const bool  useBeamConstraint 
)
protected

Definition at line 474 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References HLTObjectsMonitor_cfi::beamspot, covv_, dumpMFGeometry_cfg::delta, fill_trackinfo(), fill_weights(), get_TransientVertex(), mps_fire::i, nt, TransientVertex::originalTracks(), TransientVertex::position(), tkfirstv_, tkmap_, tkweight_, update(), HLT_2024v11_cff::useBeamConstraint, xv_, yv_, PV3DBase< T, PVType, FrameType >::z(), and zv_.

Referenced by fit().

476  {
477  // fit a single vertex using all tracks in the tracklist
478  const unsigned int nt = cluster.originalTracks().size();
479  const int max_iterations = 50;
480 
481  // initialize, vectors with size=1 here to avoid code duplication from the multivertex case in update()
482  const double zclu = cluster.position().z();
483  xv_ = {beamspot.x(zclu)};
484  yv_ = {beamspot.y(zclu)};
485  zv_ = {zclu};
486  tkfirstv_ = {0, nt};
487  covv_.clear();
488 
490  tkweight_.assign(nt, 0.);
491  tkmap_.clear();
492  tkmap_.reserve(nt);
493  for (unsigned int i = 0; i < nt; i++) {
494  tkmap_.emplace_back(i); // trivial map for single vertex fits
495  }
496 
497  float beam_weight = useBeamConstraint ? 1. : 0.;
498 
499  double delta = 0;
500  unsigned int nit = 0;
501  while ((nit == 0) || ((delta > 0.0001) && (nit < max_iterations))) {
503  delta = update(beamspot, beam_weight, false);
504  nit++;
505  }
506 
507  if ((nit >= max_iterations) && (delta > 0.1)) {
508  edm::LogWarning("AdaptiveChisquarePrimaryVertexFitter")
509  << "single vertex fit, iteration limit reached " << nit << " last delta = " << delta << std::endl
510  << " nt = " << cluster.originalTracks().size() << std::endl;
511  }
512 
513  // fill the covariance matrices
514  update(beamspot, beam_weight, true);
515 
516  // put the result into a transient vertex
517  std::vector<std::pair<unsigned int, float>> vertex_track_weights;
518  for (unsigned int i = 0; i < nt; i++) {
519  vertex_track_weights.emplace_back(i, tkweight_[i]);
520  }
521 
522  return get_TransientVertex(0, vertex_track_weights, cluster.originalTracks(), beam_weight, beamspot);
523 }
void fill_weights(const reco::BeamSpot &, const double beta=1.)
GlobalPoint position() const
T z() const
Definition: PV3DBase.h:61
void fill_trackinfo(const std::vector< reco::TransientTrack > &, const reco::BeamSpot &)
std::vector< reco::TransientTrack > const & originalTracks() const
int nt
Definition: AMPTWrapper.h:42
double update(const reco::BeamSpot &, float beam_weight, const bool fill_covariances=false)
TransientVertex get_TransientVertex(const unsigned int, const std::vector< std::pair< unsigned int, float >> &, const std::vector< reco::TransientTrack > &, const float, const reco::BeamSpot &)
Log< level::Warning, false > LogWarning

◆ remove_vertex()

void AdaptiveChisquarePrimaryVertexFitter::remove_vertex ( unsigned int  k)
protected

Definition at line 220 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References covv_, dqmdumpme::k, tkfirstv_, tkmap_, tkweight_, xv_, yv_, and zv_.

Referenced by clean().

220  {
221  // remove a vertex or rather merge it with it's neighbour
222  // used for multi-vertex fits only
223  unsigned int nv = xv_.size();
224  if (nv < 2)
225  return;
226 
227  // 1) remove the vertex from the vertex list
228  xv_.erase(xv_.begin() + k);
229  yv_.erase(yv_.begin() + k);
230  zv_.erase(zv_.begin() + k);
231  covv_.erase(covv_.begin() + k);
232 
233  // 2) adjust the track-map map
234  // 2a) remove the map entries that belong the deleted vertex
235  const unsigned int num_erased_map_entries = tkfirstv_[k + 1] - tkfirstv_[k];
236  tkmap_.erase(tkmap_.begin() + tkfirstv_[k], tkmap_.begin() + tkfirstv_[k + 1]);
237  tkweight_.erase(tkweight_.begin() + tkfirstv_[k], tkweight_.begin() + tkfirstv_[k + 1]);
238  // 2b) adjust pointers for the following vertices, including the dummy entry behind the last (now [nv-1])
239  for (unsigned int k1 = k + 1; k1 < nv + 1; k1++) {
240  tkfirstv_[k1] -= num_erased_map_entries;
241  }
242  // 2c) erase the pointer of the removed vertex
243  tkfirstv_.erase(tkfirstv_.begin() + k);
244 }

◆ track_in_vertex_chsq()

double AdaptiveChisquarePrimaryVertexFitter::track_in_vertex_chsq ( const TrackInfo ti,
const double  xvtx,
const double  yvtx,
const double  zvtx 
)
protected

Definition at line 21 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References cms::cuda::assert(), AdaptiveChisquarePrimaryVertexFitter::TrackInfo::b1, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::b2, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::H1, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::H2, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::S11, AdaptiveChisquarePrimaryVertexFitter::TrackInfo::S12, and AdaptiveChisquarePrimaryVertexFitter::TrackInfo::S22.

Referenced by fill_weights(), and get_TransientVertex().

24  {
25  double F1 = ti.b1 + xvtx * ti.H1[0] + yvtx * ti.H1[1]; // using H1[2]=0
26  double F2 = ti.b2 + xvtx * ti.H2[0] + yvtx * ti.H2[1] - zvtx; // using H2[2]=-1
27  double chsq = F1 * F1 * ti.S11 + F2 * F2 * ti.S22 + 2. * F1 * F2 * ti.S12;
28 #ifdef PVTX_DEBUG
29  assert((chsq >= 0) && " negative chi**2");
30 #endif
31  return chsq;
32 }
assert(be >=bs)

◆ update()

double AdaptiveChisquarePrimaryVertexFitter::update ( const reco::BeamSpot beamspot,
float  beam_weight,
const bool  fill_covariances = false 
)
protected

Definition at line 246 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References funct::abs(), HltBtagPostValidation_cff::c, covv_, MillePedeFileConverter_cfg::e, get_inverse_beam_covariance(), mps_fire::i, dqmiolumiharvest::j, dqmdumpme::k, MainPageGenerator::l, SiStripPI::max, nt, tkfirstv_, tkmap_, tkweight_, trackinfo_, w(), xv_, yv_, and zv_.

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

248  {
249  double rho_vtx = 0;
250  double delta_z = 0;
251  double delta_x = 0;
252  double delta_y = 0;
253  unsigned const int nt = trackinfo_.size();
254  unsigned const int nv = xv_.size();
255  if (fill_covariances) {
256  covv_.clear();
257  }
258 
259  // initial value for S, 0 or inverse of the beamspot covariance matrix
260  Error3 S0;
261  double c_beam[3] = {0, 0, 0};
262  if (beam_weight > 0) {
264  for (unsigned int j = 0; j < 3; j++) {
265  c_beam[j] = -(S0(j, 0) * beamspot.x0() + S0(j, 1) * beamspot.y0() + S0(j, 2) * beamspot.z0());
266  }
267  }
268 
269  for (unsigned int k = 0; k < nv; k++) {
270  rho_vtx = 0;
271  Error3 S(S0);
272  // sum track contributions
273  double c[3] = {c_beam[0], c_beam[1], c_beam[2]};
274  for (unsigned int j = tkfirstv_[k]; j < tkfirstv_[k + 1]; j++) {
275  const unsigned int i = tkmap_[j];
276  const auto w = tkweight_[j];
277  rho_vtx += w;
278  S += w * trackinfo_[i].C;
279  for (unsigned int l = 0; l < 3; l++) {
280  c[l] += w * trackinfo_[i].g[l];
281  }
282  }
283 
284 #ifdef PVTX_DEBUG
285  if ((fabs(S(1, 2) - S(2, 1)) > 1e-3) || (fabs(S(0, 2) - S(2, 0)) > 1e-3) || (fabs(S(0, 1) - S(1, 0)) > 1e-3) ||
286  (S(0, 0) <= 0) || (S(0, 0) <= 0) || (S(0, 0) <= 0)) {
287  edm::LogWarning("AdaptiveChisquarePrimaryVertexFitter") << "update() bad S-matrix S=" << std::endl
288  << S << std::endl;
289  edm::LogWarning("AdaptiveChisquarePrimaryVertexFitter")
290  << "n-vertex = " << nv << " n-track = " << nt << std::endl;
291  }
292 #endif
293 
294  const auto xold = xv_[k];
295  const auto yold = yv_[k];
296  const auto zold = zv_[k];
297 
298  if (S.Invert()) {
299  xv_[k] = -(S(0, 0) * c[0] + S(0, 1) * c[1] + S(0, 2) * c[2]);
300  yv_[k] = -(S(1, 0) * c[0] + S(1, 1) * c[1] + S(1, 2) * c[2]);
301  zv_[k] = -(S(2, 0) * c[0] + S(2, 1) * c[1] + S(2, 2) * c[2]);
302  if (fill_covariances) {
303  covv_.emplace_back(S);
304  }
305  } else {
306  edm::LogWarning("AdaptiveChisquarePrimaryVertexFitter") << "update() Matrix inversion failed" << S << std::endl;
307  if (fill_covariances) {
308  Error3 covv_dummy;
309  covv_dummy(0, 0) = 100.;
310  covv_dummy(1, 1) = 100.;
311  covv_dummy(2, 2) = 100.;
312  covv_.emplace_back(covv_dummy);
313  }
314  }
315 
316  if ((nt > 1) && (rho_vtx > 1.0)) {
317  delta_x = std::max(delta_x, std::abs(xv_[k] - xold));
318  delta_y = std::max(delta_y, std::abs(yv_[k] - yold));
319  delta_z = std::max(delta_z, std::abs(zv_[k] - zold));
320  }
321 
322  } // vertex loop
323 
324  return std::max(delta_z, std::max(delta_x, delta_y));
325 }
T w() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int nt
Definition: AMPTWrapper.h:42
Log< level::Warning, false > LogWarning

◆ verify()

void AdaptiveChisquarePrimaryVertexFitter::verify ( )
inlineprotected

Definition at line 31 of file AdaptiveChisquarePrimaryVertexFitter.h.

References cms::cuda::assert(), mps_fire::i, dqmiolumiharvest::j, dqmdumpme::k, nt, tkfirstv_, tkmap_, tkweight_, trackinfo_, xv_, yv_, and zv_.

31  { // DEBUG only
32  unsigned int nt = trackinfo_.size();
33  unsigned int nv = xv_.size();
34  assert((yv_.size() == nv) && "yv size");
35  assert((zv_.size() == nv) && "zv size");
36  assert((tkfirstv_.size() == (nv + 1)) && "tkfirstv size");
37  assert((tkmap_.size() == tkweight_.size()) && "tkmapsize <> tkweightssize");
38  for (unsigned int k = 0; k < nv; k++) {
39  assert((tkfirstv_[k] < tkweight_.size()) && "tkfirst[k]");
40  assert((tkfirstv_[k + 1] <= tkweight_.size()) && "tkfirst[k+1]");
41  assert((tkfirstv_[k] <= tkfirstv_[k + 1]) && "tkfirst[k/k+1]");
42  for (unsigned int j = tkfirstv_[k]; j < tkfirstv_[k + 1]; k++) {
43  assert((j < tkmap_.size()) && "illegal tkfirst entry");
44  unsigned int i = tkmap_[j];
45  assert((i < nt) && "illegal tkmap entry");
46  assert((tkweight_[i] >= 0) && "negative tkweight or nan");
47  assert((tkweight_[i] <= 1) && "tkweight > 1 or nan");
48  }
49  }
50  };
assert(be >=bs)
int nt
Definition: AMPTWrapper.h:42

◆ vertices()

std::vector< TransientVertex > AdaptiveChisquarePrimaryVertexFitter::vertices ( const std::vector< reco::TransientTrack > &  tracks,
const std::vector< TransientVertex > &  clusters,
const reco::BeamSpot beamspot,
const bool  useBeamConstraint 
)
protected

Definition at line 384 of file AdaptiveChisquarePrimaryVertexFitter.cc.

References clean(), bsc_activity_cfg::clusters, covv_, dumpMFGeometry_cfg::delta, fill_trackinfo(), fill_weights(), get_TransientVertex(), mps_fire::i, dqmiolumiharvest::j, dqmdumpme::k, make_vtx_trk_map(), FSQDQM_cfi::pvs, tkfirstv_, tkmap_, tkweight_, trackinfo_, DiMuonV_cfg::tracks, update(), HLT_2024v11_cff::useBeamConstraint, xv_, yv_, and zv_.

Referenced by fit().

388  {
389  // simultaneous fit of all vertices in the input list
390 
391  const int max_iterations = 50;
392 
393  // initialize the vertices
394  const unsigned int nv = clusters.size();
395  xv_.clear();
396  xv_.reserve(nv);
397  yv_.clear();
398  yv_.reserve(nv);
399  zv_.clear();
400  zv_.reserve(nv);
401  tkfirstv_.clear();
402  tkfirstv_.reserve(nv + 1);
403  covv_.clear();
404  covv_.reserve(nv);
405 
406  // seeds
407  for (auto &clu : clusters) {
408  const double zclu = clu.position().z();
409  xv_.emplace_back(beamspot.x(zclu));
410  yv_.emplace_back(beamspot.y(zclu));
411  zv_.emplace_back(zclu);
412  }
413 
415 
416  make_vtx_trk_map(5.); // use tracks within 5 sigma windows (if that is less than z_cutoff_)
417 
418  float beam_weight = useBeamConstraint ? 1. : 0.;
419 
420  double delta = 0;
421  unsigned int nit = 0;
422  while ((nit == 0) || ((delta > 0.0001) && (nit < max_iterations))) {
424  delta = update(beamspot, beam_weight, false);
425  nit++;
426  }
427  if ((nit >= max_iterations) && (delta > 0.01)) {
428  edm::LogWarning("AdaptiveChisquarePrimaryVertexFitter")
429  << "iteration limit reached " << nit << " last delta = " << delta << std::endl
430  << " nv = " << nv << " nt = " << tracks.size() << std::endl;
431  }
432 
433  // may need to remove collapsed vertices
434  nit = 0;
435  while ((xv_.size() > 1) && (nit < max_iterations) && (clean())) {
437  update(beamspot, beam_weight, false);
438  nit++;
439  }
440 
441  // fill the covariance matrices
442  update(beamspot, beam_weight, true);
443 
444  // assign tracks to vertices
445  std::vector<unsigned int> track_to_vertex(trackinfo_.size(), nv);
446  // for each track identify the vertex that wants it most
447  std::vector<double> maxweight(trackinfo_.size(), -1.);
448  for (unsigned int k = 0; k < nv; k++) {
449  for (unsigned int j = tkfirstv_[k]; j < tkfirstv_[k + 1]; j++) {
450  const unsigned int i = tkmap_[j];
451  if (tkweight_[j] > maxweight[i]) {
452  maxweight[i] = tkweight_[j];
453  track_to_vertex[i] = k;
454  }
455  }
456  }
457 
458  // fill the fit result into transient vertices
459  std::vector<TransientVertex> pvs;
460  for (unsigned int k = 0; k < xv_.size(); k++) {
461  std::vector<std::pair<unsigned int, float>> vertex_tracks_weights;
462  for (unsigned int j = tkfirstv_[k]; j < tkfirstv_[k + 1]; j++) {
463  unsigned int i = tkmap_[j];
464  if (track_to_vertex[i] == k) {
465  vertex_tracks_weights.emplace_back(tkmap_[j], tkweight_[j]);
466  }
467  }
468  pvs.emplace_back(get_TransientVertex(k, vertex_tracks_weights, tracks, beam_weight, beamspot));
469  }
470 
471  return pvs;
472 }
void fill_weights(const reco::BeamSpot &, const double beta=1.)
void fill_trackinfo(const std::vector< reco::TransientTrack > &, const reco::BeamSpot &)
double update(const reco::BeamSpot &, float beam_weight, const bool fill_covariances=false)
TransientVertex get_TransientVertex(const unsigned int, const std::vector< std::pair< unsigned int, float >> &, const std::vector< reco::TransientTrack > &, const float, const reco::BeamSpot &)
Log< level::Warning, false > LogWarning

Member Data Documentation

◆ chi_cutoff_

double AdaptiveChisquarePrimaryVertexFitter::chi_cutoff_
protected

◆ covv_

std::vector<Error3> AdaptiveChisquarePrimaryVertexFitter::covv_
protected

◆ min_trackweight_

double AdaptiveChisquarePrimaryVertexFitter::min_trackweight_
protected

◆ multivertexfit_

double AdaptiveChisquarePrimaryVertexFitter::multivertexfit_
protected

◆ tkfirstv_

std::vector<unsigned int> AdaptiveChisquarePrimaryVertexFitter::tkfirstv_
protected

◆ tkmap_

std::vector<unsigned int> AdaptiveChisquarePrimaryVertexFitter::tkmap_
protected

◆ tkweight_

std::vector<double> AdaptiveChisquarePrimaryVertexFitter::tkweight_
protected

◆ trackinfo_

std::vector<TrackInfo> AdaptiveChisquarePrimaryVertexFitter::trackinfo_
protected

◆ xv_

std::vector<double> AdaptiveChisquarePrimaryVertexFitter::xv_
protected

◆ yv_

std::vector<double> AdaptiveChisquarePrimaryVertexFitter::yv_
protected

◆ z_cutoff_

double AdaptiveChisquarePrimaryVertexFitter::z_cutoff_
protected

◆ zv_

std::vector<double> AdaptiveChisquarePrimaryVertexFitter::zv_
protected