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 edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
 
virtual std::string name () const
 
virtual bool qualityFilter (const TempTrajectory &) const
 
virtual bool qualityFilter (const Trajectory &) const
 
void setEvent (const edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
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 SiPixelClusterShapeCachetheCache
 
edm::EDGetTokenT
< SiPixelClusterShapeCache
theCacheToken
 
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 18 of file ClusterShapeTrajectoryFilter.h.

Constructor & Destructor Documentation

ClusterShapeTrajectoryFilter::ClusterShapeTrajectoryFilter ( const edm::ParameterSet iConfig,
edm::ConsumesCollector iC 
)

Definition at line 39 of file ClusterShapeTrajectoryFilter.cc.

39  :
41  theFilter(nullptr)
42 {}
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::EDGetTokenT< SiPixelClusterShapeCache > theCacheToken
T getParameter(std::string const &) const
const ClusterShapeHitFilter * theFilter
ClusterShapeTrajectoryFilter::~ClusterShapeTrajectoryFilter ( )
virtual

Definition at line 44 of file ClusterShapeTrajectoryFilter.cc.

45 {
46 }

Member Function Documentation

virtual std::string ClusterShapeTrajectoryFilter::name ( void  ) const
inlinevirtual
bool ClusterShapeTrajectoryFilter::qualityFilter ( const TempTrajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 233 of file ClusterShapeTrajectoryFilter.cc.

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

234 {
235  TempTrajectory t = trajectory;
236 
237  // Check if ok
238  if(toBeContinued(t)) return true;
239 
240  // Should take out last
241  if(t.measurements().size() <= 3) return false;
242 
243  return true;
244 }
const DataContainer & measurements() const
virtual bool toBeContinued(TempTrajectory &) const
size_type size() const
Definition: bqueue.h:167
bool ClusterShapeTrajectoryFilter::qualityFilter ( const Trajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 226 of file ClusterShapeTrajectoryFilter.cc.

227 {
228  return true;
229 }
void ClusterShapeTrajectoryFilter::setEvent ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
overridevirtual

Reimplemented from TrajectoryFilter.

Definition at line 48 of file ClusterShapeTrajectoryFilter.cc.

References edm::EventSetup::get(), edm::Event::getByToken(), edm::Handle< T >::product(), edm::ESHandle< class >::product(), theCache, theCacheToken, and theFilter.

48  {
50  iSetup.get<TrajectoryFilter::Record>().get("ClusterShapeHitFilter", shape);
51  theFilter = shape.product();
52 
54  iEvent.getByToken(theCacheToken, cache);
55  theCache = cache.product();
56 }
edm::EDGetTokenT< SiPixelClusterShapeCache > theCacheToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
const ClusterShapeHitFilter * theFilter
const SiPixelClusterShapeCache * theCache
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
bool ClusterShapeTrajectoryFilter::toBeContinued ( TempTrajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 130 of file ClusterShapeTrajectoryFilter.cc.

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

131 {
132  assert(theCache);
133  TempTrajectory::DataContainer tms = trajectory.measurements();
134 
136  tm = tms.rbegin(); tm!= tms.rend(); --tm)
137  {
138  const TrackingRecHit* ttRecHit = &(*((*tm).recHit()));
139 
140  if(ttRecHit->isValid())
141  {
142  const TrackingRecHit * tRecHit = ttRecHit->hit();
143 
144  TrajectoryStateOnSurface ts = (*tm).updatedState();
145  GlobalVector gdir = ts.globalDirection();
146 
152  { // pixel
153  const SiPixelRecHit* recHit =
154  dynamic_cast<const SiPixelRecHit *>(tRecHit);
155 
156  if(recHit != 0)
157  if(! theFilter->isCompatible(*recHit, gdir, *theCache))
158  {
159  LogTrace("TrajectFilter")
160  << " [TrajectFilter] fail pixel";
161  return false;
162  }
163  }
164  else if(GeomDetEnumerators::isTrackerStrip(ttRecHit->det()->subDetector()))
165  { // strip
166  if(dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit) != 0)
167  { // glued
168  const SiStripMatchedRecHit2D* recHit =
169  dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit);
170 
171  if(recHit != 0)
172  {
173  if(! theFilter->isCompatible(recHit->monoHit(), gdir))
174  {
175  LogTrace("TrajectFilter")
176  << " [TrajectFilter] fail strip matched 1st";
177  return false;
178  }
179 
180  if(! theFilter->isCompatible(recHit->stereoHit(), gdir))
181  {
182  LogTrace("TrajectFilter")
183  << " [TrajectFilter] fail strip matched 2nd";
184  return false;
185  }
186  }
187  }
188  else
189  { // single
190  if(dynamic_cast<const SiStripRecHit2D *>(tRecHit) != 0)
191  { // normal
192  const SiStripRecHit2D* recHit =
193  dynamic_cast<const SiStripRecHit2D *>(tRecHit);
194 
195  if(recHit != 0)
196  if(! theFilter->isCompatible(*recHit, gdir))
197  {
198  LogTrace("TrajectFilter")
199  << " [TrajectFilter] fail strip single";
200  return false;
201  }
202  }
203  else
204  { // projected
205  const ProjectedSiStripRecHit2D* recHit =
206  dynamic_cast<const ProjectedSiStripRecHit2D *>(tRecHit);
207 
208  if(recHit != 0)
209  if(! theFilter->isCompatible(recHit->originalHit(), gdir))
210  {
211  LogTrace("TrajectFilter")
212  << " [TrajectFilter] fail strip projected";
213  return false;
214  }
215  }
216  }
217  }
218  }
219  }
220 
221  return true;
222 }
const_iterator rend() const
Definition: bqueue.h:164
assert(m_qm.get())
const DataContainer & measurements() const
const ClusterShapeHitFilter * theFilter
bool isTrackerStrip(const GeomDetEnumerators::SubDetector m)
const GeomDet * det() const
const SiPixelClusterShapeCache * theCache
#define LogTrace(id)
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, PixelData const *pd=0) const
const_iterator rbegin() const
Definition: bqueue.h:163
SiStripRecHit2D originalHit() const
virtual TrackingRecHit const * hit() const
SiStripRecHit2D stereoHit() const
bool isValid() const
SiStripRecHit2D monoHit() const
virtual SubDetector subDetector() const
Which subdetector.
Definition: GeomDet.cc:49
GlobalVector globalDirection() const
Our base class.
Definition: SiPixelRecHit.h:23
bool ClusterShapeTrajectoryFilter::toBeContinued ( Trajectory trajectory) const
virtual

Implements TrajectoryFilter.

Definition at line 60 of file ClusterShapeTrajectoryFilter.cc.

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

61 {
63  vector<TrajectoryMeasurement> tms = trajectory.measurements();
64 
65  for(vector<TrajectoryMeasurement>::const_iterator
66  tm = tms.begin(); tm!= tms.end(); tm++)
67  {
68  const TrackingRecHit* ttRecHit = &(*((*tm).recHit()));
69 
70  if(ttRecHit->isValid())
71  {
72  const TrackingRecHit * tRecHit = ttRecHit->hit();
73 
74  TrajectoryStateOnSurface ts = (*tm).updatedState();
75  const GlobalVector gdir = ts.globalDirection();
76 
82  { // pixel
83  const SiPixelRecHit* recHit =
84  dynamic_cast<const SiPixelRecHit *>(tRecHit);
85 
86  if(recHit != 0)
87  return theFilter->isCompatible(*recHit, gdir, *theCache);
88  }
89  else if(GeomDetEnumerators::isTrackerStrip(ttRecHit->det()->subDetector()))
90  { // strip
91  if(dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit) != 0)
92  { // glued
93  const SiStripMatchedRecHit2D* recHit =
94  dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit);
95 
96  if(recHit != 0)
97  {
98  return (theFilter->isCompatible(recHit->monoHit() , gdir) &&
99  theFilter->isCompatible(recHit->stereoHit(), gdir));
100  }
101  }
102  else
103  { // single
104  if(dynamic_cast<const SiStripRecHit2D *>(tRecHit) != 0)
105  { // normal
106  const SiStripRecHit2D* recHit =
107  dynamic_cast<const SiStripRecHit2D *>(tRecHit);
108 
109  if(recHit != 0)
110  return theFilter->isCompatible(*recHit, gdir);
111  }
112  else
113  { // projected
114  const ProjectedSiStripRecHit2D* recHit =
115  dynamic_cast<const ProjectedSiStripRecHit2D *>(tRecHit);
116 
117  if(recHit != 0)
118  return theFilter->isCompatible(recHit->originalHit(), gdir);
119  }
120  }
121  }
122  }
123  }
124 
125  return true;
126 }
assert(m_qm.get())
const ClusterShapeHitFilter * theFilter
bool isTrackerStrip(const GeomDetEnumerators::SubDetector m)
DataContainer const & measurements() const
Definition: Trajectory.h:203
const GeomDet * det() const
const SiPixelClusterShapeCache * theCache
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, PixelData const *pd=0) const
SiStripRecHit2D originalHit() const
virtual TrackingRecHit const * hit() const
SiStripRecHit2D stereoHit() const
bool isValid() const
SiStripRecHit2D monoHit() const
virtual SubDetector subDetector() const
Which subdetector.
Definition: GeomDet.cc:49
GlobalVector globalDirection() const
Our base class.
Definition: SiPixelRecHit.h:23

Member Data Documentation

const SiPixelClusterShapeCache* ClusterShapeTrajectoryFilter::theCache
private

Definition at line 36 of file ClusterShapeTrajectoryFilter.h.

Referenced by setEvent().

edm::EDGetTokenT<SiPixelClusterShapeCache> ClusterShapeTrajectoryFilter::theCacheToken
private

Definition at line 35 of file ClusterShapeTrajectoryFilter.h.

Referenced by setEvent().

const ClusterShapeHitFilter* ClusterShapeTrajectoryFilter::theFilter
private

Definition at line 37 of file ClusterShapeTrajectoryFilter.h.

Referenced by setEvent().