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 Member Functions | Private Attributes
PixelClusterShapeSeedComparitor Class Reference
Inheritance diagram for PixelClusterShapeSeedComparitor:
SeedComparitor

Public Member Functions

virtual bool compatible (const SeedingHitSet &hits, const TrackingRegion &region) const override
 
virtual bool compatible (const TrajectorySeed &seed) const override
 
virtual bool compatible (const TrajectoryStateOnSurface &, SeedingHitSet::ConstRecHitPointer hit) const override
 
virtual bool compatible (const SeedingHitSet &hits, const GlobalTrajectoryParameters &helixStateAtVertex, const FastHelix &helix, const TrackingRegion &region) const override
 
virtual bool compatible (const SeedingHitSet &hits, const GlobalTrajectoryParameters &straightLineStateAtVertex, const TrackingRegion &region) const override
 
virtual void init (const edm::Event &ev, const edm::EventSetup &es) override
 
 PixelClusterShapeSeedComparitor (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
virtual ~PixelClusterShapeSeedComparitor ()
 
- Public Member Functions inherited from SeedComparitor
virtual ~SeedComparitor ()
 

Private Member Functions

bool compatibleHit (const TrackingRecHit &hit, const GlobalVector &direction) const
 

Private Attributes

const bool filterAtHelixStage_
 
edm::ESHandle
< ClusterShapeHitFilter
filterHandle_
 
std::string filterName_
 
const bool filterPixelHits_
 
const bool filterStripHits_
 
const SiPixelClusterShapeCachepixelClusterShapeCache_
 
edm::EDGetTokenT
< SiPixelClusterShapeCache
pixelClusterShapeCacheToken_
 

Detailed Description

Definition at line 25 of file ClusterShapeSeedComparitor.cc.

Constructor & Destructor Documentation

PixelClusterShapeSeedComparitor::PixelClusterShapeSeedComparitor ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 54 of file ClusterShapeSeedComparitor.cc.

References edm::ConsumesCollector::consumes(), filterPixelHits_, edm::ParameterSet::getParameter(), and pixelClusterShapeCacheToken_.

54  :
55  filterName_(cfg.getParameter<std::string>("ClusterShapeHitFilterName")),
56  pixelClusterShapeCache_(nullptr),
57  filterAtHelixStage_(cfg.getParameter<bool>("FilterAtHelixStage")),
58  filterPixelHits_(cfg.getParameter<bool>("FilterPixelHits")),
59  filterStripHits_(cfg.getParameter<bool>("FilterStripHits"))
60 {
61  if(filterPixelHits_) {
63  }
64 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< SiPixelClusterShapeCache > pixelClusterShapeCacheToken_
const SiPixelClusterShapeCache * pixelClusterShapeCache_
PixelClusterShapeSeedComparitor::~PixelClusterShapeSeedComparitor ( )
virtual

Definition at line 66 of file ClusterShapeSeedComparitor.cc.

67 {
68 }

Member Function Documentation

virtual bool PixelClusterShapeSeedComparitor::compatible ( const SeedingHitSet hits,
const TrackingRegion region 
) const
inlineoverridevirtual

Implements SeedComparitor.

Definition at line 30 of file ClusterShapeSeedComparitor.cc.

30 { return true; }
virtual bool PixelClusterShapeSeedComparitor::compatible ( const TrajectorySeed seed) const
inlineoverridevirtual

Implements SeedComparitor.

Definition at line 31 of file ClusterShapeSeedComparitor.cc.

31 { return true; }
bool PixelClusterShapeSeedComparitor::compatible ( const TrajectoryStateOnSurface tsos,
SeedingHitSet::ConstRecHitPointer  hit 
) const
overridevirtual

Implements SeedComparitor.

Definition at line 82 of file ClusterShapeSeedComparitor.cc.

References assert(), compatibleHit(), filterAtHelixStage_, TrajectoryStateOnSurface::globalDirection(), and TrajectoryStateOnSurface::isValid().

84 {
85  if (filterAtHelixStage_) return true;
86  assert(hit->isValid() && tsos.isValid());
87  return compatibleHit(*hit, tsos.globalDirection());
88 }
assert(m_qm.get())
bool compatibleHit(const TrackingRecHit &hit, const GlobalVector &direction) const
GlobalVector globalDirection() const
bool PixelClusterShapeSeedComparitor::compatible ( const SeedingHitSet hits,
const GlobalTrajectoryParameters helixStateAtVertex,
const FastHelix helix,
const TrackingRegion region 
) const
overridevirtual

Implements SeedComparitor.

Definition at line 99 of file ClusterShapeSeedComparitor.cc.

References FastHelix::circle(), compatibleHit(), filterAtHelixStage_, i, FastHelix::isValid(), GlobalTrajectoryParameters::momentum(), gen::n, perp2(), GlobalTrajectoryParameters::position(), SeedingHitSet::size(), mathSSE::sqrt(), PV3DBase< T, PVType, FrameType >::x(), FastCircle::x0(), PV3DBase< T, PVType, FrameType >::y(), FastCircle::y0(), and PV3DBase< T, PVType, FrameType >::z().

103 {
104  if (!filterAtHelixStage_) return true;
105 
106  if(!helix.isValid()) edm::LogWarning("InvalidHelix") << "PixelClusterShapeSeedComparitor helix is not valid, result is bad";
107 
108  float xc = helix.circle().x0(), yc = helix.circle().y0();
109 
110  GlobalPoint vertex = helixStateAtVertex.position();
111  GlobalVector momvtx = helixStateAtVertex.momentum();
112  float x0 = vertex.x(), y0 = vertex.y();
113  for (unsigned int i = 0, n = hits.size(); i < n; ++i) {
114  auto const & hit = *hits[i];
115  GlobalPoint pos = hit.globalPosition();
116  float x1 = pos.x(), y1 = pos.y(), dx1 = x1 - xc, dy1 = y1 - yc;
117 
118  // now figure out the proper tangent vector
119  float perpx = -dy1, perpy = dx1;
120  if (perpx * (x1-x0) + perpy * (y1 - y0) < 0) {
121  perpy = -perpy; perpx = -perpx;
122  }
123 
124  // now normalize (perpx, perpy, 1.0) to momentum (px, py, pz)
125  float perp2 = perpx*perpx + perpy*perpy;
126  float pmom2 = momvtx.x()*momvtx.x() + momvtx.y()*momvtx.y(), momz2 = momvtx.z()*momvtx.z(), mom2 = pmom2 + momz2;
127  float perpscale = sqrt(pmom2/mom2 / perp2), zscale = sqrt((1-pmom2/mom2));
128  GlobalVector gdir(perpx*perpscale, perpy*perpscale, (momvtx.z() > 0 ? zscale : -zscale));
129 
130  if (!compatibleHit(hit, gdir)) {
131  return false; // not yet
132  }
133  }
134  return true;
135 }
int i
Definition: DBlmapReader.cc:9
bool isValid() const
Definition: FastHelix.h:63
double x0() const
Definition: FastCircle.h:50
T y() const
Definition: PV3DBase.h:63
const FastCircle & circle() const
Definition: FastHelix.h:67
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
T perp2() const
Squared magnitude of transverse component.
double y0() const
Definition: FastCircle.h:52
unsigned int size() const
Definition: SeedingHitSet.h:44
bool compatibleHit(const TrackingRecHit &hit, const GlobalVector &direction) const
T x() const
Definition: PV3DBase.h:62
bool PixelClusterShapeSeedComparitor::compatible ( const SeedingHitSet hits,
const GlobalTrajectoryParameters straightLineStateAtVertex,
const TrackingRegion region 
) const
overridevirtual

Implements SeedComparitor.

Definition at line 91 of file ClusterShapeSeedComparitor.cc.

94 {
95  return true;
96 }
bool PixelClusterShapeSeedComparitor::compatibleHit ( const TrackingRecHit hit,
const GlobalVector direction 
) const
private

Definition at line 138 of file ClusterShapeSeedComparitor.cc.

References assert(), edm::hlt::Exception, filterHandle_, filterPixelHits_, filterStripHits_, TrackingRecHit::geographicalId(), SiStripMatchedRecHit2D::monoCluster(), SiStripMatchedRecHit2D::monoId(), ProjectedSiStripRecHit2D::originalHit(), pixelClusterShapeCache_, SiStripMatchedRecHit2D::stereoCluster(), SiStripMatchedRecHit2D::stereoId(), and DetId::subdetId().

Referenced by compatible().

139 {
140  if (hit.geographicalId().subdetId() <= 2) {
141  if (!filterPixelHits_) return true;
142  const SiPixelRecHit *pixhit = dynamic_cast<const SiPixelRecHit *>(&hit);
143  if (pixhit == 0) throw cms::Exception("LogicError", "Found a valid hit on the pixel detector which is not a SiPixelRecHit\n");
144  //printf("Cheching hi hit on detid %10d, local direction is x = %9.6f, y = %9.6f, z = %9.6f\n", hit.geographicalId().rawId(), direction.x(), direction.y(), direction.z());
145  return filterHandle_->isCompatible(*pixhit, direction, *pixelClusterShapeCache_);
146  } else {
147  if (!filterStripHits_) return true;
148  const std::type_info &tid = typeid(*&hit);
149  if (tid == typeid(SiStripMatchedRecHit2D)) {
150  const SiStripMatchedRecHit2D* matchedHit = dynamic_cast<const SiStripMatchedRecHit2D *>(&hit);
151  assert(matchedHit != 0);
152  return (filterHandle_->isCompatible(DetId(matchedHit->monoId()), matchedHit->monoCluster(), direction) &&
153  filterHandle_->isCompatible(DetId(matchedHit->stereoId()), matchedHit->stereoCluster(), direction));
154  } else if (tid == typeid(SiStripRecHit2D)) {
155  const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D *>(&hit);
156  assert(recHit != 0);
157  return filterHandle_->isCompatible(*recHit, direction);
158  } else if (tid == typeid(ProjectedSiStripRecHit2D)) {
159  const ProjectedSiStripRecHit2D* precHit = dynamic_cast<const ProjectedSiStripRecHit2D *>(&hit);
160  assert(precHit != 0);
161  return filterHandle_->isCompatible(precHit->originalHit(), direction);
162  } else {
163  //printf("Questo e' un %s, che ci fo?\n", tid.name());
164  return true;
165  }
166  }
167 }
unsigned int stereoId() const
SiStripCluster const & monoCluster() const
assert(m_qm.get())
edm::ESHandle< ClusterShapeHitFilter > filterHandle_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
SiStripRecHit2D originalHit() const
Definition: DetId.h:18
SiStripCluster const & stereoCluster() const
unsigned int monoId() const
const SiPixelClusterShapeCache * pixelClusterShapeCache_
DetId geographicalId() const
Our base class.
Definition: SiPixelRecHit.h:23
void PixelClusterShapeSeedComparitor::init ( const edm::Event ev,
const edm::EventSetup es 
)
overridevirtual

Implements SeedComparitor.

Definition at line 71 of file ClusterShapeSeedComparitor.cc.

References filterHandle_, filterName_, filterPixelHits_, edm::EventSetup::get(), edm::Event::getByToken(), pixelClusterShapeCache_, pixelClusterShapeCacheToken_, and edm::Handle< T >::product().

71  {
73  if(filterPixelHits_) {
77  }
78 }
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
edm::ESHandle< ClusterShapeHitFilter > filterHandle_
edm::EDGetTokenT< SiPixelClusterShapeCache > pixelClusterShapeCacheToken_
T const * product() const
Definition: Handle.h:81
const T & get() const
Definition: EventSetup.h:55
const SiPixelClusterShapeCache * pixelClusterShapeCache_

Member Data Documentation

const bool PixelClusterShapeSeedComparitor::filterAtHelixStage_
private

Definition at line 49 of file ClusterShapeSeedComparitor.cc.

Referenced by compatible().

edm::ESHandle<ClusterShapeHitFilter> PixelClusterShapeSeedComparitor::filterHandle_
private

Definition at line 46 of file ClusterShapeSeedComparitor.cc.

Referenced by compatibleHit(), and init().

std::string PixelClusterShapeSeedComparitor::filterName_
private

Definition at line 45 of file ClusterShapeSeedComparitor.cc.

Referenced by init().

const bool PixelClusterShapeSeedComparitor::filterPixelHits_
private
const bool PixelClusterShapeSeedComparitor::filterStripHits_
private

Definition at line 50 of file ClusterShapeSeedComparitor.cc.

Referenced by compatibleHit().

const SiPixelClusterShapeCache* PixelClusterShapeSeedComparitor::pixelClusterShapeCache_
private

Definition at line 48 of file ClusterShapeSeedComparitor.cc.

Referenced by compatibleHit(), and init().

edm::EDGetTokenT<SiPixelClusterShapeCache> PixelClusterShapeSeedComparitor::pixelClusterShapeCacheToken_
private

Definition at line 47 of file ClusterShapeSeedComparitor.cc.

Referenced by init(), and PixelClusterShapeSeedComparitor().