CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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, HitVector
HitVectorMap
 
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

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.

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.

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

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  }
unsigned int numberOfParameters() const
reco::FormulaEvaluator pixelEfficiencyFunction_
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
template<typename TRACK_TYPE, typename HIT_TYPE>
virtual CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::~CTPPSTimingTrackRecognition ( )
virtualdefault

Member Function Documentation

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.

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.
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 198 of file CTPPSTimingTrackRecognition.h.

References SiStripPI::max, SiStripPI::min, mps_fire::result, hit::x, CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::SpatialRange::xBegin, CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::SpatialRange::xEnd, hit::y, CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::SpatialRange::yBegin, CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::SpatialRange::yEnd, hit::z, CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::SpatialRange::zBegin, and CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::SpatialRange::zEnd.

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

References mps_fire::i, submitPVValidationJobs::params, CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::DimensionParameters::rangeBegin, and CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::DimensionParameters::rangeEnd.

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

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 229 of file CTPPSTimingTrackRecognition.h.

References validate-o2o-wbm::f, funct::pow(), mathSSE::sqrt(), and histoStyle::weight.

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

Member Data Documentation

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().

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

Definition at line 60 of file CTPPSTimingTrackRecognition.h.

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

Definition at line 61 of file CTPPSTimingTrackRecognition.h.

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

Definition at line 58 of file CTPPSTimingTrackRecognition.h.

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

Definition at line 59 of file CTPPSTimingTrackRecognition.h.

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

Definition at line 62 of file CTPPSTimingTrackRecognition.h.