CMS 3D CMS Logo

ClusterShapeSeedComparitor.cc
Go to the documentation of this file.
20 #include <cstdio>
21 #include <cassert>
22 
23 #include <iostream>
24 
26 public:
29  void init(const edm::Event &ev, const edm::EventSetup &es) override;
30  bool compatible(const SeedingHitSet &hits) const override { return true; }
32  bool compatible(const SeedingHitSet &hits,
33  const GlobalTrajectoryParameters &helixStateAtVertex,
34  const FastHelix &helix) const override;
35 
36 private:
37  bool compatibleHit(const TrackingRecHit &hit, const GlobalVector &direction) const;
38 
43  const bool filterAtHelixStage_;
45 };
46 
49  : filterName_(cfg.getParameter<std::string>("ClusterShapeHitFilterName")),
51  filterAtHelixStage_(cfg.getParameter<bool>("FilterAtHelixStage")),
52  filterPixelHits_(cfg.getParameter<bool>("FilterPixelHits")),
53  filterStripHits_(cfg.getParameter<bool>("FilterStripHits")) {
54  if (filterPixelHits_) {
56  iC.consumes<SiPixelClusterShapeCache>(cfg.getParameter<edm::InputTag>("ClusterShapeCacheSrc"));
57  }
58 }
59 
61 
64  if (filterPixelHits_) {
68  }
69 }
70 
74  return true;
75  assert(hit->isValid() && tsos.isValid());
76  return compatibleHit(*hit, tsos.globalDirection());
77 }
78 
80  const GlobalTrajectoryParameters &helixStateAtVertex,
81  const FastHelix &helix) const {
83  return true;
84 
85  if (!helix.isValid() //check still if it's a straight line, which are OK
86  && !helix.circle().isLine()) //complain if it's not even a straight line
87  edm::LogWarning("InvalidHelix") << "PixelClusterShapeSeedComparitor helix is not valid, result is bad";
88 
89  float xc = helix.circle().x0(), yc = helix.circle().y0();
90 
91  GlobalPoint vertex = helixStateAtVertex.position();
92  GlobalVector momvtx = helixStateAtVertex.momentum();
93  float x0 = vertex.x(), y0 = vertex.y();
94  for (unsigned int i = 0, n = hits.size(); i < n; ++i) {
95  auto const &hit = *hits[i];
96  GlobalPoint pos = hit.globalPosition();
97  float x1 = pos.x(), y1 = pos.y(), dx1 = x1 - xc, dy1 = y1 - yc;
98 
99  // now figure out the proper tangent vector
100  float perpx = -dy1, perpy = dx1;
101  if (perpx * (x1 - x0) + perpy * (y1 - y0) < 0) {
102  perpy = -perpy;
103  perpx = -perpx;
104  }
105 
106  // now normalize (perpx, perpy, 1.0) to momentum (px, py, pz)
107  float perp2 = perpx * perpx + perpy * perpy;
108  float pmom2 = momvtx.x() * momvtx.x() + momvtx.y() * momvtx.y(), momz2 = momvtx.z() * momvtx.z(),
109  mom2 = pmom2 + momz2;
110  float perpscale = sqrt(pmom2 / mom2 / perp2), zscale = sqrt((1 - pmom2 / mom2));
111  GlobalVector gdir(perpx * perpscale, perpy * perpscale, (momvtx.z() > 0 ? zscale : -zscale));
112 
113  if (!compatibleHit(hit, gdir)) {
114  return false; // not yet
115  }
116  }
117  return true;
118 }
119 
121  if (hit.geographicalId().subdetId() <= 2) {
122  if (!filterPixelHits_)
123  return true;
124  const SiPixelRecHit *pixhit = dynamic_cast<const SiPixelRecHit *>(&hit);
125  if (pixhit == nullptr)
126  throw cms::Exception("LogicError", "Found a valid hit on the pixel detector which is not a SiPixelRecHit\n");
127  //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());
128  return filterHandle_->isCompatible(*pixhit, direction, *pixelClusterShapeCache_);
129  } else {
130  if (!filterStripHits_)
131  return true;
132  const std::type_info &tid = typeid(*&hit);
133  if (tid == typeid(SiStripMatchedRecHit2D)) {
134  const SiStripMatchedRecHit2D *matchedHit = dynamic_cast<const SiStripMatchedRecHit2D *>(&hit);
135  assert(matchedHit != nullptr);
136  return (filterHandle_->isCompatible(DetId(matchedHit->monoId()), matchedHit->monoCluster(), direction) &&
137  filterHandle_->isCompatible(DetId(matchedHit->stereoId()), matchedHit->stereoCluster(), direction));
138  } else if (tid == typeid(SiStripRecHit2D)) {
139  const SiStripRecHit2D *recHit = dynamic_cast<const SiStripRecHit2D *>(&hit);
140  assert(recHit != nullptr);
141  return filterHandle_->isCompatible(*recHit, direction);
142  } else if (tid == typeid(ProjectedSiStripRecHit2D)) {
143  const ProjectedSiStripRecHit2D *precHit = dynamic_cast<const ProjectedSiStripRecHit2D *>(&hit);
144  assert(precHit != nullptr);
145  return filterHandle_->isCompatible(precHit->originalHit(), direction);
146  } else {
147  //printf("Questo e' un %s, che ci fo?\n", tid.name());
148  return true;
149  }
150  }
151 }
152 
153 DEFINE_EDM_PLUGIN(SeedComparitorFactory, PixelClusterShapeSeedComparitor, "PixelClusterShapeSeedComparitor");
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
bool isValid() const
Definition: FastHelix.h:57
unsigned int stereoId() const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
double x0() const
Definition: FastCircle.h:43
bool isLine() const
Definition: FastCircle.h:51
PixelClusterShapeSeedComparitor(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
SiStripCluster const & monoCluster() const
#define nullptr
T y() const
Definition: PV3DBase.h:60
bool ev
edm::ESHandle< ClusterShapeHitFilter > filterHandle_
const FastCircle & circle() const
Definition: FastHelix.h:61
BaseTrackerRecHit const * ConstRecHitPointer
Definition: SeedingHitSet.h:10
edm::EDGetTokenT< SiPixelClusterShapeCache > pixelClusterShapeCacheToken_
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
bool isCompatible(const SiPixelRecHit &recHit, const LocalVector &ldir, const SiPixelClusterShapeCache &clusterShapeCache, PixelData const *pd=0) const
bool compatible(const SeedingHitSet &hits) const override
SiStripRecHit2D originalHit() const
Definition: DetId.h:17
T const * product() const
Definition: Handle.h:69
T perp2() const
Squared magnitude of transverse component.
double y0() const
Definition: FastCircle.h:45
SiStripCluster const & stereoCluster() const
T get() const
Definition: EventSetup.h:73
unsigned int size() const
Definition: SeedingHitSet.h:41
bool compatibleHit(const TrackingRecHit &hit, const GlobalVector &direction) const
#define DEFINE_EDM_PLUGIN(factory, type, name)
unsigned int monoId() const
const SiPixelClusterShapeCache * pixelClusterShapeCache_
DetId geographicalId() const
T x() const
Definition: PV3DBase.h:59
void init(const edm::Event &ev, const edm::EventSetup &es) override
GlobalVector globalDirection() const
Our base class.
Definition: SiPixelRecHit.h:23