CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes
TotemRPLocalTrackFitterAlgorithm Class Reference

Algorithm for fitting tracks through a single RP. More...

#include <TotemRPLocalTrackFitterAlgorithm.h>

Classes

struct  RPDetCoordinateAlgebraObjs
 

Public Member Functions

bool fitTrack (const edm::DetSetVector< TotemRPRecHit > &hits, double z_0, const CTPPSGeometry &tot_geom, TotemRPLocalTrack &fitted_track)
 performs the track fit, returns true if successful More...
 
void reset ()
 Resets the reconstruction-data cache. More...
 
 TotemRPLocalTrackFitterAlgorithm (const edm::ParameterSet &conf)
 

Private Member Functions

RPDetCoordinateAlgebraObjsgetDetAlgebraData (unsigned int det_id, const CTPPSGeometry &tot_rp_geom)
 
void multiplyByDiagonalInPlace (TMatrixD &mt, const TVectorD &diag)
 A matrix multiplication shorthand. More...
 
RPDetCoordinateAlgebraObjs prepareReconstAlgebraData (unsigned int det_id, const CTPPSGeometry &tot_rp_geom)
 Build the reconstruction data. More...
 

Static Private Member Functions

static TVector3 convert3vector (const CTPPSGeometry::Vector &v)
 

Private Attributes

std::unordered_map< unsigned int, RPDetCoordinateAlgebraObjsdet_data_map_
 A cache of reconstruction data. Must be reset every time the geometry chagnges. More...
 
RPTopology rp_topology_
 

Detailed Description

Algorithm for fitting tracks through a single RP.

Definition at line 31 of file TotemRPLocalTrackFitterAlgorithm.h.

Constructor & Destructor Documentation

◆ TotemRPLocalTrackFitterAlgorithm()

TotemRPLocalTrackFitterAlgorithm::TotemRPLocalTrackFitterAlgorithm ( const edm::ParameterSet conf)

Definition at line 22 of file TotemRPLocalTrackFitterAlgorithm.cc.

Member Function Documentation

◆ convert3vector()

static TVector3 TotemRPLocalTrackFitterAlgorithm::convert3vector ( const CTPPSGeometry::Vector v)
inlinestaticprivate

Definition at line 73 of file TotemRPLocalTrackFitterAlgorithm.h.

◆ fitTrack()

bool TotemRPLocalTrackFitterAlgorithm::fitTrack ( const edm::DetSetVector< TotemRPRecHit > &  hits,
double  z_0,
const CTPPSGeometry tot_geom,
TotemRPLocalTrack fitted_track 
)

performs the track fit, returns true if successful

Definition at line 63 of file TotemRPLocalTrackFitterAlgorithm.cc.

67  {
68  fitted_track.setValid(false);
69 
70  // bind hits with their algebra objects
71  struct HitWithAlg {
72  unsigned int detId;
73  const TotemRPRecHit *hit;
74  RPDetCoordinateAlgebraObjs *alg;
75  };
76 
77  vector<HitWithAlg> applicable_hits;
78 
79  for (auto &ds : hits) {
80  unsigned int detId = ds.detId();
81 
82  for (auto &h : ds) {
83  RPDetCoordinateAlgebraObjs *alg = getDetAlgebraData(detId, tot_geom);
84  if (alg->available_)
85  applicable_hits.push_back({detId, &h, alg});
86  }
87  }
88 
89  if (applicable_hits.size() < 5)
90  return false;
91 
92  TMatrixD H(applicable_hits.size(), 4);
93  TVectorD V(applicable_hits.size());
94  TVectorD V_inv(applicable_hits.size());
95  TVectorD U(applicable_hits.size());
96 
97  for (unsigned int i = 0; i < applicable_hits.size(); ++i) {
98  RPDetCoordinateAlgebraObjs *alg_obj = applicable_hits[i].alg;
99 
100  H(i, 0) = alg_obj->readout_direction_.X();
101  H(i, 1) = alg_obj->readout_direction_.Y();
102  double delta_z = alg_obj->centre_of_det_global_position_.Z() - z_0;
103  H(i, 2) = alg_obj->readout_direction_.X() * delta_z;
104  H(i, 3) = alg_obj->readout_direction_.Y() * delta_z;
105  double var = applicable_hits[i].hit->sigma();
106  var *= var;
107  V[i] = var;
108  V_inv[i] = 1.0 / var;
109  U[i] = applicable_hits[i].hit->position() - alg_obj->rec_u_0_;
110  }
111 
112  TMatrixD H_T_V_inv(TMatrixD::kTransposed, H);
113  multiplyByDiagonalInPlace(H_T_V_inv, V_inv);
114  TMatrixD V_a(H_T_V_inv);
115  TMatrixD V_a_mult(V_a, TMatrixD::kMult, H);
116  try {
117  V_a_mult.Invert();
118  } catch (cms::Exception &e) {
119  LogError("TotemRPLocalTrackFitterAlgorithm") << "Error in TotemRPLocalTrackFitterAlgorithm::fitTrack > "
120  << "Fit matrix is singular. Skipping.";
121  return false;
122  }
123 
124  TMatrixD u_to_a(V_a_mult, TMatrixD::kMult, H_T_V_inv);
125  TVectorD a(U);
126  a *= u_to_a;
127 
128  fitted_track.setZ0(z_0);
129  fitted_track.setParameterVector(a);
130  fitted_track.setCovarianceMatrix(V_a_mult);
131 
132  double Chi_2 = 0;
133  for (unsigned int i = 0; i < applicable_hits.size(); ++i) {
134  RPDetCoordinateAlgebraObjs *alg_obj = applicable_hits[i].alg;
135  TVector2 readout_dir = alg_obj->readout_direction_;
136  double det_z = alg_obj->centre_of_det_global_position_.Z();
137  double sigma_str = applicable_hits[i].hit->sigma();
138  double sigma_str_2 = sigma_str * sigma_str;
139  TVector2 fited_det_xy_point = fitted_track.trackPoint(det_z);
140  double U_readout = applicable_hits[i].hit->position() - alg_obj->rec_u_0_;
141  double U_fited = (readout_dir *= fited_det_xy_point);
142  double residual = U_fited - U_readout;
143  TMatrixD V_T_Cov_X_Y(1, 2);
144  V_T_Cov_X_Y(0, 0) = readout_dir.X();
145  V_T_Cov_X_Y(0, 1) = readout_dir.Y();
146  TMatrixD V_T_Cov_X_Y_mult(V_T_Cov_X_Y, TMatrixD::kMult, fitted_track.trackPointInterpolationCovariance(det_z));
147  double fit_strip_var = V_T_Cov_X_Y_mult(0, 0) * readout_dir.X() + V_T_Cov_X_Y_mult(0, 1) * readout_dir.Y();
148  double pull_normalization = sqrt(sigma_str_2 - fit_strip_var);
149  double pull = residual / pull_normalization;
150 
151  Chi_2 += residual * residual / sigma_str_2;
152 
154  *(applicable_hits[i].hit), TVector3(fited_det_xy_point.X(), fited_det_xy_point.Y(), det_z), residual, pull);
155  fitted_track.addHit(applicable_hits[i].detId, hit_point);
156  }
157 
158  fitted_track.setChiSquared(Chi_2);
159  fitted_track.setValid(true);
160  return true;

References a, TotemRPLocalTrack::addHit(), TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::available_, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::centre_of_det_global_position_, MillePedeFileConverter_cfg::e, h, class-composition::H, hfClusterShapes_cfi::hits, mps_fire::i, JetComb::kMult, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::readout_direction_, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::rec_u_0_, TotemRPLocalTrack::setChiSquared(), TotemRPLocalTrack::setCovarianceMatrix(), TotemRPLocalTrack::setParameterVector(), TotemRPLocalTrack::setValid(), TotemRPLocalTrack::setZ0(), mathSSE::sqrt(), TotemRPLocalTrack::trackPoint(), TotemRPLocalTrack::trackPointInterpolationCovariance(), mitigatedMETSequence_cff::U, cms::cuda::V, and trigObjTnPSource_cfi::var.

Referenced by TotemRPLocalTrackFitter::produce().

◆ getDetAlgebraData()

TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs * TotemRPLocalTrackFitterAlgorithm::getDetAlgebraData ( unsigned int  det_id,
const CTPPSGeometry tot_rp_geom 
)
private

Returns the reconstruction data for the chosen detector from the cache DetReconstructionDataMap. If it is not yet in the cache, calls PrepareReconstAlgebraData to make it.

Definition at line 50 of file TotemRPLocalTrackFitterAlgorithm.cc.

52  {
53  auto it = det_data_map_.find(det_id);
54  if (it != det_data_map_.end()) {
55  return &(it->second);
56  } else {
57  det_data_map_[det_id] = prepareReconstAlgebraData(det_id, tot_rp_geom);
58  return &det_data_map_[det_id];
59  }

◆ multiplyByDiagonalInPlace()

void TotemRPLocalTrackFitterAlgorithm::multiplyByDiagonalInPlace ( TMatrixD &  mt,
const TVectorD &  diag 
)
private

A matrix multiplication shorthand.

Definition at line 164 of file TotemRPLocalTrackFitterAlgorithm.cc.

165  {
166  for (int i = 0; i < mt.GetNrows(); ++i) {
167  for (int j = 0; j < mt.GetNcols(); ++j) {
168  mt[i][j] *= diag[j];
169  }
170  }

References mps_fire::i, dqmiolumiharvest::j, and TtSemiLepEvtBuilder_cfi::mt.

◆ prepareReconstAlgebraData()

TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs TotemRPLocalTrackFitterAlgorithm::prepareReconstAlgebraData ( unsigned int  det_id,
const CTPPSGeometry tot_rp_geom 
)
private

Build the reconstruction data.

Definition at line 31 of file TotemRPLocalTrackFitterAlgorithm.cc.

32  {
33  RPDetCoordinateAlgebraObjs det_algebra_obj;
34 
35  det_algebra_obj.centre_of_det_global_position_ = convert3vector(tot_rp_geom.sensorTranslation(det_id));
36 
37  TVector3 rd_dir = convert3vector(tot_rp_geom.localToGlobalDirection(det_id, rp_topology_.GetStripReadoutAxisDir()));
38 
39  TVector2 v(rd_dir.X(), rd_dir.Y());
40  det_algebra_obj.readout_direction_ = v.Unit();
41  det_algebra_obj.rec_u_0_ = 0.0;
42  det_algebra_obj.available_ = true;
43  det_algebra_obj.rec_u_0_ =
44  -(det_algebra_obj.readout_direction_ * det_algebra_obj.centre_of_det_global_position_.XYvector());
45 
46  return det_algebra_obj;

References TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::available_, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::centre_of_det_global_position_, CTPPSGeometry::localToGlobalDirection(), TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::readout_direction_, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::rec_u_0_, CTPPSGeometry::sensorTranslation(), and findQualityFiles::v.

◆ reset()

void TotemRPLocalTrackFitterAlgorithm::reset ( void  )

Resets the reconstruction-data cache.

Definition at line 26 of file TotemRPLocalTrackFitterAlgorithm.cc.

Referenced by TotemRPLocalTrackFitter::produce().

Member Data Documentation

◆ det_data_map_

std::unordered_map<unsigned int, RPDetCoordinateAlgebraObjs> TotemRPLocalTrackFitterAlgorithm::det_data_map_
private

A cache of reconstruction data. Must be reset every time the geometry chagnges.

Definition at line 59 of file TotemRPLocalTrackFitterAlgorithm.h.

◆ rp_topology_

RPTopology TotemRPLocalTrackFitterAlgorithm::rp_topology_
private

Definition at line 61 of file TotemRPLocalTrackFitterAlgorithm.h.

TotemRPLocalTrack::FittedRecHit
Definition: TotemRPLocalTrack.h:39
TotemRPRecHit
Reconstructed hit in TOTEM RP.
Definition: TotemRPRecHit.h:17
class-composition.H
H
Definition: class-composition.py:31
cms::cuda::V
cudaStream_t T uint32_t const T *__restrict__ const uint32_t *__restrict__ uint32_t int cudaStream_t V
Definition: HistoContainer.h:99
mps_fire.i
i
Definition: mps_fire.py:355
TotemRPLocalTrackFitterAlgorithm::multiplyByDiagonalInPlace
void multiplyByDiagonalInPlace(TMatrixD &mt, const TVectorD &diag)
A matrix multiplication shorthand.
Definition: TotemRPLocalTrackFitterAlgorithm.cc:164
TotemRPLocalTrack::setChiSquared
void setChiSquared(double &chiSquared)
Definition: TotemRPLocalTrack.h:116
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
TotemRPLocalTrack::addHit
void addHit(unsigned int detId, const FittedRecHit &hit)
Definition: TotemRPLocalTrack.h:82
h
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
Definition: L1TUtmAlgorithmRcd.h:4
TotemRPLocalTrack::trackPoint
TVector2 trackPoint(double z) const
returns (x, y) vector
Definition: TotemRPLocalTrack.h:123
findQualityFiles.v
v
Definition: findQualityFiles.py:179
trigObjTnPSource_cfi.var
var
Definition: trigObjTnPSource_cfi.py:21
TotemRPLocalTrackFitterAlgorithm::det_data_map_
std::unordered_map< unsigned int, RPDetCoordinateAlgebraObjs > det_data_map_
A cache of reconstruction data. Must be reset every time the geometry chagnges.
Definition: TotemRPLocalTrackFitterAlgorithm.h:59
TotemRPLocalTrack::trackPointInterpolationCovariance
TMatrixD trackPointInterpolationCovariance(double z) const
Definition: TotemRPLocalTrack.cc:13
h
RPTopology::GetStripReadoutAxisDir
const Vector & GetStripReadoutAxisDir() const
Definition: RPTopology.h:29
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
TotemRPLocalTrackFitterAlgorithm::rp_topology_
RPTopology rp_topology_
Definition: TotemRPLocalTrackFitterAlgorithm.h:61
CTPPSGeometry::localToGlobalDirection
Vector localToGlobalDirection(unsigned int id, const Vector &) const
Definition: CTPPSGeometry.cc:170
TotemRPLocalTrack::setZ0
void setZ0(double z0)
Definition: TotemRPLocalTrack.h:95
mitigatedMETSequence_cff.U
U
Definition: mitigatedMETSequence_cff.py:36
CTPPSGeometry::sensorTranslation
Vector sensorTranslation(unsigned int id) const
Definition: CTPPSGeometry.cc:182
edm::LogError
Definition: MessageLogger.h:183
a
double a
Definition: hdecay.h:119
TotemRPLocalTrackFitterAlgorithm::getDetAlgebraData
RPDetCoordinateAlgebraObjs * getDetAlgebraData(unsigned int det_id, const CTPPSGeometry &tot_rp_geom)
Definition: TotemRPLocalTrackFitterAlgorithm.cc:50
JetComb::kMult
Definition: TtSemiLepJetComb.h:35
TotemRPLocalTrackFitterAlgorithm::convert3vector
static TVector3 convert3vector(const CTPPSGeometry::Vector &v)
Definition: TotemRPLocalTrackFitterAlgorithm.h:73
TotemRPLocalTrack::setValid
void setValid(bool valid)
Definition: TotemRPLocalTrack.h:135
TotemRPLocalTrack::setParameterVector
void setParameterVector(const TVectorD &track_params_vector)
Definition: TotemRPLocalTrack.cc:59
TotemRPLocalTrackFitterAlgorithm::prepareReconstAlgebraData
RPDetCoordinateAlgebraObjs prepareReconstAlgebraData(unsigned int det_id, const CTPPSGeometry &tot_rp_geom)
Build the reconstruction data.
Definition: TotemRPLocalTrackFitterAlgorithm.cc:31
TtSemiLepEvtBuilder_cfi.mt
mt
Definition: TtSemiLepEvtBuilder_cfi.py:47
cms::Exception
Definition: Exception.h:70
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
TotemRPLocalTrack::setCovarianceMatrix
void setCovarianceMatrix(const TMatrixD &par_covariance_matrix)
Definition: TotemRPLocalTrack.cc:78
hit
Definition: SiStripHitEffFromCalibTree.cc:88
MillePedeFileConverter_cfg.e
e
Definition: MillePedeFileConverter_cfg.py:37