CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
ClusterShapeTrajectoryFilter Class Reference

#include <ClusterShapeTrajectoryFilter.h>

Inheritance diagram for ClusterShapeTrajectoryFilter:
TrajectoryFilter

Public Member Functions

 ClusterShapeTrajectoryFilter (const ClusterShapeHitFilter *f)
 
virtual std::string name () const
 
virtual bool qualityFilter (const TempTrajectory &) const
 
virtual bool qualityFilter (const Trajectory &) const
 
virtual bool toBeContinued (TempTrajectory &) const
 
virtual bool toBeContinued (Trajectory &) const
 
virtual ~ClusterShapeTrajectoryFilter ()
 
- Public Member Functions inherited from TrajectoryFilter
virtual bool operator() (TempTrajectory &t) const
 
virtual bool operator() (Trajectory &t) const
 
virtual ~TrajectoryFilter ()
 

Private Attributes

const ClusterShapeHitFiltertheFilter
 

Additional Inherited Members

- Public Types inherited from TrajectoryFilter
typedef CkfComponentsRecord Record
 
- Static Public Attributes inherited from TrajectoryFilter
static const bool qualityFilterIfNotContributing =true
 
static const bool toBeContinuedIfNotContributing =true
 

Detailed Description

Definition at line 16 of file ClusterShapeTrajectoryFilter.h.

Constructor & Destructor Documentation

ClusterShapeTrajectoryFilter::ClusterShapeTrajectoryFilter ( const ClusterShapeHitFilter f)
inline

Definition at line 21 of file ClusterShapeTrajectoryFilter.h.

21 :theFilter(f){}
const ClusterShapeHitFilter * theFilter
ClusterShapeTrajectoryFilter::~ClusterShapeTrajectoryFilter ( )
virtual

Definition at line 35 of file ClusterShapeTrajectoryFilter.cc.

36 {
37 }

Member Function Documentation

virtual std::string ClusterShapeTrajectoryFilter::name ( void  ) const
inlinevirtual

Implements TrajectoryFilter.

Definition at line 31 of file ClusterShapeTrajectoryFilter.h.

31 { return "ClusterShapeTrajectoryFilter"; }
bool ClusterShapeTrajectoryFilter::qualityFilter ( const TempTrajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 206 of file ClusterShapeTrajectoryFilter.cc.

References TempTrajectory::measurements(), cmsutils::bqueue< T >::size(), and lumiQTWidget::t.

207 {
208  TempTrajectory t = trajectory;
209 
210  // Check if ok
211  if(toBeContinued(t)) return true;
212 
213  // Should take out last
214  if(t.measurements().size() <= 3) return false;
215 
216  return true;
217 }
const DataContainer & measurements() const
virtual bool toBeContinued(TempTrajectory &) const
size_type size() const
Definition: bqueue.h:146
bool ClusterShapeTrajectoryFilter::qualityFilter ( const Trajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 199 of file ClusterShapeTrajectoryFilter.cc.

200 {
201  return true;
202 }
bool ClusterShapeTrajectoryFilter::toBeContinued ( TempTrajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 107 of file ClusterShapeTrajectoryFilter.cc.

References TransientTrackingRecHit::det(), TrajectoryStateOnSurface::globalDirection(), TransientTrackingRecHit::hit(), TrackingRecHit::isValid(), LogTrace, TempTrajectory::measurements(), SiStripMatchedRecHit2D::monoHit(), ProjectedSiStripRecHit2D::originalHit(), GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, cmsutils::bqueue< T >::rbegin(), cmsutils::bqueue< T >::rend(), SiStripMatchedRecHit2D::stereoHit(), and GeomDet::subDetector().

108 {
109  TempTrajectory::DataContainer tms = trajectory.measurements();
110 
112  tm = tms.rbegin(); tm!= tms.rend(); --tm)
113  {
114  const TransientTrackingRecHit* ttRecHit = &(*((*tm).recHit()));
115 
116  if(ttRecHit->isValid())
117  {
118  const TrackingRecHit * tRecHit = ttRecHit->hit();
119 
120  TrajectoryStateOnSurface ts = (*tm).updatedState();
121  GlobalVector gdir = ts.globalDirection();
122 
123  if(ttRecHit->det()->subDetector() == GeomDetEnumerators::PixelBarrel ||
125  { // pixel
126  const SiPixelRecHit* recHit =
127  dynamic_cast<const SiPixelRecHit *>(tRecHit);
128 
129  if(recHit != 0)
130  if(! theFilter->isCompatible(*recHit, gdir))
131  {
132  LogTrace("TrajectFilter")
133  << " [TrajectFilter] fail pixel";
134  return false;
135  }
136  }
137  else
138  { // strip
139  if(dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit) != 0)
140  { // glued
141  const SiStripMatchedRecHit2D* recHit =
142  dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit);
143 
144  if(recHit != 0)
145  {
146  if(! theFilter->isCompatible(recHit->monoHit(), gdir))
147  {
148  LogTrace("TrajectFilter")
149  << " [TrajectFilter] fail strip matched 1st";
150  return false;
151  }
152 
153  if(! theFilter->isCompatible(recHit->stereoHit(), gdir))
154  {
155  LogTrace("TrajectFilter")
156  << " [TrajectFilter] fail strip matched 2nd";
157  return false;
158  }
159  }
160  }
161  else
162  { // single
163  if(dynamic_cast<const SiStripRecHit2D *>(tRecHit) != 0)
164  { // normal
165  const SiStripRecHit2D* recHit =
166  dynamic_cast<const SiStripRecHit2D *>(tRecHit);
167 
168  if(recHit != 0)
169  if(! theFilter->isCompatible(*recHit, gdir))
170  {
171  LogTrace("TrajectFilter")
172  << " [TrajectFilter] fail strip single";
173  return false;
174  }
175  }
176  else
177  { // projected
178  const ProjectedSiStripRecHit2D* recHit =
179  dynamic_cast<const ProjectedSiStripRecHit2D *>(tRecHit);
180 
181  if(recHit != 0)
182  if(! theFilter->isCompatible(recHit->originalHit(), gdir))
183  {
184  LogTrace("TrajectFilter")
185  << " [TrajectFilter] fail strip projected";
186  return false;
187  }
188  }
189  }
190  }
191  }
192  }
193 
194  return true;
195 }
const_iterator rend() const
Definition: bqueue.h:145
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, PixelData const *pd=0) const
const DataContainer & measurements() const
virtual SubDetector subDetector() const =0
Which subdetector.
virtual const TrackingRecHit * hit() const =0
const ClusterShapeHitFilter * theFilter
#define LogTrace(id)
SiStripRecHit2D stereoHit() const
iterator rbegin()
Definition: bqueue.h:143
bool isValid() const
SiStripRecHit2D monoHit() const
const GeomDet * det() const
The GomeDet* can be zero for InvalidTransientRecHits and for TConstraintRecHit2Ds.
const SiStripRecHit2D & originalHit() const
GlobalVector globalDirection() const
Our base class.
Definition: SiPixelRecHit.h:22
bool ClusterShapeTrajectoryFilter::toBeContinued ( Trajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 41 of file ClusterShapeTrajectoryFilter.cc.

References TransientTrackingRecHit::det(), TrajectoryStateOnSurface::globalDirection(), TransientTrackingRecHit::hit(), TrackingRecHit::isValid(), Trajectory::measurements(), SiStripMatchedRecHit2D::monoHit(), ProjectedSiStripRecHit2D::originalHit(), GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, SiStripMatchedRecHit2D::stereoHit(), and GeomDet::subDetector().

42 {
43  vector<TrajectoryMeasurement> tms = trajectory.measurements();
44 
45  for(vector<TrajectoryMeasurement>::const_iterator
46  tm = tms.begin(); tm!= tms.end(); tm++)
47  {
48  const TransientTrackingRecHit* ttRecHit = &(*((*tm).recHit()));
49 
50  if(ttRecHit->isValid())
51  {
52  const TrackingRecHit * tRecHit = ttRecHit->hit();
53 
54  TrajectoryStateOnSurface ts = (*tm).updatedState();
55  const GlobalVector gdir = ts.globalDirection();
56 
57  if(ttRecHit->det()->subDetector() == GeomDetEnumerators::PixelBarrel ||
59  { // pixel
60  const SiPixelRecHit* recHit =
61  dynamic_cast<const SiPixelRecHit *>(tRecHit);
62 
63  if(recHit != 0)
64  return theFilter->isCompatible(*recHit, gdir);
65  }
66  else
67  { // strip
68  if(dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit) != 0)
69  { // glued
70  const SiStripMatchedRecHit2D* recHit =
71  dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit);
72 
73  if(recHit != 0)
74  {
75  return (theFilter->isCompatible(recHit->monoHit() , gdir) &&
76  theFilter->isCompatible(recHit->stereoHit(), gdir));
77  }
78  }
79  else
80  { // single
81  if(dynamic_cast<const SiStripRecHit2D *>(tRecHit) != 0)
82  { // normal
83  const SiStripRecHit2D* recHit =
84  dynamic_cast<const SiStripRecHit2D *>(tRecHit);
85 
86  if(recHit != 0)
87  return theFilter->isCompatible(*recHit, gdir);
88  }
89  else
90  { // projected
91  const ProjectedSiStripRecHit2D* recHit =
92  dynamic_cast<const ProjectedSiStripRecHit2D *>(tRecHit);
93 
94  if(recHit != 0)
95  return theFilter->isCompatible(recHit->originalHit(), gdir);
96  }
97  }
98  }
99  }
100  }
101 
102  return true;
103 }
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, PixelData const *pd=0) const
virtual SubDetector subDetector() const =0
Which subdetector.
virtual const TrackingRecHit * hit() const =0
const ClusterShapeHitFilter * theFilter
DataContainer const & measurements() const
Definition: Trajectory.h:203
SiStripRecHit2D stereoHit() const
bool isValid() const
SiStripRecHit2D monoHit() const
const GeomDet * det() const
The GomeDet* can be zero for InvalidTransientRecHits and for TConstraintRecHit2Ds.
const SiStripRecHit2D & originalHit() const
GlobalVector globalDirection() const
Our base class.
Definition: SiPixelRecHit.h:22

Member Data Documentation

const ClusterShapeHitFilter* ClusterShapeTrajectoryFilter::theFilter
private

Definition at line 35 of file ClusterShapeTrajectoryFilter.h.