CMS 3D CMS Logo

List of all members | Classes | Public Member Functions | Private Attributes
PixelHitMatcher Class Reference

#include <PixelHitMatcher.h>

Classes

struct  BarrelMeasurementEstimator
 
struct  ForwardMeasurementEstimator
 

Public Member Functions

std::vector< SeedWithInfooperator() (const std::vector< const TrajectorySeedCollection * > &seedsV, const GlobalPoint &xmeas, const GlobalPoint &vprim, float energy, int charge) const
 
 PixelHitMatcher (float phi1min, float phi1max, float phi2minB, float phi2maxB, float phi2minF, float phi2maxF, float z2minB, float z2maxB, float r2minF, float r2maxF, float rMinI, float rMaxI, bool useRecoVertex)
 
void set1stLayer (float dummyphi1min, float dummyphi1max)
 
void set1stLayerZRange (float zmin1, float zmax1)
 
void set2ndLayer (float dummyphi2minB, float dummyphi2maxB, float dummyphi2minF, float dummyphi2maxF)
 
void setES (const MagneticField *, const TrackerGeometry *trackerGeometry)
 

Private Attributes

BarrelMeasurementEstimator meas1stBLayer
 
ForwardMeasurementEstimator meas1stFLayer
 
BarrelMeasurementEstimator meas2ndBLayer
 
ForwardMeasurementEstimator meas2ndFLayer
 
std::unique_ptr< PropagatorWithMaterialprop1stLayer
 
std::unique_ptr< PropagatorWithMaterialprop2ndLayer
 
const MagneticFieldtheMagField
 
const TrackerGeometrytheTrackerGeometry
 
const bool useRecoVertex_
 

Detailed Description

Description: Class to match an ECAL cluster to the pixel hits. Two compatible hits in the pixel layers are required.

Implementation: future redesign

Definition at line 44 of file PixelHitMatcher.h.

Constructor & Destructor Documentation

PixelHitMatcher::PixelHitMatcher ( float  phi1min,
float  phi1max,
float  phi2minB,
float  phi2maxB,
float  phi2minF,
float  phi2maxF,
float  z2minB,
float  z2maxB,
float  r2minF,
float  r2maxF,
float  rMinI,
float  rMaxI,
bool  useRecoVertex 
)

Definition at line 61 of file PixelHitMatcher.cc.

74  : //zmin1 and zmax1 are dummy at this moment, set from beamspot later
75  meas1stBLayer{phi1min, phi1max, 0., 0.},
76  meas2ndBLayer{phi2minB, phi2maxB, z2minB, z2maxB},
77  meas1stFLayer{phi1min, phi1max, 0., 0., rMinI, rMaxI},
78  meas2ndFLayer{phi2minF, phi2maxF, r2minF, r2maxF, rMinI, rMaxI},
79  prop1stLayer(nullptr),
80  prop2ndLayer(nullptr),
BarrelMeasurementEstimator meas2ndBLayer
BarrelMeasurementEstimator meas1stBLayer
std::unique_ptr< PropagatorWithMaterial > prop2ndLayer
std::unique_ptr< PropagatorWithMaterial > prop1stLayer
ForwardMeasurementEstimator meas2ndFLayer
ForwardMeasurementEstimator meas1stFLayer
const bool useRecoVertex_

Member Function Documentation

std::vector< SeedWithInfo > PixelHitMatcher::operator() ( const std::vector< const TrajectorySeedCollection * > &  seedsV,
const GlobalPoint xmeas,
const GlobalPoint vprim,
float  energy,
int  charge 
) const

Definition at line 112 of file PixelHitMatcher.cc.

References funct::cos(), TrackerGeometry::dets(), HLT_2018_cff::distance, dt, GeomDet::gdetIndex(), FTSFromVertexToPointFactory::get(), hcalSimParameters_cfi::he, hfClusterShapes_cfi::hits, mps_fire::i, globals_cff::id1, globals_cff::id2, meas1stBLayer, meas1stFLayer, meas2ndBLayer, meas2ndFLayer, FreeTrajectoryState::momentum(), PV3DBase< T, PVType, FrameType >::perp(), FreeTrajectoryState::position(), prop1stLayer, prop2ndLayer, mps_fire::result, SimDataFormats::CaloAnalysis::sc, SurveyInfoScenario_cff::seed, mathSSE::sqrt(), DetId::subdetId(), GeomDet::surface(), theMagField, theTrackerGeometry, useRecoVertex_, bphysicsOniaDQM_cfi::vertex, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), PV3DBase< T, PVType, FrameType >::z(), and HLT_FULL_cff::zVertex.

116  {
117  auto xmeas_r = xmeas.perp();
118 
119  const float phicut = std::cos(2.5);
120 
123  TrajectoryStateOnSurface tsos(fts, *bpb(fts.position(), fts.momentum()));
124 
125  std::vector<SeedWithInfo> result;
126  unsigned int allSeedsSize = 0;
127  for (auto const sc : seedsV)
128  allSeedsSize += sc->size();
129 
131 
132  auto ndets = theTrackerGeometry->dets().size();
133 
134  int iTsos[ndets];
135  for (auto &i : iTsos)
136  i = -1;
137  std::vector<TrajectoryStateOnSurface> vTsos;
138  vTsos.reserve(allSeedsSize);
139 
140  for (const auto seeds : seedsV) {
141  for (const auto &seed : *seeds) {
142  std::vector<GlobalPoint> hitGpMap;
143  if (seed.nHits() > 9) {
144  edm::LogWarning("GsfElectronAlgo|UnexpectedSeed") << "We cannot deal with seeds having more than 9 hits.";
145  continue;
146  }
147 
148  const TrajectorySeed::range &hits = seed.recHits();
149  // cache the global points
150 
151  for (auto it = hits.first; it != hits.second; ++it) {
152  hitGpMap.emplace_back(it->globalPosition());
153  }
154 
155  //iterate on the hits
156  auto he = hits.second - 1;
157  for (auto it1 = hits.first; it1 < he; ++it1) {
158  if (!it1->isValid())
159  continue;
160  auto idx1 = std::distance(hits.first, it1);
161  const DetId id1 = it1->geographicalId();
162  const GeomDet *geomdet1 = it1->det();
163 
164  auto ix1 = geomdet1->gdetIndex();
165 
166  /* VI: this generates regression (other cut is just in phi). in my opinion it is safe and makes sense
167  auto away = geomdet1->position().basicVector().dot(xmeas.basicVector()) <0;
168  if (away) continue;
169  */
170 
171  const GlobalPoint &hit1Pos = hitGpMap[idx1];
172  auto dt = hit1Pos.x() * xmeas.x() + hit1Pos.y() * xmeas.y();
173  if (dt < 0)
174  continue;
175  if (dt < phicut * (xmeas_r * hit1Pos.perp()))
176  continue;
177 
178  if (iTsos[ix1] < 0) {
179  iTsos[ix1] = vTsos.size();
180  vTsos.push_back(prop1stLayer->propagate(tsos, geomdet1->surface()));
181  }
182  auto tsos1 = &vTsos[iTsos[ix1]];
183 
184  if (!tsos1->isValid())
185  continue;
186  bool est = id1.subdetId() % 2 ? meas1stBLayer(vprim, *tsos1, hit1Pos) : meas1stFLayer(vprim, *tsos1, hit1Pos);
187  if (!est)
188  continue;
189  EleRelPointPair pp1(hit1Pos, tsos1->globalParameters().position(), vprim);
190  const math::XYZPoint relHit1Pos(hit1Pos - vprim), relTSOSPos(tsos1->globalParameters().position() - vprim);
191  const int subDet1 = id1.subdetId();
192  const float dRz1 = (id1.subdetId() % 2 ? pp1.dZ() : pp1.dPerp());
193  const float dPhi1 = pp1.dPhi();
194  // setup our vertex
195  double zVertex;
196  if (!useRecoVertex_) {
197  // we don't know the z vertex position, get it from linear extrapolation
198  // compute the z vertex from the cluster point and the found pixel hit
199  const double pxHit1z = hit1Pos.z();
200  const double pxHit1x = hit1Pos.x();
201  const double pxHit1y = hit1Pos.y();
202  const double r1diff =
203  std::sqrt((pxHit1x - vprim.x()) * (pxHit1x - vprim.x()) + (pxHit1y - vprim.y()) * (pxHit1y - vprim.y()));
204  const double r2diff =
205  std::sqrt((xmeas.x() - pxHit1x) * (xmeas.x() - pxHit1x) + (xmeas.y() - pxHit1y) * (xmeas.y() - pxHit1y));
206  zVertex = pxHit1z - r1diff * (xmeas.z() - pxHit1z) / r2diff;
207  } else {
208  // here use rather the reco vertex z position
209  zVertex = vprim.z();
210  }
211  GlobalPoint vertex(vprim.x(), vprim.y(), zVertex);
213  // now find the matching hit
214  for (auto it2 = it1 + 1; it2 != hits.second; ++it2) {
215  if (!it2->isValid())
216  continue;
217  auto idx2 = std::distance(hits.first, it2);
218  const DetId id2 = it2->geographicalId();
219  const GeomDet *geomdet2 = it2->det();
220  const auto det_key = std::make_pair(geomdet2->gdetIndex(), hit1Pos);
221  const TrajectoryStateOnSurface *tsos2;
222  auto tsos2_itr = mapTsos2Fast.find(det_key);
223  if (tsos2_itr != mapTsos2Fast.end()) {
224  tsos2 = &(tsos2_itr->second);
225  } else {
226  auto empl_result = mapTsos2Fast.emplace(det_key, prop2ndLayer->propagate(fts2, geomdet2->surface()));
227  tsos2 = &(empl_result.first->second);
228  }
229  if (!tsos2->isValid())
230  continue;
231  const GlobalPoint &hit2Pos = hitGpMap[idx2];
232  bool est2 =
233  id2.subdetId() % 2 ? meas2ndBLayer(vertex, *tsos2, hit2Pos) : meas2ndFLayer(vertex, *tsos2, hit2Pos);
234  if (est2) {
235  EleRelPointPair pp2(hit2Pos, tsos2->globalParameters().position(), vertex);
236  const int subDet2 = id2.subdetId();
237  const float dRz2 = (subDet2 % 2 == 1) ? pp2.dZ() : pp2.dPerp();
238  const float dPhi2 = pp2.dPhi();
239  const unsigned char hitsMask = (1 << idx1) | (1 << idx2);
240  result.push_back({seed, hitsMask, subDet2, dRz2, dPhi2, subDet1, dRz1, dPhi1});
241  }
242  } // inner loop on hits
243  } // outer loop on hits
244  } // loop on seeds
245  } //loop on vector of seeds
246 
247  return result;
248 }
float dt
Definition: AMPTWrapper.h:136
BarrelMeasurementEstimator meas2ndBLayer
std::pair< const_iterator, const_iterator > range
static FreeTrajectoryState get(MagneticField const &magField, GlobalPoint const &xmeas, GlobalPoint const &xvert, float momentum, TrackCharge charge)
T perp() const
Definition: PV3DBase.h:69
BarrelMeasurementEstimator meas1stBLayer
int gdetIndex() const
Definition: GeomDet.h:87
T y() const
Definition: PV3DBase.h:60
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
std::unique_ptr< PropagatorWithMaterial > prop2ndLayer
std::unique_ptr< PropagatorWithMaterial > prop1stLayer
const DetContainer & dets() const override
Returm a vector of all GeomDet (including all GeomDetUnits)
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:48
const MagneticField * theMagField
GlobalVector momentum() const
Definition: DetId.h:17
GlobalPoint position() const
std::unordered_map< std::pair< int, GlobalPoint >, T, HashIntGlobalPointPair > IntGlobalPointPairUnorderedMap
Definition: utils.h:20
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
ForwardMeasurementEstimator meas2ndFLayer
ForwardMeasurementEstimator meas1stFLayer
const bool useRecoVertex_
T x() const
Definition: PV3DBase.h:59
const TrackerGeometry * theTrackerGeometry
void PixelHitMatcher::set1stLayer ( float  dummyphi1min,
float  dummyphi1max 
)
void PixelHitMatcher::set1stLayerZRange ( float  zmin1,
float  zmax1 
)
void PixelHitMatcher::set2ndLayer ( float  dummyphi2minB,
float  dummyphi2maxB,
float  dummyphi2minF,
float  dummyphi2maxF 
)
void PixelHitMatcher::setES ( const MagneticField magField,
const TrackerGeometry trackerGeometry 
)

Definition at line 104 of file PixelHitMatcher.cc.

References alongMomentum, EgHLTOffHistBins_cfi::mass, oppositeToMomentum, prop1stLayer, prop2ndLayer, theMagField, and theTrackerGeometry.

Referenced by ElectronSeedGenerator::setupES().

104  {
105  theMagField = magField;
106  theTrackerGeometry = trackerGeometry;
107  float mass = .000511; // electron propagation
108  prop1stLayer = std::make_unique<PropagatorWithMaterial>(oppositeToMomentum, mass, theMagField);
109  prop2ndLayer = std::make_unique<PropagatorWithMaterial>(alongMomentum, mass, theMagField);
110 }
std::unique_ptr< PropagatorWithMaterial > prop2ndLayer
std::unique_ptr< PropagatorWithMaterial > prop1stLayer
const MagneticField * theMagField
const TrackerGeometry * theTrackerGeometry

Member Data Documentation

BarrelMeasurementEstimator PixelHitMatcher::meas1stBLayer
private

Definition at line 93 of file PixelHitMatcher.h.

Referenced by operator()(), set1stLayer(), and set1stLayerZRange().

ForwardMeasurementEstimator PixelHitMatcher::meas1stFLayer
private

Definition at line 95 of file PixelHitMatcher.h.

Referenced by operator()(), set1stLayer(), and set1stLayerZRange().

BarrelMeasurementEstimator PixelHitMatcher::meas2ndBLayer
private

Definition at line 94 of file PixelHitMatcher.h.

Referenced by operator()(), and set2ndLayer().

ForwardMeasurementEstimator PixelHitMatcher::meas2ndFLayer
private

Definition at line 96 of file PixelHitMatcher.h.

Referenced by operator()(), and set2ndLayer().

std::unique_ptr<PropagatorWithMaterial> PixelHitMatcher::prop1stLayer
private

Definition at line 97 of file PixelHitMatcher.h.

Referenced by operator()(), and setES().

std::unique_ptr<PropagatorWithMaterial> PixelHitMatcher::prop2ndLayer
private

Definition at line 98 of file PixelHitMatcher.h.

Referenced by operator()(), and setES().

const MagneticField* PixelHitMatcher::theMagField
private

Definition at line 99 of file PixelHitMatcher.h.

Referenced by operator()(), and setES().

const TrackerGeometry* PixelHitMatcher::theTrackerGeometry
private

Definition at line 100 of file PixelHitMatcher.h.

Referenced by operator()(), and setES().

const bool PixelHitMatcher::useRecoVertex_
private

Definition at line 101 of file PixelHitMatcher.h.

Referenced by operator()().