CMS 3D CMS Logo

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

#include <CAHitTripletGenerator.h>

Classes

class  QuantityDependsPt
 
class  QuantityDependsPtEval
 

Public Types

typedef LayerHitMapCache LayerCacheType
 
typedef OrderedHitSeeds ResultType
 

Public Member Functions

 CAHitTripletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
 CAHitTripletGenerator (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)
 
 ~CAHitTripletGenerator ()=default
 

Static Public Member Functions

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

Static Public Attributes

static constexpr unsigned int minLayers = 3
 

Private Attributes

const float caHardPtCut = 0.f
 
CACut caPhiCut
 
CACut caThetaCut
 
const float extraHitRPhitolerance
 
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 CAHitTripletGenerator.h.

Member Typedef Documentation

◆ LayerCacheType

Definition at line 33 of file CAHitTripletGenerator.h.

◆ ResultType

Definition at line 36 of file CAHitTripletGenerator.h.

Constructor & Destructor Documentation

◆ CAHitTripletGenerator() [1/2]

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

Definition at line 39 of file CAHitTripletGenerator.h.

39 : CAHitTripletGenerator(cfg, iC) {}
CAHitTripletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)

◆ CAHitTripletGenerator() [2/2]

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

Definition at line 28 of file CAHitTripletGenerator.cc.

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

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

◆ ~CAHitTripletGenerator()

CAHitTripletGenerator::~CAHitTripletGenerator ( )
default

Member Function Documentation

◆ fillDescriptions()

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

Definition at line 46 of file CAHitTripletGenerator.cc.

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

46  {
47  desc.add<double>("extraHitRPhitolerance", 0.06);
48  desc.add<bool>("useBendingCorrection", false);
49  desc.add<double>("CAThetaCut", 0.00125);
50  desc.add<double>("CAPhiCut", 0.1);
51 
52  edm::ParameterSetDescription validatorCACut;
53  validatorCACut.add<string>("seedingLayers", "BPix1+BPix2+BPix3");
54  validatorCACut.add<double>("cut", 0.00125);
55  std::vector<edm::ParameterSet> defaultCACutVector;
56  edm::ParameterSet defaultCACut;
57  defaultCACut.addParameter<string>("seedingLayers", "");
58  defaultCACut.addParameter<double>("cut", -1.);
59  defaultCACutVector.push_back(defaultCACut);
60  desc.addVPSet("CAThetaCut_byTriplets", validatorCACut, defaultCACutVector);
61  desc.addVPSet("CAPhiCut_byTriplets", validatorCACut, defaultCACutVector);
62 
63  desc.add<double>("CAHardPtCut", 0);
64 
65  edm::ParameterSetDescription descMaxChi2;
66  descMaxChi2.add<double>("pt1", 0.8);
67  descMaxChi2.add<double>("pt2", 2);
68  descMaxChi2.add<double>("value1", 50);
69  descMaxChi2.add<double>("value2", 8);
70  descMaxChi2.add<bool>("enabled", true);
71  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
72 
73  edm::ParameterSetDescription descComparitor;
74  descComparitor.add<std::string>("ComponentName", "none");
75  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
76  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
77 }
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
ParameterDescriptionBase * add(U const &iLabel, T const &value)

◆ fillDescriptionsLabel()

static const char* CAHitTripletGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 45 of file CAHitTripletGenerator.h.

45 { return "caHitTriplet"; }

◆ hitNtuplets()

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

Definition at line 160 of file CAHitTripletGenerator.cc.

References funct::abs(), caHardPtCut, caPhiCut, caThetaCut, hltPixelTracks_cff::chi2, RZLine::chi2(), PixelRecoUtilities::curvature(), ThirdHitPredictionFromCircle::curvature(), relativeConstraints::error, CAHitTripletGenerator::QuantityDependsPt::evaluator(), extraHitRPhitolerance, CellularAutomaton::findTriplets(), g, CellularAutomaton::getAllCells(), mps_fire::i, l1ctLayer2EG_cff::id, PixelRecoUtilities::inversePt(), PixelPluginsPhase0_cfi::isBarrel, edm::isNotFinite(), hgcalTopologyTester_cfi::layers, mps_splice::line, maxChi2, PixelSubdetector::PixelBarrel, point, DiDispStaMuonMonitor_cfi::pt, ihd::RegionLayerSets< T >::region(), HLT_2022v14_cff::region, mps_fire::result, CACut::setCutValuesByLayerIds(), funct::sqr(), mathSSE::sqrt(), theComparitor, theField, useBendingCorrection, CAHitTripletGenerator::QuantityDependsPtEval::value(), x, and y.

162  {
163  CAGraph g;
164 
165  std::vector<const HitDoublets*> hitDoublets;
166 
167  std::vector<CACell::CAntuplet> foundTriplets;
168 
169  int index = 0;
170  for (const auto& regionLayerPairs : regionDoublets) {
171  const TrackingRegion& region = regionLayerPairs.region();
172  hitDoublets.clear();
173  foundTriplets.clear();
174 
175  if (index == 0) {
176  createGraphStructure(layers, g);
179  } else {
180  clearGraphStructure(layers, g);
181  }
182  fillGraph(layers, regionLayerPairs, g, hitDoublets);
183  CellularAutomaton ca(g);
184  ca.findTriplets(hitDoublets, foundTriplets, region, caThetaCut, caPhiCut, caHardPtCut);
185 
186  auto& allCells = ca.getAllCells();
187 
188  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(*theField);
189 
190  // re-used thoughout, need to be vectors because of RZLine interface
191  std::array<float, 3> bc_r;
192  std::array<float, 3> bc_z;
193  std::array<float, 3> bc_errZ2;
194  std::array<GlobalPoint, 3> gps;
195  std::array<GlobalError, 3> ges;
196  std::array<bool, 3> barrels;
197 
198  unsigned int numberOfFoundTriplets = foundTriplets.size();
199  for (unsigned int tripletId = 0; tripletId < numberOfFoundTriplets; ++tripletId) {
200  OrderedHitTriplet tmpTriplet(allCells[foundTriplets[tripletId][0]].getInnerHit(),
201  allCells[foundTriplets[tripletId][0]].getOuterHit(),
202  allCells[foundTriplets[tripletId][1]].getOuterHit());
203 
204  auto isBarrel = [](const unsigned id) -> bool { return id == PixelSubdetector::PixelBarrel; };
205  for (unsigned int i = 0; i < 2; ++i) {
206  auto const& ahit = allCells[foundTriplets[tripletId][i]].getInnerHit();
207  gps[i] = ahit->globalPosition();
208  ges[i] = ahit->globalPositionError();
209  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
210  }
211 
212  auto const& ahit = allCells[foundTriplets[tripletId][1]].getOuterHit();
213  gps[2] = ahit->globalPosition();
214  ges[2] = ahit->globalPositionError();
215  barrels[2] = isBarrel(ahit->geographicalId().subdetId());
216 
217  PixelRecoLineRZ line(gps[0], gps[2]);
218  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2], extraHitRPhitolerance);
219  const float curvature = predictionRPhi.curvature(ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
220  const float abscurv = std::abs(curvature);
221  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
222  float chi2 = std::numeric_limits<float>::quiet_NaN();
223  // TODO: Do we have any use case to not use bending correction?
224 
225  if (useBendingCorrection) {
226  // Following PixelFitterByConformalMappingAndLine
227  const float simpleCot = (gps.back().z() - gps.front().z()) / (gps.back().perp() - gps.front().perp());
228  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, *theField);
229  for (int i = 0; i < 3; ++i) {
230  const GlobalPoint& point = gps[i];
231  const GlobalError& error = ges[i];
232  bc_r[i] = sqrt(sqr(point.x() - region.origin().x()) + sqr(point.y() - region.origin().y()));
234  bc_z[i] = point.z() - region.origin().z();
235  bc_errZ2[i] = (barrels[i]) ? error.czz() : error.rerr(point) * sqr(simpleCot);
236  }
237  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
238  chi2 = rzLine.chi2();
239  } else {
240  RZLine rzLine(gps, ges, barrels);
241  chi2 = rzLine.chi2();
242  }
243 
244  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2) {
245  continue;
246  }
247 
248  if (theComparitor) {
249  if (!theComparitor->compatible(tmpTriplet)) {
250  continue;
251  }
252  }
253  result[index].emplace_back(tmpTriplet);
254  }
255  index++;
256  }
257 }
const QuantityDependsPt maxChi2
void setCutValuesByLayerIds(CAGraph &caLayers)
Definition: CACut.h:43
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
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)
const MagneticField * theField
QuantityDependsPtEval evaluator(const MagneticField &field) const
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T inversePt(T curvature, const MagneticField &field)
Definition: RZLine.h:12
Square< F >::type sqr(const F &f)
Definition: Square.h:14
std::unique_ptr< SeedComparitor > theComparitor
*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 CAHitTripletGenerator::initEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Definition at line 79 of file CAHitTripletGenerator.cc.

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

79  {
80  if (theComparitor)
81  theComparitor->init(ev, es);
83 }
const MagneticField * theField
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > theFieldToken
std::unique_ptr< SeedComparitor > theComparitor

Member Data Documentation

◆ caHardPtCut

const float CAHitTripletGenerator::caHardPtCut = 0.f
private

Definition at line 133 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

◆ caPhiCut

CACut CAHitTripletGenerator::caPhiCut
private

Definition at line 132 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

◆ caThetaCut

CACut CAHitTripletGenerator::caThetaCut
private

Definition at line 131 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

◆ extraHitRPhitolerance

const float CAHitTripletGenerator::extraHitRPhitolerance
private

Definition at line 126 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

◆ maxChi2

const QuantityDependsPt CAHitTripletGenerator::maxChi2
private

Definition at line 128 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

◆ minLayers

constexpr unsigned int CAHitTripletGenerator::minLayers = 3
static

Definition at line 35 of file CAHitTripletGenerator.h.

◆ theComparitor

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

Definition at line 56 of file CAHitTripletGenerator.h.

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

◆ theField

const MagneticField* CAHitTripletGenerator::theField = nullptr
private

Definition at line 59 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets(), and initEvent().

◆ theFieldToken

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

Definition at line 58 of file CAHitTripletGenerator.h.

Referenced by initEvent().

◆ theLayerCache

LayerCacheType CAHitTripletGenerator::theLayerCache
private

Definition at line 54 of file CAHitTripletGenerator.h.

◆ useBendingCorrection

const bool CAHitTripletGenerator::useBendingCorrection
private

Definition at line 129 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().