CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes
CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE > Class Template Referenceabstract

#include <CTPPSTimingTrackRecognition.h>

Classes

struct  DimensionParameters
 
struct  SpatialRange
 Structure representing a 3D range in space. More...
 

Public Member Functions

virtual void addHit (const HIT_TYPE &recHit)=0
 Add new hit to the set from which the tracks are reconstructed. More...
 
virtual void clear ()
 Reset internal state of a class instance. More...
 
 CTPPSTimingTrackRecognition (const edm::ParameterSet &iConfig)
 
virtual int produceTracks (edm::DetSet< TRACK_TYPE > &tracks)=0
 Produce a collection of tracks, given its hits collection. More...
 
virtual ~CTPPSTimingTrackRecognition ()=default
 

Protected Types

typedef std::vector< HIT_TYPE > HitVector
 
typedef std::unordered_map< int, HitVectorHitVectorMap
 
typedef std::vector< TRACK_TYPE > TrackVector
 

Protected Member Functions

SpatialRange getHitSpatialRange (const HitVector &hits)
 
void producePartialTracks (const HitVector &hits, const DimensionParameters &param, float(*getHitCenter)(const HIT_TYPE &), float(*getHitRangeWidth)(const HIT_TYPE &), void(*setTrackCenter)(TRACK_TYPE &, float), void(*setTrackSigma)(TRACK_TYPE &, float), TrackVector &result)
 
bool timeEval (const HitVector &hits, float &meanTime, float &timeSigma) const
 

Protected Attributes

HitVectorMap hitVectorMap_
 RecHit vectors that should be processed separately while reconstructing tracks. More...
 
reco::FormulaEvaluator pixelEfficiencyFunction_
 
const float resolution_
 
const float sigma_
 
const float threshold_
 
const float thresholdFromMaximum_
 
const float tolerance_
 

Detailed Description

template<typename TRACK_TYPE, typename HIT_TYPE>
class CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >

Class intended to perform general CTPPS timing detectors track recognition, as well as construction of specialized classes (for now CTPPSDiamond and TotemTiming local tracks).

Definition at line 31 of file CTPPSTimingTrackRecognition.h.

Member Typedef Documentation

◆ HitVector

template<typename TRACK_TYPE, typename HIT_TYPE>
typedef std::vector<HIT_TYPE> CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::HitVector
protected

Definition at line 66 of file CTPPSTimingTrackRecognition.h.

◆ HitVectorMap

template<typename TRACK_TYPE, typename HIT_TYPE>
typedef std::unordered_map<int, HitVector> CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::HitVectorMap
protected

Definition at line 67 of file CTPPSTimingTrackRecognition.h.

◆ TrackVector

template<typename TRACK_TYPE, typename HIT_TYPE>
typedef std::vector<TRACK_TYPE> CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::TrackVector
protected

Definition at line 65 of file CTPPSTimingTrackRecognition.h.

Constructor & Destructor Documentation

◆ CTPPSTimingTrackRecognition()

template<typename TRACK_TYPE, typename HIT_TYPE>
CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::CTPPSTimingTrackRecognition ( const edm::ParameterSet iConfig)
inline

Definition at line 33 of file CTPPSTimingTrackRecognition.h.

34  : threshold_(iConfig.getParameter<double>("threshold")),
35  thresholdFromMaximum_(iConfig.getParameter<double>("thresholdFromMaximum")),
36  resolution_(iConfig.getParameter<double>("resolution")),
37  sigma_(iConfig.getParameter<double>("sigma")),
38  tolerance_(iConfig.getParameter<double>("tolerance")),
39  pixelEfficiencyFunction_(iConfig.getParameter<std::string>("pixelEfficiencyFunction")) {
41  throw cms::Exception("CTPPSTimingTrackRecognition")
42  << "Invalid number of parameters to the pixel efficiency function! "
44  }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
unsigned int numberOfParameters() const
reco::FormulaEvaluator pixelEfficiencyFunction_

◆ ~CTPPSTimingTrackRecognition()

template<typename TRACK_TYPE, typename HIT_TYPE>
virtual CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::~CTPPSTimingTrackRecognition ( )
virtualdefault

Member Function Documentation

◆ addHit()

template<typename TRACK_TYPE, typename HIT_TYPE>
virtual void CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::addHit ( const HIT_TYPE &  recHit)
pure virtual

Add new hit to the set from which the tracks are reconstructed.

Implemented in CTPPSDiamondTrackRecognition, and TotemTimingTrackRecognition.

◆ clear()

template<typename TRACK_TYPE, typename HIT_TYPE>
virtual void CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::clear ( void  )
inlinevirtual

Reset internal state of a class instance.

Reimplemented in CTPPSDiamondTrackRecognition.

Definition at line 50 of file CTPPSTimingTrackRecognition.h.

Referenced by CTPPSDiamondTrackRecognition::clear().

50 { hitVectorMap_.clear(); }
HitVectorMap hitVectorMap_
RecHit vectors that should be processed separately while reconstructing tracks.

◆ getHitSpatialRange()

template<class TRACK_TYPE , class HIT_TYPE >
CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::SpatialRange CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::getHitSpatialRange ( const HitVector hits)
inlineprotected

Retrieve the bounds of a 3D range in which all hits from given collection are contained.

Parameters
[in]hitshits collection to retrieve the range from

Definition at line 196 of file CTPPSTimingTrackRecognition.h.

196  {
197  bool initialized = false;
198  SpatialRange result = {};
199 
200  for (const auto& hit : hits) {
201  const float xBegin = hit.x() - 0.5f * hit.xWidth(), xEnd = hit.x() + 0.5f * hit.xWidth();
202  const float yBegin = hit.y() - 0.5f * hit.yWidth(), yEnd = hit.y() + 0.5f * hit.yWidth();
203  const float zBegin = hit.z() - 0.5f * hit.zWidth(), zEnd = hit.z() + 0.5f * hit.zWidth();
204 
205  if (!initialized) {
206  result.xBegin = xBegin;
207  result.xEnd = xEnd;
208  result.yBegin = yBegin;
209  result.yEnd = yEnd;
210  result.zBegin = zBegin;
211  result.zEnd = zEnd;
212  initialized = true;
213  continue;
214  }
215  result.xBegin = std::min(result.xBegin, xBegin);
216  result.xEnd = std::max(result.xEnd, xEnd);
217  result.yBegin = std::min(result.yBegin, yBegin);
218  result.yEnd = std::max(result.yEnd, yEnd);
219  result.zBegin = std::min(result.zBegin, zBegin);
220  result.zEnd = std::max(result.zEnd, zEnd);
221  }
222 
223  return result;
224 }

◆ producePartialTracks()

template<class TRACK_TYPE, class HIT_TYPE>
void CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::producePartialTracks ( const HitVector hits,
const DimensionParameters param,
float(*)(const HIT_TYPE &)  getHitCenter,
float(*)(const HIT_TYPE &)  getHitRangeWidth,
void(*)(TRACK_TYPE &, float)  setTrackCenter,
void(*)(TRACK_TYPE &, float)  setTrackSigma,
TrackVector result 
)
inlineprotected

Produce all partial tracks from given set with regard to single dimension.

Parameters
[in]hitsvector of hits from which the tracks are created
[in]paramdescribe all parameters used by 1D track recognition algorithm
[in]getHitCenterfunction extracting hit's center in the dimension that the partial tracks relate to
[in]getHitRangeWidthanalogue to getHitCenter, but extracts hit's width in specific dimension
[in]setTrackCenterfunction used to set track's position in considered dimension
[in]setTrackSigmafunction used to set track's sigma in considered dimension
[out]resultvector to which produced tracks are appended

Definition at line 119 of file CTPPSTimingTrackRecognition.h.

126  {
127  const float invResolution = 1. / resolution_;
128  const float profileRangeMargin = sigma_ * 3.;
129  const float profileRangeBegin = param.rangeBegin - profileRangeMargin;
130  const float profileRangeEnd = param.rangeEnd + profileRangeMargin;
131 
132  std::vector<float> hitProfile((profileRangeEnd - profileRangeBegin) * invResolution + 1, 0.);
133  // extra component to make sure that the profile drops below the threshold at range's end
134  *hitProfile.rbegin() = -1.f;
135 
136  // Creates hit profile
137  for (auto const& hit : hits) {
138  const float center = getHitCenter(hit), rangeWidth = getHitRangeWidth(hit);
139  std::vector<double> params{center, rangeWidth, sigma_};
140  for (unsigned int i = 0; i < hitProfile.size(); ++i)
141  hitProfile[i] +=
142  pixelEfficiencyFunction_.evaluate(std::vector<double>{profileRangeBegin + i * resolution_}, params);
143  }
144 
145  bool underThreshold = true;
146  float rangeMaximum = -1.0f;
147  bool trackRangeFound = false;
148  int trackRangeBegin = 0, trackRangeEnd = 0;
149 
150  // Searches for tracks in the hit profile
151  for (unsigned int i = 0; i < hitProfile.size(); i++) {
152  if (hitProfile[i] > rangeMaximum)
153  rangeMaximum = hitProfile[i];
154 
155  // Going above the threshold
156  if (underThreshold && hitProfile[i] > threshold_) {
157  underThreshold = false;
158  trackRangeBegin = i;
159  rangeMaximum = hitProfile[i];
160  }
161 
162  // Going under the threshold
163  else if (!underThreshold && hitProfile[i] <= threshold_) {
164  underThreshold = true;
165  trackRangeEnd = i;
166  trackRangeFound = true;
167  }
168 
169  // Finds all tracks within the track range
170  if (trackRangeFound) {
171  float trackThreshold = rangeMaximum - thresholdFromMaximum_;
172  int trackBegin;
173  bool underTrackThreshold = true;
174 
175  for (int j = trackRangeBegin; j <= trackRangeEnd; j++) {
176  if (underTrackThreshold && hitProfile[j] > trackThreshold) {
177  underTrackThreshold = false;
178  trackBegin = j;
179  } else if (!underTrackThreshold && hitProfile[j] <= trackThreshold) {
180  underTrackThreshold = true;
181  TRACK_TYPE track;
182  float leftMargin = profileRangeBegin + resolution_ * trackBegin;
183  float rightMargin = profileRangeBegin + resolution_ * j;
184  setTrackCenter(track, 0.5f * (leftMargin + rightMargin));
185  setTrackSigma(track, 0.5f * (rightMargin - leftMargin));
186  result.push_back(track);
187  }
188  }
189  trackRangeFound = false;
190  }
191  }
192 }
double evaluate(V const &iVariables, P const &iParameters) const
double f[11][100]
reco::FormulaEvaluator pixelEfficiencyFunction_

◆ produceTracks()

template<typename TRACK_TYPE, typename HIT_TYPE>
virtual int CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::produceTracks ( edm::DetSet< TRACK_TYPE > &  tracks)
pure virtual

Produce a collection of tracks, given its hits collection.

Implemented in CTPPSDiamondTrackRecognition, and TotemTimingTrackRecognition.

◆ timeEval()

template<class TRACK_TYPE , class HIT_TYPE >
bool CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::timeEval ( const HitVector hits,
float &  meanTime,
float &  timeSigma 
) const
inlineprotected

Evaluate the time + associated uncertainty for a given track

Note
General remarks:
  • track's time = weighted mean of all hit times with time precision as weight,
  • track's time sigma = uncertainty of the weighted mean
  • hit is ignored if the time precision is equal to 0

Definition at line 227 of file CTPPSTimingTrackRecognition.h.

229  {
230  float mean_num = 0.f, mean_denom = 0.f;
231  bool valid_hits = false;
232  for (const auto& hit : hits) {
233  if (hit.tPrecision() <= 0.)
234  continue;
235  const float weight = std::pow(hit.tPrecision(), -2);
236  mean_num += weight * hit.time();
237  mean_denom += weight;
238  valid_hits = true; // at least one valid hit to account for
239  }
240  mean_time = valid_hits ? (mean_num / mean_denom) : 0.f;
241  time_sigma = valid_hits ? std::sqrt(1.f / mean_denom) : 0.f;
242  return valid_hits;
243 }
Definition: weight.py:1
T sqrt(T t)
Definition: SSEVec.h:23
double f[11][100]
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

Member Data Documentation

◆ hitVectorMap_

template<typename TRACK_TYPE, typename HIT_TYPE>
HitVectorMap CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::hitVectorMap_
protected

RecHit vectors that should be processed separately while reconstructing tracks.

Definition at line 70 of file CTPPSTimingTrackRecognition.h.

Referenced by CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >::clear().

◆ pixelEfficiencyFunction_

template<typename TRACK_TYPE, typename HIT_TYPE>
reco::FormulaEvaluator CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::pixelEfficiencyFunction_
protected

◆ resolution_

template<typename TRACK_TYPE, typename HIT_TYPE>
const float CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::resolution_
protected

Definition at line 60 of file CTPPSTimingTrackRecognition.h.

◆ sigma_

template<typename TRACK_TYPE, typename HIT_TYPE>
const float CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::sigma_
protected

Definition at line 61 of file CTPPSTimingTrackRecognition.h.

◆ threshold_

template<typename TRACK_TYPE, typename HIT_TYPE>
const float CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::threshold_
protected

Definition at line 58 of file CTPPSTimingTrackRecognition.h.

◆ thresholdFromMaximum_

template<typename TRACK_TYPE, typename HIT_TYPE>
const float CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::thresholdFromMaximum_
protected

Definition at line 59 of file CTPPSTimingTrackRecognition.h.

◆ tolerance_

template<typename TRACK_TYPE, typename HIT_TYPE>
const float CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::tolerance_
protected

Definition at line 62 of file CTPPSTimingTrackRecognition.h.