CMS 3D CMS Logo

TrajSeedMatcher.cc
Go to the documentation of this file.
2 
7 
11 
14 
17 
19 
20 constexpr float TrajSeedMatcher::kElectronMass_;
21 
22 namespace {
23  auto makeMatchingCuts(std::vector<edm::ParameterSet> const& cutsPSets) {
24  std::vector<std::unique_ptr<TrajSeedMatcher::MatchingCuts> > matchingCuts;
25 
26  for (const auto& cutPSet : cutsPSets) {
27  int version = cutPSet.getParameter<int>("version");
28  switch (version) {
29  case 1:
30  matchingCuts.emplace_back(std::make_unique<TrajSeedMatcher::MatchingCutsV1>(cutPSet));
31  break;
32  case 2:
33  matchingCuts.emplace_back(std::make_unique<TrajSeedMatcher::MatchingCutsV2>(cutPSet));
34  break;
35  default:
36  throw cms::Exception("InvalidConfig") << " Error TrajSeedMatcher::TrajSeedMatcher pixel match cuts version "
37  << version << " not recognised" << std::endl;
38  }
39  }
40 
41  return matchingCuts;
42  }
43 
44  TrajSeedMatcher::SCHitMatch makeSCHitMatch(const GlobalPoint& vtxPos,
45  const TrajectoryStateOnSurface& trajState,
46  const TrackingRecHit& hit,
47  float et,
48  float eta,
49  float phi,
50  int charge,
51  int nrClus) {
52  EleRelPointPair pointPair(hit.globalPosition(), trajState.globalParameters().position(), vtxPos);
53  float dRZ = hit.geographicalId().subdetId() == PixelSubdetector::PixelBarrel ? pointPair.dZ() : pointPair.dPerp();
54  return {hit.geographicalId(), hit.globalPosition(), dRZ, pointPair.dPhi(), hit, et, eta, phi, charge, nrClus};
55  }
56 
57  const std::vector<TrajSeedMatcher::MatchInfo> makeMatchInfoVector(
58  std::vector<TrajSeedMatcher::SCHitMatch> const& posCharge,
59  std::vector<TrajSeedMatcher::SCHitMatch> const& negCharge) {
60  std::vector<TrajSeedMatcher::MatchInfo> matchInfos;
61  size_t nrHitsMax = std::max(posCharge.size(), negCharge.size());
62  for (size_t hitNr = 0; hitNr < nrHitsMax; hitNr++) {
63  DetId detIdPos = hitNr < posCharge.size() ? posCharge[hitNr].detId : DetId(0);
64  float dRZPos = hitNr < posCharge.size() ? posCharge[hitNr].dRZ : std::numeric_limits<float>::max();
65  float dPhiPos = hitNr < posCharge.size() ? posCharge[hitNr].dPhi : std::numeric_limits<float>::max();
66 
67  DetId detIdNeg = hitNr < negCharge.size() ? negCharge[hitNr].detId : DetId(0);
68  float dRZNeg = hitNr < negCharge.size() ? negCharge[hitNr].dRZ : std::numeric_limits<float>::max();
69  float dPhiNeg = hitNr < negCharge.size() ? negCharge[hitNr].dPhi : std::numeric_limits<float>::max();
70 
71  if (detIdPos != detIdNeg && (detIdPos.rawId() != 0 && detIdNeg.rawId() != 0)) {
72  cms::Exception("LogicError") << " error in " << __FILE__ << ", " << __LINE__
73  << " hits to be combined have different detIDs, this should not be possible and "
74  "nothing good will come of it";
75  }
76  DetId detId = detIdPos.rawId() != 0 ? detIdPos : detIdNeg;
77  matchInfos.push_back({detId, dRZPos, dRZNeg, dPhiPos, dPhiNeg});
78  }
79  return matchInfos;
80  }
81 }; // namespace
82 
84  : magFieldToken{cc.esConsumes()},
85  paramMagFieldToken{cc.esConsumes(pset.getParameter<edm::ESInputTag>("paramMagField"))},
86  navSchoolToken{cc.esConsumes(pset.getParameter<edm::ESInputTag>("navSchool"))},
87  detLayerGeomToken{cc.esConsumes(pset.getParameter<edm::ESInputTag>("detLayerGeom"))},
88  useRecoVertex{pset.getParameter<bool>("useRecoVertex")},
89  enableHitSkipping{pset.getParameter<bool>("enableHitSkipping")},
90  requireExactMatchCount{pset.getParameter<bool>("requireExactMatchCount")},
91  useParamMagFieldIfDefined{pset.getParameter<bool>("useParamMagFieldIfDefined")},
92  minNrHits{pset.getParameter<std::vector<unsigned int> >("minNrHits")},
93  minNrHitsValidLayerBins{pset.getParameter<std::vector<int> >("minNrHitsValidLayerBins")},
94  matchingCuts{makeMatchingCuts(pset.getParameter<std::vector<edm::ParameterSet> >("matchingCuts"))} {
95  if (minNrHitsValidLayerBins.size() + 1 != minNrHits.size()) {
96  throw cms::Exception("InvalidConfig")
97  << " TrajSeedMatcher::TrajSeedMatcher minNrHitsValidLayerBins should be 1 less than minNrHits when its "
98  << minNrHitsValidLayerBins.size() << " vs " << minNrHits.size();
99  }
100 }
101 
103  math::XYZPoint const& vprim,
105  edm::EventSetup const& iSetup,
107  : seeds_{seeds},
108  vprim_(vprim.x(), vprim.y(), vprim.z()),
109  cfg_{cfg},
110  magField_{iSetup.getData(cfg_.magFieldToken)},
111  magFieldParam_{iSetup.getData(cfg_.paramMagFieldToken)},
112  measTkEvt_{measTkEvt},
113  navSchool_{iSetup.getData(cfg_.navSchoolToken)},
114  detLayerGeom_{iSetup.getData(cfg_.detLayerGeomToken)},
115  forwardPropagator_(alongMomentum, kElectronMass_, &magField_),
116  backwardPropagator_(oppositeToMomentum, kElectronMass_, &magField_) {}
117 
120  desc.add<bool>("useRecoVertex", false);
121  desc.add<bool>("enableHitSkipping", false);
122  desc.add<bool>("requireExactMatchCount", true);
123  desc.add<bool>("useParamMagFieldIfDefined", true);
124  desc.add<edm::ESInputTag>("paramMagField", edm::ESInputTag{"", "ParabolicMf"});
125  desc.add<edm::ESInputTag>("navSchool", edm::ESInputTag{"", "SimpleNavigationSchool"});
126  desc.add<edm::ESInputTag>("detLayerGeom", edm::ESInputTag{"", "hltESPGlobalDetLayerGeometry"});
127  desc.add<std::vector<int> >("minNrHitsValidLayerBins", {4});
128  desc.add<std::vector<unsigned int> >("minNrHits", {2, 3});
129 
131  auto cutDescCases = 1 >> (edm::ParameterDescription<double>("dPhiMax", 0.04, true) and
132  edm::ParameterDescription<double>("dRZMax", 0.09, true) and
133  edm::ParameterDescription<double>("dRZMaxLowEtThres", 20., true) and
134  edm::ParameterDescription<std::vector<double> >("dRZMaxLowEtEtaBins", {1., 1.5}, true) and
135  edm::ParameterDescription<std::vector<double> >("dRZMaxLowEt", {0.09, 0.15, 0.09}, true)) or
136  2 >> (edm::ParameterDescription<std::vector<double> >("dPhiMaxHighEt", {0.003}, true) and
137  edm::ParameterDescription<std::vector<double> >("dPhiMaxHighEtThres", {0.0}, true) and
138  edm::ParameterDescription<std::vector<double> >("dPhiMaxLowEtGrad", {0.0}, true) and
139  edm::ParameterDescription<std::vector<double> >("dRZMaxHighEt", {0.005}, true) and
140  edm::ParameterDescription<std::vector<double> >("dRZMaxHighEtThres", {30}, true) and
141  edm::ParameterDescription<std::vector<double> >("dRZMaxLowEtGrad", {-0.002}, true) and
142  edm::ParameterDescription<std::vector<double> >("etaBins", {}, true));
143  cutsDesc.ifValue(edm::ParameterDescription<int>("version", 1, true), std::move(cutDescCases));
144 
146  defaults.addParameter<double>("dPhiMax", 0.04);
147  defaults.addParameter<double>("dRZMax", 0.09);
148  defaults.addParameter<double>("dRZMaxLowEtThres", 0.09);
149  defaults.addParameter<std::vector<double> >("dRZMaxLowEtEtaBins", std::vector<double>{1., 1.5});
150  defaults.addParameter<std::vector<double> >("dRZMaxLowEt", std::vector<double>{0.09, 0.09, 0.09});
151  defaults.addParameter<int>("version", 1);
152  desc.addVPSet("matchingCuts", cutsDesc, std::vector<edm::ParameterSet>{defaults, defaults, defaults});
153  return desc;
154 }
155 
156 std::vector<TrajSeedMatcher::SeedWithInfo> TrajSeedMatcher::operator()(const GlobalPoint& candPos, const float energy) {
157  clearCache();
158 
159  std::vector<SeedWithInfo> matchedSeeds;
160 
161  //these are super expensive functions
162  TrajectoryStateOnSurface scTrajStateOnSurfNeg = makeTrajStateOnSurface(candPos, energy, -1);
163  TrajectoryStateOnSurface scTrajStateOnSurfPos = makeTrajStateOnSurface(candPos, energy, 1);
164 
165  for (const auto& seed : seeds_) {
166  std::vector<SCHitMatch> matchedHitsNeg = processSeed(seed, candPos, energy, scTrajStateOnSurfNeg);
167  std::vector<SCHitMatch> matchedHitsPos = processSeed(seed, candPos, energy, scTrajStateOnSurfPos);
168 
169  int nrValidLayersPos = 0;
170  int nrValidLayersNeg = 0;
171  if (matchedHitsNeg.size() >= 2) {
172  nrValidLayersNeg = getNrValidLayersAlongTraj(matchedHitsNeg[0], matchedHitsNeg[1], candPos, energy, -1);
173  }
174  if (matchedHitsPos.size() >= 2) {
175  nrValidLayersPos = getNrValidLayersAlongTraj(matchedHitsPos[0], matchedHitsPos[1], candPos, energy, +1);
176  }
177 
178  int nrValidLayers = std::max(nrValidLayersNeg, nrValidLayersPos);
179  size_t nrHitsRequired = getNrHitsRequired(nrValidLayers);
180  bool matchCountPasses;
182  // If the input seed collection is not cross-cleaned, an exact match is necessary to
183  // prevent redundant seeds.
184  matchCountPasses = matchedHitsNeg.size() == nrHitsRequired || matchedHitsPos.size() == nrHitsRequired;
185  } else {
186  matchCountPasses = matchedHitsNeg.size() >= nrHitsRequired || matchedHitsPos.size() >= nrHitsRequired;
187  }
188  if (matchCountPasses) {
189  matchedSeeds.push_back({seed, makeMatchInfoVector(matchedHitsPos, matchedHitsNeg), nrValidLayers});
190  }
191  }
192  return matchedSeeds;
193 }
194 
195 std::vector<TrajSeedMatcher::SCHitMatch> TrajSeedMatcher::processSeed(const TrajectorySeed& seed,
196  const GlobalPoint& candPos,
197  const float energy,
198  const TrajectoryStateOnSurface& initialTrajState) {
199  //next try passing these variables in once...
200  const float candEta = candPos.eta();
201  const float candEt = energy * std::sin(candPos.theta());
202  const int charge = initialTrajState.charge();
203 
204  std::vector<SCHitMatch> matches;
205  FreeTrajectoryState firstMatchFreeTraj;
206  GlobalPoint prevHitPos;
208  const auto nCuts = cfg_.matchingCuts.size();
209  for (size_t iHit = 0;
210  matches.size() < nCuts && iHit < seed.nHits() && (cfg_.enableHitSkipping || iHit == matches.size());
211  iHit++) {
212  auto const& recHit = *(seed.recHits().begin() + iHit);
213 
214  if (!recHit.isValid()) {
215  continue;
216  }
217 
218  const bool doFirstMatch = matches.empty();
219 
220  auto const& trajState = doFirstMatch
221  ? getTrajStateFromVtx(recHit, initialTrajState, backwardPropagator_)
222  : getTrajStateFromPoint(recHit, firstMatchFreeTraj, prevHitPos, forwardPropagator_);
223  if (!trajState.isValid()) {
224  continue;
225  }
226 
227  auto const& vtxForMatchObject = doFirstMatch ? vprim_ : vertex;
228  auto match = makeSCHitMatch(vtxForMatchObject, trajState, recHit, candEt, candEta, candPos.phi(), charge, 1);
229 
230  if ((*cfg_.matchingCuts[matches.size()])(match)) {
231  matches.push_back(match);
232  if (doFirstMatch) {
233  //now we can figure out the z vertex
234  double zVertex = cfg_.useRecoVertex ? vprim_.z() : getZVtxFromExtrapolation(vprim_, match.hitPos, candPos);
236  firstMatchFreeTraj = ftsFromVertexToPoint(match.hitPos, vertex, energy, charge);
237  }
238  prevHitPos = match.hitPos;
239  }
240  }
241  return matches;
242 }
243 
244 // compute the z vertex from the candidate position and the found pixel hit
246  const GlobalPoint& hitPos,
247  const GlobalPoint& candPos) {
248  auto sq = [](float x) { return x * x; };
249  auto calRDiff = [sq](const GlobalPoint& p1, const GlobalPoint& p2) {
250  return std::sqrt(sq(p2.x() - p1.x()) + sq(p2.y() - p1.y()));
251  };
252  const double r1Diff = calRDiff(primeVtxPos, hitPos);
253  const double r2Diff = calRDiff(hitPos, candPos);
254  return hitPos.z() - r1Diff * (candPos.z() - hitPos.z()) / r2Diff;
255 }
256 
258  const TrajectoryStateOnSurface& initialState,
260  auto& trajStateFromVtxCache =
262 
263  auto key = hit.det()->gdetIndex();
264  auto res = trajStateFromVtxCache.find(key);
265  if (res != trajStateFromVtxCache.end())
266  return res->second;
267  else { //doesnt exist, need to make it
268  //FIXME: check for efficiency
269  auto val = trajStateFromVtxCache.emplace(key, propagator.propagate(initialState, hit.det()->surface()));
270  return val.first->second;
271  }
272 }
273 
275  const FreeTrajectoryState& initialState,
276  const GlobalPoint& point,
278  auto& trajStateFromPointCache =
280 
281  auto key = std::make_pair(hit.det()->gdetIndex(), point);
282  auto res = trajStateFromPointCache.find(key);
283  if (res != trajStateFromPointCache.end())
284  return res->second;
285  else { //doesnt exist, need to make it
286  //FIXME: check for efficiency
287  auto val = trajStateFromPointCache.emplace(key, propagator.propagate(initialState, hit.det()->surface()));
288  return val.first->second;
289  }
290 }
291 
293  const float energy,
294  const int charge) const {
295  auto freeTS = ftsFromVertexToPoint(pos, vprim_, energy, charge);
296  return TrajectoryStateOnSurface(freeTS, *PerpendicularBoundPlaneBuilder{}(freeTS.position(), freeTS.momentum()));
297 }
298 
304 }
305 
307  const SCHitMatch& hit1, const SCHitMatch& hit2, const GlobalPoint& candPos, const float energy, const int charge) {
308  double zVertex = cfg_.useRecoVertex ? vprim_.z() : getZVtxFromExtrapolation(vprim_, hit1.hitPos, candPos);
310 
311  auto firstMatchFreeTraj = ftsFromVertexToPoint(hit1.hitPos, vertex, energy, charge);
312  auto const& secondHitTraj = getTrajStateFromPoint(hit2.hit, firstMatchFreeTraj, hit1.hitPos, forwardPropagator_);
313  return getNrValidLayersAlongTraj(hit2.hit.geographicalId(), secondHitTraj);
314 }
315 
316 int TrajSeedMatcher::getNrValidLayersAlongTraj(const DetId& hitId, const TrajectoryStateOnSurface& hitTrajState) const {
317  const DetLayer* detLayer = detLayerGeom_.idToLayer(hitId);
318  if (detLayer == nullptr)
319  return 0;
320 
321  const FreeTrajectoryState& hitFreeState = *hitTrajState.freeState();
322  auto const inLayers = navSchool_.compatibleLayers(*detLayer, hitFreeState, oppositeToMomentum);
323  const auto outLayers = navSchool_.compatibleLayers(*detLayer, hitFreeState, alongMomentum);
324 
325  int nrValidLayers = 1; //because our current hit is also valid and wont be included in the count otherwise
326  int nrPixInLayers = 0;
327  int nrPixOutLayers = 0;
328  for (auto layer : inLayers) {
329  if (GeomDetEnumerators::isTrackerPixel(layer->subDetector())) {
330  nrPixInLayers++;
331  if (layerHasValidHits(*layer, hitTrajState, backwardPropagator_))
332  nrValidLayers++;
333  }
334  }
335  for (auto layer : outLayers) {
336  if (GeomDetEnumerators::isTrackerPixel(layer->subDetector())) {
337  nrPixOutLayers++;
338  if (layerHasValidHits(*layer, hitTrajState, forwardPropagator_))
339  nrValidLayers++;
340  }
341  }
342  return nrValidLayers;
343 }
344 
346  const TrajectoryStateOnSurface& hitSurState,
347  const Propagator& propToLayerFromState) const {
348  //FIXME: do not hardcode with werid magic numbers stolen from ancient tracking code
349  //its taken from https://cmssdt.cern.ch/dxr/CMSSW/source/RecoTracker/TrackProducer/interface/TrackProducerBase.icc#165
350  //which inspires this code
351  Chi2MeasurementEstimator estimator(30., -3.0, 0.5, 2.0, 0.5, 1.e12); // same as defauts....
352 
353  const std::vector<GeometricSearchDet::DetWithState>& detWithState =
354  layer.compatibleDets(hitSurState, propToLayerFromState, estimator);
355  if (detWithState.empty())
356  return false;
357  else {
358  DetId id = detWithState.front().first->geographicalId();
360  if (measDet.isActive())
361  return true;
362  else
363  return false;
364  }
365 }
366 
367 size_t TrajSeedMatcher::getNrHitsRequired(const int nrValidLayers) const {
368  for (size_t binNr = 0; binNr < cfg_.minNrHitsValidLayerBins.size(); binNr++) {
369  if (nrValidLayers < cfg_.minNrHitsValidLayerBins[binNr])
370  return cfg_.minNrHits[binNr];
371  }
372  return cfg_.minNrHits.back();
373 }
374 
376  : dPhiMax_(pset.getParameter<double>("dPhiMax")),
377  dRZMax_(pset.getParameter<double>("dRZMax")),
378  dRZMaxLowEtThres_(pset.getParameter<double>("dRZMaxLowEtThres")),
379  dRZMaxLowEtEtaBins_(pset.getParameter<std::vector<double> >("dRZMaxLowEtEtaBins")),
380  dRZMaxLowEt_(pset.getParameter<std::vector<double> >("dRZMaxLowEt")) {
381  if (dRZMaxLowEtEtaBins_.size() + 1 != dRZMaxLowEt_.size()) {
382  throw cms::Exception("InvalidConfig") << " dRZMaxLowEtEtaBins should be 1 less than dRZMaxLowEt when its "
383  << dRZMaxLowEtEtaBins_.size() << " vs " << dRZMaxLowEt_.size();
384  }
385 }
386 
388  if (dPhiMax_ >= 0 && std::abs(scHitMatch.dPhi) > dPhiMax_)
389  return false;
390 
391  const float dRZMax = getDRZCutValue(scHitMatch.et, scHitMatch.eta);
392  if (dRZMax_ >= 0 && std::abs(scHitMatch.dRZ) > dRZMax)
393  return false;
394 
395  return true;
396 }
397 
398 float TrajSeedMatcher::MatchingCutsV1::getDRZCutValue(const float scEt, const float scEta) const {
399  if (scEt >= dRZMaxLowEtThres_)
400  return dRZMax_;
401  else {
402  const float absEta = std::abs(scEta);
403  for (size_t etaNr = 0; etaNr < dRZMaxLowEtEtaBins_.size(); etaNr++) {
404  if (absEta < dRZMaxLowEtEtaBins_[etaNr])
405  return dRZMaxLowEt_[etaNr];
406  }
407  return dRZMaxLowEt_.back();
408  }
409 }
410 
412  : dPhiHighEt_(pset.getParameter<std::vector<double> >("dPhiMaxHighEt")),
413  dPhiHighEtThres_(pset.getParameter<std::vector<double> >("dPhiMaxHighEtThres")),
414  dPhiLowEtGrad_(pset.getParameter<std::vector<double> >("dPhiMaxLowEtGrad")),
415  dRZHighEt_(pset.getParameter<std::vector<double> >("dRZMaxHighEt")),
416  dRZHighEtThres_(pset.getParameter<std::vector<double> >("dRZMaxHighEtThres")),
417  dRZLowEtGrad_(pset.getParameter<std::vector<double> >("dRZMaxLowEtGrad")),
418  etaBins_(pset.getParameter<std::vector<double> >("etaBins")) {
419  auto binSizeCheck = [](size_t sizeEtaBins, const std::vector<double>& vec, const std::string& name) {
420  if (vec.size() != sizeEtaBins + 1) {
421  throw cms::Exception("InvalidConfig")
422  << " when constructing TrajSeedMatcher::MatchingCutsV2 " << name << " has " << vec.size()
423  << " bins, it should be equal to #bins of etaBins+1" << sizeEtaBins + 1;
424  }
425  };
426  binSizeCheck(etaBins_.size(), dPhiHighEt_, "dPhiMaxHighEt");
427  binSizeCheck(etaBins_.size(), dPhiHighEtThres_, "dPhiMaxHighEtThres");
428  binSizeCheck(etaBins_.size(), dPhiLowEtGrad_, "dPhiMaxLowEtGrad");
429  binSizeCheck(etaBins_.size(), dRZHighEt_, "dRZMaxHighEt");
430  binSizeCheck(etaBins_.size(), dRZHighEtThres_, "dRZMaxHighEtThres");
431  binSizeCheck(etaBins_.size(), dRZLowEtGrad_, "dRZMaxLowEtGrad");
432 }
433 
435  size_t binNr = getBinNr(scHitMatch.eta);
436  float dPhiMax = getCutValue(scHitMatch.et, dPhiHighEt_[binNr], dPhiHighEtThres_[binNr], dPhiLowEtGrad_[binNr]);
437  if (dPhiMax >= 0 && std::abs(scHitMatch.dPhi) > dPhiMax)
438  return false;
439  float dRZMax = getCutValue(scHitMatch.et, dRZHighEt_[binNr], dRZHighEtThres_[binNr], dRZLowEtGrad_[binNr]);
440  if (dRZMax >= 0 && std::abs(scHitMatch.dRZ) > dRZMax)
441  return false;
442 
443  return true;
444 }
445 
446 //eta bins is exactly 1 smaller than the vectors which will be accessed by this bin nr
448  const float absEta = std::abs(eta);
449  for (size_t etaNr = 0; etaNr < etaBins_.size(); etaNr++) {
450  if (absEta < etaBins_[etaNr])
451  return etaNr;
452  }
453  return etaBins_.size();
454 }
std::vector< SCHitMatch > processSeed(const TrajectorySeed &seed, const GlobalPoint &candPos, const float energy, const TrajectoryStateOnSurface &initialTrajState)
MatchingCutsV1(const edm::ParameterSet &pset)
size_t getNrHitsRequired(const int nrValidLayers) const
IntGlobalPointPairUnorderedMap< TrajectoryStateOnSurface > trajStateFromPointNegChargeCache_
virtual const DetLayer * idToLayer(const DetId &detId) const
MatchingCutsV2(const edm::ParameterSet &pset)
float getDRZCutValue(const float scEt, const float scEta) const
std::vector< double > dRZLowEtGrad_
const TrajectoryStateOnSurface & getTrajStateFromVtx(const TrackingRecHit &hit, const TrajectoryStateOnSurface &initialState, const PropagatorWithMaterial &propagator)
Configuration(const edm::ParameterSet &pset, edm::ConsumesCollector &&cc)
std::vector< double > dPhiLowEtGrad_
T z() const
Definition: PV3DBase.h:61
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
IntGlobalPointPairUnorderedMap< TrajectoryStateOnSurface > trajStateFromPointPosChargeCache_
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
T eta() const
Definition: PV3DBase.h:73
const GlobalTrajectoryParameters & globalParameters() const
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
MeasurementDetWithData idToDet(const DetId &id) const
Previous MeasurementDetSystem interface.
size_t getBinNr(float eta) const
const std::vector< double > dRZMaxLowEt_
static float getZVtxFromExtrapolation(const GlobalPoint &primeVtxPos, const GlobalPoint &hitPos, const GlobalPoint &candPos)
Definition: Electron.h:6
const TrajectoryStateOnSurface & getTrajStateFromPoint(const TrackingRecHit &hit, const FreeTrajectoryState &initialState, const GlobalPoint &point, const PropagatorWithMaterial &propagator)
constexpr std::array< uint8_t, layerIndexSize > layer
const std::vector< std::unique_ptr< MatchingCuts > > matchingCuts
std::vector< TrajSeedMatcher::SeedWithInfo > operator()(const GlobalPoint &candPos, const float energy)
std::vector< double > etaBins_
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
TrajectorySeedCollection const & seeds_
PropagatorWithMaterial forwardPropagator_
Configuration const & cfg_
static edm::ParameterSetDescription makePSetDescription()
bool operator()(const SCHitMatch &scHitMatch) const override
std::vector< const DetLayer * > compatibleLayers(const DetLayer &detLayer, Args &&... args) const
Returns all layers compatible.
const std::vector< int > minNrHitsValidLayerBins
std::vector< TrajectorySeed > TrajectorySeedCollection
T sqrt(T t)
Definition: SSEVec.h:19
TrackCharge charge() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
NavigationSchool const & navSchool_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::unordered_map< int, TrajectoryStateOnSurface > trajStateFromVtxPosChargeCache_
int getNrValidLayersAlongTraj(const SCHitMatch &hit1, const SCHitMatch &hit2, const GlobalPoint &candPos, const float energy, const int charge)
std::vector< double > dPhiHighEt_
DetLayerGeometry const & detLayerGeom_
bool operator()(const SCHitMatch &scHitMatch) const override
const TrackingRecHit & hit
TrajSeedMatcher(TrajectorySeedCollection const &seeds, math::XYZPoint const &vprim, Configuration const &cfg, edm::EventSetup const &iSetup, MeasurementTrackerEvent const &measTkEvt)
Definition: DetId.h:17
const std::vector< double > dRZMaxLowEtEtaBins_
DetId geographicalId() const
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
MeasurementTrackerEvent const & measTkEvt_
const GlobalPoint vprim_
std::vector< double > dPhiHighEtThres_
auto ftsFromVertexToPoint(GlobalPoint const &point, GlobalPoint const &vertex, float energy, int charge) const
PropagatorWithMaterial backwardPropagator_
FreeTrajectoryState const * freeState(bool withErrors=true) const
static constexpr float kElectronMass_
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
std::unordered_map< int, TrajectoryStateOnSurface > trajStateFromVtxNegChargeCache_
bool layerHasValidHits(const DetLayer &layer, const TrajectoryStateOnSurface &hitSurState, const Propagator &propToLayerFromState) const
const std::vector< unsigned int > minNrHits
def move(src, dest)
Definition: eostools.py:511
bool isTrackerPixel(GeomDetEnumerators::SubDetector m)
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
std::vector< double > dRZHighEtThres_
std::vector< double > dRZHighEt_
TrajectoryStateOnSurface makeTrajStateOnSurface(const GlobalPoint &pos, const float energy, const int charge) const