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 TotemRPGeometry &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 TotemRPGeometry &tot_rp_geom)
 
void multiplyByDiagonalInPlace (TMatrixD &mt, const TVectorD &diag)
 A matrix multiplication shorthand. More...
 
RPDetCoordinateAlgebraObjs prepareReconstAlgebraData (unsigned int det_id, const TotemRPGeometry &tot_rp_geom)
 Build the reconstruction data. More...
 

Static Private Member Functions

static TVector3 convert3vector (const CLHEP::Hep3Vector &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 32 of file TotemRPLocalTrackFitterAlgorithm.h.

Constructor & Destructor Documentation

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

Definition at line 23 of file TotemRPLocalTrackFitterAlgorithm.cc.

24 {
25 }

Member Function Documentation

static TVector3 TotemRPLocalTrackFitterAlgorithm::convert3vector ( const CLHEP::Hep3Vector &  v)
inlinestaticprivate

Definition at line 67 of file TotemRPLocalTrackFitterAlgorithm.h.

68  {
69  return TVector3(v.x(),v.y(),v.z()) ;
70  }
bool TotemRPLocalTrackFitterAlgorithm::fitTrack ( const edm::DetSetVector< TotemRPRecHit > &  hits,
double  z_0,
const TotemRPGeometry tot_geom,
TotemRPLocalTrack fitted_track 
)

performs the track fit, returns true if successful

Definition at line 79 of file TotemRPLocalTrackFitterAlgorithm.cc.

References a, TotemRPLocalTrack::addHit(), TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::available_, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::centre_of_det_global_position_, MillePedeFileConverter_cfg::e, TotemRPLocalTrack::getTrackPoint(), h, class-composition::H, 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::trackPointInterpolationCovariance(), mitigatedMETSequence_cff::U, and JetChargeProducer_cfi::var.

Referenced by TotemRPLocalTrackFitter::produce().

81 {
82  fitted_track.setValid(false);
83 
84  // bind hits with their algebra objects
85  struct HitWithAlg
86  {
87  unsigned int detId;
88  const TotemRPRecHit *hit;
89  RPDetCoordinateAlgebraObjs *alg;
90  };
91 
92  vector<HitWithAlg> applicable_hits;
93 
94  for (auto &ds : hits)
95  {
96  unsigned int detId = ds.detId();
97 
98  for (auto &h : ds)
99  {
100  RPDetCoordinateAlgebraObjs * alg = getDetAlgebraData(detId, tot_geom);
101  if (alg->available_)
102  applicable_hits.push_back({ detId, &h, alg});
103  }
104  }
105 
106  if (applicable_hits.size() < 5)
107  return false;
108 
109  TMatrixD H(applicable_hits.size(), 4);
110  TVectorD V(applicable_hits.size());
111  TVectorD V_inv(applicable_hits.size());
112  TVectorD U(applicable_hits.size());
113 
114  for(unsigned int i = 0; i < applicable_hits.size(); ++i)
115  {
116  RPDetCoordinateAlgebraObjs *alg_obj = applicable_hits[i].alg;
117 
118  H(i,0) = alg_obj->readout_direction_.X();
119  H(i,1) = alg_obj->readout_direction_.Y();
120  double delta_z = alg_obj->centre_of_det_global_position_.Z()-z_0;
121  H(i,2) = alg_obj->readout_direction_.X()*delta_z;
122  H(i,3) = alg_obj->readout_direction_.Y()*delta_z;
123  double var = applicable_hits[i].hit->getSigma();
124  var*=var;
125  V[i] = var;
126  V_inv[i] = 1.0/var;
127  U[i] = applicable_hits[i].hit->getPosition() - alg_obj->rec_u_0_;
128  }
129 
130  TMatrixD H_T_V_inv(TMatrixD::kTransposed, H);
131  multiplyByDiagonalInPlace(H_T_V_inv, V_inv);
132  TMatrixD V_a(H_T_V_inv);
133  TMatrixD V_a_mult(V_a, TMatrixD::kMult, H);
134  try
135  {
136  V_a_mult.Invert();
137  }
138  catch (cms::Exception &e)
139  {
140  LogError("TotemRPLocalTrackFitterAlgorithm") << "Error in TotemRPLocalTrackFitterAlgorithm::fitTrack > "
141  << "Fit matrix is singular. Skipping.";
142  return false;
143  }
144 
145  TMatrixD u_to_a(V_a_mult, TMatrixD::kMult, H_T_V_inv);
146  TVectorD a(U);
147  a *= u_to_a;
148 
149  fitted_track.setZ0(z_0);
150  fitted_track.setParameterVector(a);
151  fitted_track.setCovarianceMatrix(V_a_mult);
152 
153  double Chi_2 = 0;
154  for(unsigned int i=0; i<applicable_hits.size(); ++i)
155  {
156  RPDetCoordinateAlgebraObjs *alg_obj = applicable_hits[i].alg;
157  TVector2 readout_dir = alg_obj->readout_direction_;
158  double det_z = alg_obj->centre_of_det_global_position_.Z();
159  double sigma_str = applicable_hits[i].hit->getSigma();
160  double sigma_str_2 = sigma_str*sigma_str;
161  TVector2 fited_det_xy_point = fitted_track.getTrackPoint(det_z);
162  double U_readout = applicable_hits[i].hit->getPosition() - alg_obj->rec_u_0_;
163  double U_fited = (readout_dir*=fited_det_xy_point);
164  double residual = U_fited - U_readout;
165  TMatrixD V_T_Cov_X_Y(1,2);
166  V_T_Cov_X_Y(0,0) = readout_dir.X();
167  V_T_Cov_X_Y(0,1) = readout_dir.Y();
168  TMatrixD V_T_Cov_X_Y_mult(V_T_Cov_X_Y, TMatrixD::kMult, fitted_track.trackPointInterpolationCovariance(det_z));
169  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();
170  double pull_normalization = sqrt(sigma_str_2 - fit_strip_var);
171  double pull = residual/pull_normalization;
172 
173  Chi_2 += residual*residual / sigma_str_2;
174 
175  TotemRPLocalTrack::FittedRecHit hit_point(*(applicable_hits[i].hit), TVector3(fited_det_xy_point.X(),
176  fited_det_xy_point.Y(), det_z), residual, pull);
177  fitted_track.addHit(applicable_hits[i].detId, hit_point);
178  }
179 
180  fitted_track.setChiSquared(Chi_2);
181  fitted_track.setValid(true);
182  return true;
183 }
void setChiSquared(double &chiSquared)
int i
Definition: DBlmapReader.cc:9
void multiplyByDiagonalInPlace(TMatrixD &mt, const TVectorD &diag)
A matrix multiplication shorthand.
TVector2 getTrackPoint(double z) const
returns (x, y) vector
void addHit(unsigned int detId, const FittedRecHit &hit)
FWCore Framework interface EventSetupRecordImplementation h
Helper function to determine trigger accepts.
TMatrixD trackPointInterpolationCovariance(double z) const
void setZ0(double z0)
T sqrt(T t)
Definition: SSEVec.h:18
Reconstructed hit in TOTEM RP.
Definition: TotemRPRecHit.h:18
void setValid(bool valid)
void setCovarianceMatrix(const TMatrixD &par_covariance_matrix)
double a
Definition: hdecay.h:121
void setParameterVector(const TVectorD &track_params_vector)
RPDetCoordinateAlgebraObjs * getDetAlgebraData(unsigned int det_id, const TotemRPGeometry &tot_rp_geom)
TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs * TotemRPLocalTrackFitterAlgorithm::getDetAlgebraData ( unsigned int  det_id,
const TotemRPGeometry 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 65 of file TotemRPLocalTrackFitterAlgorithm.cc.

66 {
67  auto it = det_data_map_.find(det_id);
68  if (it != det_data_map_.end())
69  {
70  return &(it->second);
71  } else {
72  det_data_map_[det_id] = prepareReconstAlgebraData(det_id, tot_rp_geom);
73  return &det_data_map_[det_id];
74  }
75 }
std::unordered_map< unsigned int, RPDetCoordinateAlgebraObjs > det_data_map_
A cache of reconstruction data. Must be reset every time the geometry chagnges.
RPDetCoordinateAlgebraObjs prepareReconstAlgebraData(unsigned int det_id, const TotemRPGeometry &tot_rp_geom)
Build the reconstruction data.
void TotemRPLocalTrackFitterAlgorithm::multiplyByDiagonalInPlace ( TMatrixD &  mt,
const TVectorD &  diag 
)
private

A matrix multiplication shorthand.

Definition at line 187 of file TotemRPLocalTrackFitterAlgorithm.cc.

References i, and j.

188 {
189  for(int i=0; i<mt.GetNrows(); ++i)
190  {
191  for(int j=0; j<mt.GetNcols(); ++j)
192  {
193  mt[i][j] *= diag[j];
194  }
195  }
196 }
int i
Definition: DBlmapReader.cc:9
int j
Definition: DBlmapReader.cc:9
TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs TotemRPLocalTrackFitterAlgorithm::prepareReconstAlgebraData ( unsigned int  det_id,
const TotemRPGeometry tot_rp_geom 
)
private

Build the reconstruction data.

Definition at line 38 of file TotemRPLocalTrackFitterAlgorithm.cc.

References TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::available_, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::centre_of_det_global_position_, TotemRPGeometry::GetDetTranslation(), TotemRPGeometry::LocalToGlobalDirection(), TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::readout_direction_, TotemRPLocalTrackFitterAlgorithm::RPDetCoordinateAlgebraObjs::rec_u_0_, and findQualityFiles::v.

39 {
40  RPDetCoordinateAlgebraObjs det_algebra_obj;
41 
42  det_algebra_obj.centre_of_det_global_position_ = convert3vector(tot_rp_geom.GetDetTranslation(det_id));
43 
44  HepMC::ThreeVector rp_topology_stripaxis = rp_topology_.GetStripReadoutAxisDir();
45  CLHEP::Hep3Vector rp_topology_stripaxis_clhep;
46 
47  rp_topology_stripaxis_clhep.setX(rp_topology_stripaxis.x());
48  rp_topology_stripaxis_clhep.setY(rp_topology_stripaxis.y());
49  rp_topology_stripaxis_clhep.setZ(rp_topology_stripaxis.z());
50 
51  TVector3 rd_dir = convert3vector( tot_rp_geom.LocalToGlobalDirection(det_id, rp_topology_stripaxis_clhep ) );
52 
53  TVector2 v(rd_dir.X(), rd_dir.Y());
54  det_algebra_obj.readout_direction_ = v.Unit();
55  det_algebra_obj.rec_u_0_ = 0.0;
56  det_algebra_obj.available_ = true;
57  det_algebra_obj.rec_u_0_ = - (det_algebra_obj.readout_direction_ * det_algebra_obj.centre_of_det_global_position_.XYvector());
58 
59  return det_algebra_obj;
60 }
const HepMC::ThreeVector & GetStripReadoutAxisDir() const
Definition: RPTopology.h:33
CLHEP::Hep3Vector GetDetTranslation(unsigned int id) const
CLHEP::Hep3Vector LocalToGlobalDirection(unsigned int id, const CLHEP::Hep3Vector dir) const
static TVector3 convert3vector(const CLHEP::Hep3Vector &v)
void TotemRPLocalTrackFitterAlgorithm::reset ( void  )

Resets the reconstruction-data cache.

Definition at line 30 of file TotemRPLocalTrackFitterAlgorithm.cc.

Referenced by TotemRPLocalTrackFitter::produce().

31 {
32  det_data_map_.clear();
33 }
std::unordered_map< unsigned int, RPDetCoordinateAlgebraObjs > det_data_map_
A cache of reconstruction data. Must be reset every time the geometry chagnges.

Member Data Documentation

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 53 of file TotemRPLocalTrackFitterAlgorithm.h.

RPTopology TotemRPLocalTrackFitterAlgorithm::rp_topology_
private

Definition at line 55 of file TotemRPLocalTrackFitterAlgorithm.h.