CMS 3D CMS Logo

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)
 
std::string name () const override
 
bool qualityFilter (const TempTrajectory &) const override
 
bool qualityFilter (const Trajectory &) const override
 
void setEvent (const edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
bool toBeContinued (TempTrajectory &) const override
 
bool toBeContinued (Trajectory &) const override
 
 ~ClusterShapeTrajectoryFilter () override
 
- 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< SiPixelClusterShapeCachetheCacheToken
 
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 ( )
override

Definition at line 44 of file ClusterShapeTrajectoryFilter.cc.

45 {
46 }

Member Function Documentation

std::string ClusterShapeTrajectoryFilter::name ( void  ) const
inlineoverridevirtual

Implements TrajectoryFilter.

Definition at line 32 of file ClusterShapeTrajectoryFilter.h.

Referenced by config.CFG::__str__(), and validation.Sample::digest().

32 { return "ClusterShapeTrajectoryFilter"; }
bool ClusterShapeTrajectoryFilter::qualityFilter ( const TempTrajectory trajectory) const
overridevirtual

Implements TrajectoryFilter.

Definition at line 235 of file ClusterShapeTrajectoryFilter.cc.

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

Referenced by qualityFilter(), and toBeContinued().

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

Implements TrajectoryFilter.

Definition at line 228 of file ClusterShapeTrajectoryFilter.cc.

References qualityFilter().

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

Reimplemented from TrajectoryFilter.

Definition at line 48 of file ClusterShapeTrajectoryFilter.cc.

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

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:517
const ClusterShapeHitFilter * theFilter
const SiPixelClusterShapeCache * theCache
T const * product() const
Definition: Handle.h:74
def cache(function)
Definition: utilities.py:3
T get() const
Definition: EventSetup.h:71
T const * product() const
Definition: ESHandle.h:86
bool ClusterShapeTrajectoryFilter::toBeContinued ( TempTrajectory trajectory) const
overridevirtual

Implements TrajectoryFilter.

Definition at line 131 of file ClusterShapeTrajectoryFilter.cc.

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

Referenced by qualityFilter(), setEvent(), and toBeContinued().

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

Implements TrajectoryFilter.

Definition at line 60 of file ClusterShapeTrajectoryFilter.cc.

References TrackingRecHit::det(), TrajectoryStateOnSurface::globalDirection(), TrackingRecHit::hit(), ClusterShapeHitFilter::isCompatible(), GeomDetEnumerators::isTrackerStrip(), TrackingRecHit::isValid(), Trajectory::measurements(), SiStripMatchedRecHit2D::monoHit(), ProjectedSiStripRecHit2D::originalHit(), GeomDetEnumerators::P1PXB, GeomDetEnumerators::P1PXEC, GeomDetEnumerators::P2PXB, GeomDetEnumerators::P2PXEC, GeomDetEnumerators::PixelBarrel, GeomDetEnumerators::PixelEndcap, rpcPointValidation_cfi::recHit, SiStripMatchedRecHit2D::stereoHit(), GeomDet::subDetector(), theCache, theFilter, and toBeContinued().

61 {
62  assert(theCache);
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 
83  { // pixel
84  const SiPixelRecHit* recHit =
85  dynamic_cast<const SiPixelRecHit *>(tRecHit);
86 
87  if(recHit != nullptr)
88  return theFilter->isCompatible(*recHit, gdir, *theCache);
89  }
90  else if(GeomDetEnumerators::isTrackerStrip(ttRecHit->det()->subDetector()))
91  { // strip
92  if(dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit) != nullptr)
93  { // glued
94  const SiStripMatchedRecHit2D* recHit =
95  dynamic_cast<const SiStripMatchedRecHit2D *>(tRecHit);
96 
97  if(recHit != nullptr)
98  {
99  return (theFilter->isCompatible(recHit->monoHit() , gdir) &&
100  theFilter->isCompatible(recHit->stereoHit(), gdir));
101  }
102  }
103  else
104  { // single
105  if(dynamic_cast<const SiStripRecHit2D *>(tRecHit) != nullptr)
106  { // normal
107  const SiStripRecHit2D* recHit =
108  dynamic_cast<const SiStripRecHit2D *>(tRecHit);
109 
110  if(recHit != nullptr)
111  return theFilter->isCompatible(*recHit, gdir);
112  }
113  else
114  { // projected
115  const ProjectedSiStripRecHit2D* recHit =
116  dynamic_cast<const ProjectedSiStripRecHit2D *>(tRecHit);
117 
118  if(recHit != nullptr)
119  return theFilter->isCompatible(recHit->originalHit(), gdir);
120  }
121  }
122  }
123  }
124  }
125 
126  return true;
127 }
const ClusterShapeHitFilter * theFilter
DataContainer const & measurements() const
Definition: Trajectory.h:196
const GeomDet * det() const
const SiPixelClusterShapeCache * theCache
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, PixelData const *pd=0) const
bool isTrackerStrip(GeomDetEnumerators::SubDetector m)
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:44
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(), and toBeContinued().

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(), and toBeContinued().