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 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 constexpr unsigned int minLayers = 4
 

Private Attributes

const float caHardPtCut = 0.f
 
CACut caPhiCut
 
CACut caThetaCut
 
const float extraHitRPhitolerance
 
const bool fitFastCircle
 
const bool fitFastCircleChi2Cut
 
const QuantityDependsPt maxChi2
 
std::unique_ptr< SeedComparitortheComparitor
 
const MagneticFieldtheField = nullptr
 
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecordtheFieldToken
 
LayerCacheType theLayerCache
 
const bool useBendingCorrection
 

Detailed Description

Definition at line 31 of file CAHitQuadrupletGenerator.h.

Member Typedef Documentation

◆ LayerCacheType

Definition at line 33 of file CAHitQuadrupletGenerator.h.

◆ ResultType

Definition at line 36 of file CAHitQuadrupletGenerator.h.

Constructor & Destructor Documentation

◆ CAHitQuadrupletGenerator() [1/2]

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

Definition at line 39 of file CAHitQuadrupletGenerator.h.

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

◆ CAHitQuadrupletGenerator() [2/2]

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

Definition at line 27 of file CAHitQuadrupletGenerator.cc.

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

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

◆ ~CAHitQuadrupletGenerator()

CAHitQuadrupletGenerator::~CAHitQuadrupletGenerator ( )
default

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 47 of file CAHitQuadrupletGenerator.cc.

References edm::ParameterSetDescription::add(), edm::ParameterSet::addParameter(), submitPVResolutionJobs::desc, and AlCaHLTBitMon_QueryRunRegistry::string.

47  {
48  desc.add<double>("extraHitRPhitolerance", 0.1);
49  desc.add<bool>("fitFastCircle", false);
50  desc.add<bool>("fitFastCircleChi2Cut", false);
51  desc.add<bool>("useBendingCorrection", false);
52  desc.add<double>("CAThetaCut", 0.00125);
53  desc.add<double>("CAPhiCut", 10);
54 
55  edm::ParameterSetDescription validatorCACut;
56  validatorCACut.add<string>("seedingLayers", "BPix1+BPix2+BPix3");
57  validatorCACut.add<double>("cut", 0.00125);
58  std::vector<edm::ParameterSet> defaultCACutVector;
59  edm::ParameterSet defaultCACut;
60  defaultCACut.addParameter<string>("seedingLayers", "");
61  defaultCACut.addParameter<double>("cut", -1.);
62  defaultCACutVector.push_back(defaultCACut);
63  desc.addVPSet("CAThetaCut_byTriplets", validatorCACut, defaultCACutVector);
64  desc.addVPSet("CAPhiCut_byTriplets", validatorCACut, defaultCACutVector);
65 
66  desc.add<double>("CAHardPtCut", 0);
67  desc.addOptional<bool>("CAOnlyOneLastHitPerLayerFilter")
68  ->setComment(
69  "Deprecated and has no effect. To be fully removed later when the parameter is no longer used in HLT "
70  "configurations.");
71  edm::ParameterSetDescription descMaxChi2;
72  descMaxChi2.add<double>("pt1", 0.2);
73  descMaxChi2.add<double>("pt2", 1.5);
74  descMaxChi2.add<double>("value1", 500);
75  descMaxChi2.add<double>("value2", 50);
76  descMaxChi2.add<bool>("enabled", true);
77  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
78 
79  edm::ParameterSetDescription descComparitor;
80  descComparitor.add<std::string>("ComponentName", "none");
81  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
82  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
83 }
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:136
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ fillDescriptionsLabel()

static const char* CAHitQuadrupletGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 46 of file CAHitQuadrupletGenerator.h.

46 { return "caHitQuadrupletDefault"; }

◆ hitNtuplets()

void CAHitQuadrupletGenerator::hitNtuplets ( const IntermediateHitDoublets regionDoublets,
std::vector< OrderedHitSeeds > &  result,
const SeedingLayerSetsHits layers 
)

Definition at line 162 of file CAHitQuadrupletGenerator.cc.

References funct::abs(), DummyCfis::c, caHardPtCut, caPhiCut, caThetaCut, isoTrack_cff::chi2, RZLine::chi2(), CellularAutomaton::createAndConnectCells(), PixelRecoUtilities::curvature(), relativeConstraints::error, CAHitQuadrupletGenerator::QuantityDependsPt::evaluator(), CellularAutomaton::evolve(), extraHitRPhitolerance, CellularAutomaton::findNtuplets(), fitFastCircle, fitFastCircleChi2Cut, g, CellularAutomaton::getAllCells(), mps_fire::i, EcalPhiSymFlatTableProducers_cfi::id, PixelRecoUtilities::inversePt(), PixelPluginsPhase0_cfi::isBarrel, edm::isNotFinite(), hgcalTBTopologyTester_cfi::layers, maxChi2, PixelSubdetector::PixelBarrel, point, DiDispStaMuonMonitor_cfi::pt, nano_mu_digi_cff::region, ihd::RegionLayerSets< T >::region(), mps_fire::result, CACut::setCutValuesByLayerIds(), funct::sqr(), mathSSE::sqrt(), theComparitor, theField, useBendingCorrection, CAHitQuadrupletGenerator::QuantityDependsPtEval::value(), x, and y.

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

◆ initEvent()

void CAHitQuadrupletGenerator::initEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Definition at line 85 of file CAHitQuadrupletGenerator.cc.

References makeMEIFBenchmarkPlots::ev, edm::EventSetup::getData(), theComparitor, theField, and theFieldToken.

85  {
86  if (theComparitor)
87  theComparitor->init(ev, es);
89 }
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
std::unique_ptr< SeedComparitor > theComparitor
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > theFieldToken
const MagneticField * theField

Member Data Documentation

◆ caHardPtCut

const float CAHitQuadrupletGenerator::caHardPtCut = 0.f
private

Definition at line 136 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

◆ caPhiCut

CACut CAHitQuadrupletGenerator::caPhiCut
private

Definition at line 135 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

◆ caThetaCut

CACut CAHitQuadrupletGenerator::caThetaCut
private

Definition at line 134 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

◆ extraHitRPhitolerance

const float CAHitQuadrupletGenerator::extraHitRPhitolerance
private

Definition at line 127 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

◆ fitFastCircle

const bool CAHitQuadrupletGenerator::fitFastCircle
private

Definition at line 130 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

◆ fitFastCircleChi2Cut

const bool CAHitQuadrupletGenerator::fitFastCircleChi2Cut
private

Definition at line 131 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

◆ maxChi2

const QuantityDependsPt CAHitQuadrupletGenerator::maxChi2
private

Definition at line 129 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

◆ minLayers

constexpr unsigned int CAHitQuadrupletGenerator::minLayers = 4
static

Definition at line 35 of file CAHitQuadrupletGenerator.h.

◆ theComparitor

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

Definition at line 57 of file CAHitQuadrupletGenerator.h.

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

◆ theField

const MagneticField* CAHitQuadrupletGenerator::theField = nullptr
private

Definition at line 60 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets(), and initEvent().

◆ theFieldToken

const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> CAHitQuadrupletGenerator::theFieldToken
private

Definition at line 59 of file CAHitQuadrupletGenerator.h.

Referenced by initEvent().

◆ theLayerCache

LayerCacheType CAHitQuadrupletGenerator::theLayerCache
private

Definition at line 55 of file CAHitQuadrupletGenerator.h.

◆ useBendingCorrection

const bool CAHitQuadrupletGenerator::useBendingCorrection
private

Definition at line 132 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().