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