CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes
CAHitQuadrupletGenerator Class Reference

#include <CAHitQuadrupletGenerator.h>

Classes

class  QuantityDependsPt
 
class  QuantityDependsPtEval
 

Public Types

typedef LayerHitMapCache LayerCacheType
 
typedef OrderedHitSeeds ResultType
 

Public Member Functions

 CAHitQuadrupletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
 CAHitQuadrupletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
void hitNtuplets (const IntermediateHitDoublets &regionDoublets, std::vector< OrderedHitSeeds > &result, const edm::EventSetup &es, const SeedingLayerSetsHits &layers)
 
void initEvent (const edm::Event &ev, const edm::EventSetup &es)
 
 ~CAHitQuadrupletGenerator ()=default
 

Static Public Member Functions

static void fillDescriptions (edm::ParameterSetDescription &desc)
 
static const char * fillDescriptionsLabel ()
 

Static Public Attributes

static unsigned int minLayers = 4
 

Private Attributes

const float caHardPtCut = 0.f
 
const float caPhiCut = 0.1f
 
const float caThetaCut = 0.00125f
 
const float extraHitRPhitolerance
 
const bool fitFastCircle
 
const bool fitFastCircleChi2Cut
 
const QuantityDependsPt maxChi2
 
std::unique_ptr< SeedComparitortheComparitor
 
LayerCacheType theLayerCache
 
const bool useBendingCorrection
 

Detailed Description

Definition at line 29 of file CAHitQuadrupletGenerator.h.

Member Typedef Documentation

Definition at line 31 of file CAHitQuadrupletGenerator.h.

Definition at line 34 of file CAHitQuadrupletGenerator.h.

Constructor & Destructor Documentation

CAHitQuadrupletGenerator::CAHitQuadrupletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 38 of file CAHitQuadrupletGenerator.h.

References looper::cfg, and fillDescriptions().

38 : CAHitQuadrupletGenerator(cfg, iC) {}
CAHitQuadrupletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
CAHitQuadrupletGenerator::CAHitQuadrupletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 27 of file CAHitQuadrupletGenerator.cc.

References timingPdfMaker::get, edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, and theComparitor.

28  : extraHitRPhitolerance(cfg.getParameter<double>(
29  "extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
30  maxChi2(cfg.getParameter<edm::ParameterSet>("maxChi2")),
31  fitFastCircle(cfg.getParameter<bool>("fitFastCircle")),
32  fitFastCircleChi2Cut(cfg.getParameter<bool>("fitFastCircleChi2Cut")),
33  useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
34  caThetaCut(cfg.getParameter<double>("CAThetaCut")),
35  caPhiCut(cfg.getParameter<double>("CAPhiCut")),
36  caHardPtCut(cfg.getParameter<double>("CAHardPtCut")) {
37  edm::ParameterSet comparitorPSet = cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
38  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
39  if (comparitorName != "none") {
40  theComparitor = SeedComparitorFactory::get()->create(comparitorName, comparitorPSet, iC);
41  }
42 }
T getParameter(std::string const &) const
std::unique_ptr< SeedComparitor > theComparitor
const QuantityDependsPt maxChi2
CAHitQuadrupletGenerator::~CAHitQuadrupletGenerator ( )
default

Member Function Documentation

void CAHitQuadrupletGenerator::fillDescriptions ( edm::ParameterSetDescription desc)
static

Definition at line 44 of file CAHitQuadrupletGenerator.cc.

References edm::ParameterSetDescription::add(), edm::ParameterSetDescription::addOptional(), and AlCaHLTBitMon_QueryRunRegistry::string.

44  {
45  desc.add<double>("extraHitRPhitolerance", 0.1);
46  desc.add<bool>("fitFastCircle", false);
47  desc.add<bool>("fitFastCircleChi2Cut", false);
48  desc.add<bool>("useBendingCorrection", false);
49  desc.add<double>("CAThetaCut", 0.00125);
50  desc.add<double>("CAPhiCut", 10);
51  desc.add<double>("CAHardPtCut", 0);
52  desc.addOptional<bool>("CAOnlyOneLastHitPerLayerFilter")
53  ->setComment(
54  "Deprecated and has no effect. To be fully removed later when the parameter is no longer used in HLT "
55  "configurations.");
56  edm::ParameterSetDescription descMaxChi2;
57  descMaxChi2.add<double>("pt1", 0.2);
58  descMaxChi2.add<double>("pt2", 1.5);
59  descMaxChi2.add<double>("value1", 500);
60  descMaxChi2.add<double>("value2", 50);
61  descMaxChi2.add<bool>("enabled", true);
62  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
63 
64  edm::ParameterSetDescription descComparitor;
65  descComparitor.add<std::string>("ComponentName", "none");
66  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
67  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
68 }
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static const char* CAHitQuadrupletGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 44 of file CAHitQuadrupletGenerator.h.

References ev, hgcalTopologyTester_cfi::layers, and mps_fire::result.

44 { return "caHitQuadruplet"; }
void CAHitQuadrupletGenerator::hitNtuplets ( const IntermediateHitDoublets regionDoublets,
std::vector< OrderedHitSeeds > &  result,
const edm::EventSetup es,
const SeedingLayerSetsHits layers 
)

Definition at line 146 of file CAHitQuadrupletGenerator.cc.

References funct::abs(), HltBtagPostValidation_cff::c, caHardPtCut, caPhiCut, caThetaCut, hltPixelTracks_cff::chi2, FastCircleFit::chi2(), RZLine::chi2(), CellularAutomaton::createAndConnectCells(), PixelRecoUtilities::curvature(), GlobalErrorBase< T, ErrorWeightType >::czz(), relativeConstraints::error, CAHitQuadrupletGenerator::QuantityDependsPt::evaluator(), CellularAutomaton::evolve(), extraHitRPhitolerance, CellularAutomaton::findNtuplets(), fitFastCircle, fitFastCircleChi2Cut, g, CellularAutomaton::getAllCells(), mps_fire::i, triggerObjects_cff::id, PixelRecoUtilities::inversePt(), PixelPluginsPhase0_cfi::isBarrel, edm::isNotFinite(), maxChi2, TrackingRegion::origin(), PixelSubdetector::PixelBarrel, point, DiDispStaMuonMonitor_cfi::pt, HLT_2018_cff::region, GlobalErrorBase< T, ErrorWeightType >::rerr(), funct::sqr(), mathSSE::sqrt(), theComparitor, useBendingCorrection, CAHitQuadrupletGenerator::QuantityDependsPtEval::value(), x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

149  {
150  CAGraph g;
151 
152  std::vector<const HitDoublets*> hitDoublets;
153 
154  const int numberOfHitsInNtuplet = 4;
155  std::vector<CACell::CAntuplet> foundQuadruplets;
156 
157  int index = 0;
158  for (const auto& regionLayerPairs : regionDoublets) {
159  const TrackingRegion& region = regionLayerPairs.region();
160  hitDoublets.clear();
161  foundQuadruplets.clear();
162  if (index == 0) {
163  createGraphStructure(layers, g);
164  } else {
165  clearGraphStructure(layers, g);
166  }
167 
168  fillGraph(layers, regionLayerPairs, g, hitDoublets);
169 
170  CellularAutomaton ca(g);
171 
172  ca.createAndConnectCells(hitDoublets, region, caThetaCut, caPhiCut, caHardPtCut);
173 
174  ca.evolve(numberOfHitsInNtuplet);
175 
176  ca.findNtuplets(foundQuadruplets, numberOfHitsInNtuplet);
177 
178  auto& allCells = ca.getAllCells();
179 
180  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
181 
182  // re-used thoughout
183  std::array<float, 4> bc_r;
184  std::array<float, 4> bc_z;
185  std::array<float, 4> bc_errZ2;
186  std::array<GlobalPoint, 4> gps;
187  std::array<GlobalError, 4> ges;
188  std::array<bool, 4> barrels;
189 
190  unsigned int numberOfFoundQuadruplets = foundQuadruplets.size();
191 
192  // Loop over quadruplets
193  for (unsigned int quadId = 0; quadId < numberOfFoundQuadruplets; ++quadId) {
194  auto isBarrel = [](const unsigned id) -> bool { return id == PixelSubdetector::PixelBarrel; };
195  for (unsigned int i = 0; i < 3; ++i) {
196  auto const& ahit = allCells[foundQuadruplets[quadId][i]].getInnerHit();
197  gps[i] = ahit->globalPosition();
198  ges[i] = ahit->globalPositionError();
199  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
200  }
201 
202  auto const& ahit = allCells[foundQuadruplets[quadId][2]].getOuterHit();
203  gps[3] = ahit->globalPosition();
204  ges[3] = ahit->globalPositionError();
205  barrels[3] = isBarrel(ahit->geographicalId().subdetId());
206  // TODO:
207  // - if we decide to always do the circle fit for 4 hits, we don't
208  // need ThirdHitPredictionFromCircle for the curvature; then we
209  // could remove extraHitRPhitolerance configuration parameter
210  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2], extraHitRPhitolerance);
211  const float curvature = predictionRPhi.curvature(ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
212  const float abscurv = std::abs(curvature);
213  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
214  if (theComparitor) {
215  SeedingHitSet tmpTriplet(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
216  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
217  allCells[foundQuadruplets[quadId][2]].getOuterHit());
218 
219  if (!theComparitor->compatible(tmpTriplet)) {
220  continue;
221  }
222  }
223 
224  float chi2 = std::numeric_limits<float>::quiet_NaN();
225  // TODO: Do we have any use case to not use bending correction?
226  if (useBendingCorrection) {
227  // Following PixelFitterByConformalMappingAndLine
228  const float simpleCot = (gps.back().z() - gps.front().z()) / (gps.back().perp() - gps.front().perp());
229  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
230  for (int i = 0; i < 4; ++i) {
231  const GlobalPoint& point = gps[i];
232  const GlobalError& error = ges[i];
233  bc_r[i] = sqrt(sqr(point.x() - region.origin().x()) + sqr(point.y() - region.origin().y()));
234  bc_r[i] += pixelrecoutilities::LongitudinalBendingCorrection(pt, es)(bc_r[i]);
235  bc_z[i] = point.z() - region.origin().z();
236  bc_errZ2[i] = (barrels[i]) ? error.czz() : error.rerr(point) * sqr(simpleCot);
237  }
238  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
239  chi2 = rzLine.chi2();
240  } else {
241  RZLine rzLine(gps, ges, barrels);
242  chi2 = rzLine.chi2();
243  }
244  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2) {
245  continue;
246  }
247  // TODO: Do we have any use case to not use circle fit? Maybe
248  // HLT where low-pT inefficiency is not a problem?
249  if (fitFastCircle) {
250  FastCircleFit c(gps, ges);
251  chi2 += c.chi2();
252  if (edm::isNotFinite(chi2))
253  continue;
254  if (fitFastCircleChi2Cut && chi2 > thisMaxChi2)
255  continue;
256  }
257  result[index].emplace_back(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
258  allCells[foundQuadruplets[quadId][1]].getInnerHit(),
259  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
260  allCells[foundQuadruplets[quadId][2]].getOuterHit());
261  }
262  index++;
263  }
264 }
GlobalPoint const & origin() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
T y() const
Definition: PV3DBase.h:60
std::unique_ptr< SeedComparitor > theComparitor
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 g
Definition: Activities.doc:4
T inversePt(T curvature, const edm::EventSetup &iSetup)
T curvature(T InversePt, const edm::EventSetup &iSetup)
T sqrt(T t)
Definition: SSEVec.h:19
T z() const
Definition: PV3DBase.h:61
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: RZLine.h:12
const QuantityDependsPt maxChi2
T rerr(const GlobalPoint &aPoint) const
Square< F >::type sqr(const F &f)
Definition: Square.h:14
QuantityDependsPtEval evaluator(const edm::EventSetup &es) const
T x() const
Definition: PV3DBase.h:59
*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
void CAHitQuadrupletGenerator::initEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Member Data Documentation

const float CAHitQuadrupletGenerator::caHardPtCut = 0.f
private

Definition at line 131 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::caPhiCut = 0.1f
private

Definition at line 130 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::caThetaCut = 0.00125f
private

Definition at line 129 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::extraHitRPhitolerance
private

Definition at line 122 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const bool CAHitQuadrupletGenerator::fitFastCircle
private

Definition at line 125 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const bool CAHitQuadrupletGenerator::fitFastCircleChi2Cut
private

Definition at line 126 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const QuantityDependsPt CAHitQuadrupletGenerator::maxChi2
private

Definition at line 124 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

unsigned int CAHitQuadrupletGenerator::minLayers = 4
static

Definition at line 33 of file CAHitQuadrupletGenerator.h.

std::unique_ptr<SeedComparitor> CAHitQuadrupletGenerator::theComparitor
private

Definition at line 56 of file CAHitQuadrupletGenerator.h.

Referenced by CAHitQuadrupletGenerator(), hitNtuplets(), and initEvent().

LayerCacheType CAHitQuadrupletGenerator::theLayerCache
private

Definition at line 54 of file CAHitQuadrupletGenerator.h.

const bool CAHitQuadrupletGenerator::useBendingCorrection
private

Definition at line 127 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().