CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
CTPPSDiamondTrackRecognition Class Reference

Class performing smart reconstruction for PPS Diamond Detectors. More...

#include <CTPPSDiamondTrackRecognition.h>

Inheritance diagram for CTPPSDiamondTrackRecognition:
CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >

Public Member Functions

void addHit (const CTPPSDiamondRecHit &recHit) override
 Feed a new hit to the tracks recognition algorithm. More...
 
void clear () override
 Reset internal state of a class instance. More...
 
 CTPPSDiamondTrackRecognition (const edm::ParameterSet &iConfig)
 
int produceTracks (edm::DetSet< CTPPSDiamondLocalTrack > &tracks) override
 Produce a collection of tracks for the current station, given its hits collection. More...
 
- Public Member Functions inherited from CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >
 CTPPSTimingTrackRecognition (const edm::ParameterSet &iConfig)
 
virtual ~CTPPSTimingTrackRecognition ()=default
 

Private Attributes

bool excludeSingleEdgeHits_
 
std::unordered_map< int, int > mhMap_
 

Additional Inherited Members

- Protected Types inherited from CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >
typedef std::vector< CTPPSDiamondRecHitHitVector
 
typedef std::unordered_map< int, HitVectorHitVectorMap
 
typedef std::vector< CTPPSDiamondLocalTrackTrackVector
 
- Protected Member Functions inherited from CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >
SpatialRange getHitSpatialRange (const HitVector &hits)
 
void producePartialTracks (const HitVector &hits, const DimensionParameters &param, float(*getHitCenter)(const CTPPSDiamondRecHit &), float(*getHitRangeWidth)(const CTPPSDiamondRecHit &), void(*setTrackCenter)(CTPPSDiamondLocalTrack &, float), void(*setTrackSigma)(CTPPSDiamondLocalTrack &, float), TrackVector &result)
 
bool timeEval (const HitVector &hits, float &meanTime, float &timeSigma) const
 
- Protected Attributes inherited from CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >
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

Class performing smart reconstruction for PPS Diamond Detectors.

Date
Jan 2017

Definition at line 26 of file CTPPSDiamondTrackRecognition.h.

Constructor & Destructor Documentation

CTPPSDiamondTrackRecognition::CTPPSDiamondTrackRecognition ( const edm::ParameterSet iConfig)

Member Function Documentation

void CTPPSDiamondTrackRecognition::addHit ( const CTPPSDiamondRecHit recHit)
overridevirtual

Feed a new hit to the tracks recognition algorithm.

Implements CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >.

Definition at line 32 of file CTPPSDiamondTrackRecognition.cc.

References excludeSingleEdgeHits_, CTPPSDiamondRecHit::getOOTIndex(), CTPPSDiamondRecHit::getToT(), and CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >::hitVectorMap_.

Referenced by CTPPSDiamondLocalTrackFitter::produce().

33 {
34  if ( excludeSingleEdgeHits_ && recHit.getToT() <= 0. )
35  return;
36  // store hit parameters
37  hitVectorMap_[recHit.getOOTIndex()].emplace_back( recHit );
38 }
int getOOTIndex() const
HitVectorMap hitVectorMap_
RecHit vectors that should be processed separately while reconstructing tracks.
float getToT() const
void CTPPSDiamondTrackRecognition::clear ( void  )
overridevirtual

Reset internal state of a class instance.

Reimplemented from CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >.

Definition at line 23 of file CTPPSDiamondTrackRecognition.cc.

References CTPPSTimingTrackRecognition< TRACK_TYPE, HIT_TYPE >::clear(), and mhMap_.

Referenced by CTPPSDiamondLocalTrackFitter::produce().

24 {
26  mhMap_.clear();
27 }
virtual void clear()
Reset internal state of a class instance.
std::unordered_map< int, int > mhMap_
int CTPPSDiamondTrackRecognition::produceTracks ( edm::DetSet< CTPPSDiamondLocalTrack > &  tracks)
overridevirtual

Produce a collection of tracks for the current station, given its hits collection.

Implements CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >.

Definition at line 43 of file CTPPSDiamondTrackRecognition.cc.

References CTPPSDiamondLocalTrack::containsHit(), excludeSingleEdgeHits_, f, CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >::getHitSpatialRange(), hfClusterShapes_cfi::hits, CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >::hitVectorMap_, mhMap_, position, CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >::producePartialTracks(), edm::DetSet< T >::push_back(), CTPPSTimingLocalTrack::setPosition(), CTPPSTimingLocalTrack::setPositionSigma(), CTPPSTimingLocalTrack::setT(), CTPPSTimingLocalTrack::setTSigma(), CTPPSTimingLocalTrack::setValid(), CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >::timeEval(), CTPPSTimingTrackRecognition< CTPPSDiamondLocalTrack, CTPPSDiamondRecHit >::tolerance_, HiIsolationCommonParameters_cff::track, and x.

Referenced by CTPPSDiamondLocalTrackFitter::produce().

44 {
45  int numberOfTracks = 0;
46  DimensionParameters param;
47 
48  auto getX = []( const CTPPSDiamondRecHit& hit ){ return hit.getX(); };
49  auto getXWidth = []( const CTPPSDiamondRecHit& hit ){ return hit.getXWidth(); };
50  auto setX = []( CTPPSDiamondLocalTrack& track, float x ){ track.setPosition( math::XYZPoint( x, 0., 0. ) ); };
51  auto setXSigma = []( CTPPSDiamondLocalTrack& track, float sigma ){ track.setPositionSigma( math::XYZPoint( sigma, 0., 0. ) ); };
52 
53  for ( const auto& hitBatch: hitVectorMap_ ) {
54  // separate the tracking for each bunch crossing
55  const auto& oot = hitBatch.first;
56  const auto& hits = hitBatch.second;
57 
58  auto hitRange = getHitSpatialRange( hits );
59 
60  TrackVector xPartTracks;
61 
62  // produce tracks in x dimension
63  param.rangeBegin = hitRange.xBegin;
64  param.rangeEnd = hitRange.xEnd;
65  producePartialTracks( hits, param, getX, getXWidth, setX, setXSigma, xPartTracks );
66 
67  if ( xPartTracks.empty() )
68  continue;
69 
70  const float yRangeCenter = 0.5f*( hitRange.yBegin + hitRange.yEnd );
71  const float zRangeCenter = 0.5f*( hitRange.zBegin + hitRange.zEnd );
72  const float ySigma = 0.5f*( hitRange.yEnd - hitRange.yBegin );
73  const float zSigma = 0.5f*( hitRange.zEnd - hitRange.zBegin );
74 
75  for ( const auto& xTrack: xPartTracks ) {
76  math::XYZPoint position( xTrack.getX0(), yRangeCenter, zRangeCenter );
77  math::XYZPoint positionSigma( xTrack.getX0Sigma(), ySigma, zSigma );
78 
79  const int multipleHits = ( mhMap_.find(oot) != mhMap_.end() )
80  ? mhMap_[oot]
81  : 0;
82  CTPPSDiamondLocalTrack newTrack( position, positionSigma, 0.f, 0.f, oot, multipleHits );
83 
84  // find contributing hits
85  HitVector componentHits;
86  for ( const auto& hit : hits )
87  if ( newTrack.containsHit( hit, tolerance_ ) && ( !excludeSingleEdgeHits_ || hit.getToT() > 0. ) )
88  componentHits.emplace_back( hit );
89  // compute timing information
90  float mean_time = 0.f, time_sigma = 0.f;
91  bool valid_hits = timeEval( componentHits, mean_time, time_sigma );
92  newTrack.setValid( valid_hits );
93  newTrack.setT( mean_time );
94  newTrack.setTSigma( time_sigma );
95 
96  tracks.push_back( newTrack );
97  }
98  }
99 
100  return numberOfTracks;
101 }
void push_back(const T &t)
Definition: DetSet.h:68
void producePartialTracks(const HitVector &hits, const DimensionParameters &param, float(*getHitCenter)(const CTPPSDiamondRecHit &), float(*getHitRangeWidth)(const CTPPSDiamondRecHit &), void(*setTrackCenter)(CTPPSDiamondLocalTrack &, float), void(*setTrackSigma)(CTPPSDiamondLocalTrack &, float), TrackVector &result)
Reconstructed hit in diamond detectors.
bool timeEval(const HitVector &hits, float &meanTime, float &timeSigma) const
double f[11][100]
void setPosition(const math::XYZPoint &pos0)
HitVectorMap hitVectorMap_
RecHit vectors that should be processed separately while reconstructing tracks.
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
static int position[264][3]
Definition: ReadPGInfo.cc:509
std::unordered_map< int, int > mhMap_
void setPositionSigma(const math::XYZPoint &pos0_sigma)

Member Data Documentation

bool CTPPSDiamondTrackRecognition::excludeSingleEdgeHits_
private

Definition at line 39 of file CTPPSDiamondTrackRecognition.h.

Referenced by addHit(), and produceTracks().

std::unordered_map<int,int> CTPPSDiamondTrackRecognition::mhMap_
private

Definition at line 38 of file CTPPSDiamondTrackRecognition.h.

Referenced by clear(), and produceTracks().