CMS 3D CMS Logo

CAHitQuadrupletGenerator.cc
Go to the documentation of this file.
1 #include <functional>
2 
12 
13 #include "CAGraph.h"
14 #include "CellularAutomaton.h"
15 
16 namespace {
17  template <typename T>
18  T sqr(T x) {
19  return x * x;
20  }
21 } // namespace
22 
23 using namespace std;
24 
25 constexpr unsigned int CAHitQuadrupletGenerator::minLayers;
26 
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  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 }
45 
47  desc.add<double>("extraHitRPhitolerance", 0.1);
48  desc.add<bool>("fitFastCircle", false);
49  desc.add<bool>("fitFastCircleChi2Cut", false);
50  desc.add<bool>("useBendingCorrection", false);
51  desc.add<double>("CAThetaCut", 0.00125);
52  desc.add<double>("CAPhiCut", 10);
53 
54  edm::ParameterSetDescription validatorCACut;
55  validatorCACut.add<string>("seedingLayers", "BPix1+BPix2+BPix3");
56  validatorCACut.add<double>("cut", 0.00125);
57  std::vector<edm::ParameterSet> defaultCACutVector;
58  edm::ParameterSet defaultCACut;
59  defaultCACut.addParameter<string>("seedingLayers", "");
60  defaultCACut.addParameter<double>("cut", -1.);
61  defaultCACutVector.push_back(defaultCACut);
62  desc.addVPSet("CAThetaCut_byTriplets", validatorCACut, defaultCACutVector);
63  desc.addVPSet("CAPhiCut_byTriplets", validatorCACut, defaultCACutVector);
64 
65  desc.add<double>("CAHardPtCut", 0);
66  desc.addOptional<bool>("CAOnlyOneLastHitPerLayerFilter")
67  ->setComment(
68  "Deprecated and has no effect. To be fully removed later when the parameter is no longer used in HLT "
69  "configurations.");
70  edm::ParameterSetDescription descMaxChi2;
71  descMaxChi2.add<double>("pt1", 0.2);
72  descMaxChi2.add<double>("pt2", 1.5);
73  descMaxChi2.add<double>("value1", 500);
74  descMaxChi2.add<double>("value2", 50);
75  descMaxChi2.add<bool>("enabled", true);
76  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
77 
78  edm::ParameterSetDescription descComparitor;
79  descComparitor.add<std::string>("ComponentName", "none");
80  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
81  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
82 }
83 
85  if (theComparitor)
86  theComparitor->init(ev, es);
87 }
88 namespace {
89  void createGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
90  for (unsigned int i = 0; i < layers.size(); i++) {
91  for (unsigned int j = 0; j < 4; ++j) {
92  auto vertexIndex = 0;
93  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j].name());
94  if (foundVertex == g.theLayers.end()) {
95  g.theLayers.emplace_back(layers[i][j].name(), layers[i][j].detLayer()->seqNum(), layers[i][j].hits().size());
96  vertexIndex = g.theLayers.size() - 1;
97  } else {
98  vertexIndex = foundVertex - g.theLayers.begin();
99  }
100  if (j == 0) {
101  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(), vertexIndex) == g.theRootLayers.end()) {
102  g.theRootLayers.emplace_back(vertexIndex);
103  }
104  }
105  }
106  }
107  }
108  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
109  g.theLayerPairs.clear();
110  for (unsigned int i = 0; i < g.theLayers.size(); i++) {
111  g.theLayers[i].theInnerLayers.clear();
112  g.theLayers[i].theInnerLayerPairs.clear();
113  g.theLayers[i].theOuterLayers.clear();
114  g.theLayers[i].theOuterLayerPairs.clear();
115  for (auto& v : g.theLayers[i].isOuterHitOfCell)
116  v.clear();
117  }
118  }
119  void fillGraph(const SeedingLayerSetsHits& layers,
120  const IntermediateHitDoublets::RegionLayerSets& regionLayerPairs,
121  CAGraph& g,
122  std::vector<const HitDoublets*>& hitDoublets) {
123  for (unsigned int i = 0; i < layers.size(); i++) {
124  for (unsigned int j = 0; j < 4; ++j) {
125  auto vertexIndex = 0;
126  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j].name());
127  if (foundVertex == g.theLayers.end()) {
128  vertexIndex = g.theLayers.size() - 1;
129  } else {
130  vertexIndex = foundVertex - g.theLayers.begin();
131  }
132 
133  if (j > 0) {
134  auto innerVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j - 1].name());
135 
136  CALayerPair tmpInnerLayerPair(innerVertex - g.theLayers.begin(), vertexIndex);
137 
138  if (std::find(g.theLayerPairs.begin(), g.theLayerPairs.end(), tmpInnerLayerPair) == g.theLayerPairs.end()) {
139  auto found = std::find_if(regionLayerPairs.begin(),
140  regionLayerPairs.end(),
142  return pair.innerLayerIndex() == layers[i][j - 1].index() &&
143  pair.outerLayerIndex() == layers[i][j].index();
144  });
145  if (found != regionLayerPairs.end()) {
146  hitDoublets.emplace_back(&(found->doublets()));
147  g.theLayerPairs.push_back(tmpInnerLayerPair);
148  g.theLayers[vertexIndex].theInnerLayers.push_back(innerVertex - g.theLayers.begin());
149  innerVertex->theOuterLayers.push_back(vertexIndex);
150  g.theLayers[vertexIndex].theInnerLayerPairs.push_back(g.theLayerPairs.size() - 1);
151  innerVertex->theOuterLayerPairs.push_back(g.theLayerPairs.size() - 1);
152  }
153  }
154  }
155  }
156  }
157  }
158 } // namespace
159 
161  std::vector<OrderedHitSeeds>& result,
162  const edm::EventSetup& es,
163  const SeedingLayerSetsHits& layers) {
164  CAGraph g;
165 
166  std::vector<const HitDoublets*> hitDoublets;
167 
168  const int numberOfHitsInNtuplet = 4;
169  std::vector<CACell::CAntuplet> foundQuadruplets;
170 
171  int index = 0;
172  for (const auto& regionLayerPairs : regionDoublets) {
173  const TrackingRegion& region = regionLayerPairs.region();
174  hitDoublets.clear();
175  foundQuadruplets.clear();
176  if (index == 0) {
177  createGraphStructure(layers, g);
180  } else {
181  clearGraphStructure(layers, g);
182  }
183 
184  fillGraph(layers, regionLayerPairs, g, hitDoublets);
185 
186  CellularAutomaton ca(g);
187 
189 
190  ca.evolve(numberOfHitsInNtuplet);
191 
192  ca.findNtuplets(foundQuadruplets, numberOfHitsInNtuplet);
193 
194  auto& allCells = ca.getAllCells();
195 
196  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
197 
198  // re-used thoughout
199  std::array<float, 4> bc_r;
200  std::array<float, 4> bc_z;
201  std::array<float, 4> bc_errZ2;
202  std::array<GlobalPoint, 4> gps;
203  std::array<GlobalError, 4> ges;
204  std::array<bool, 4> barrels;
205 
206  unsigned int numberOfFoundQuadruplets = foundQuadruplets.size();
207 
208  // Loop over quadruplets
209  for (unsigned int quadId = 0; quadId < numberOfFoundQuadruplets; ++quadId) {
210  auto isBarrel = [](const unsigned id) -> bool { return id == PixelSubdetector::PixelBarrel; };
211  for (unsigned int i = 0; i < 3; ++i) {
212  auto const& ahit = allCells[foundQuadruplets[quadId][i]].getInnerHit();
213  gps[i] = ahit->globalPosition();
214  ges[i] = ahit->globalPositionError();
215  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
216  }
217 
218  auto const& ahit = allCells[foundQuadruplets[quadId][2]].getOuterHit();
219  gps[3] = ahit->globalPosition();
220  ges[3] = ahit->globalPositionError();
221  barrels[3] = isBarrel(ahit->geographicalId().subdetId());
222  // TODO:
223  // - if we decide to always do the circle fit for 4 hits, we don't
224  // need ThirdHitPredictionFromCircle for the curvature; then we
225  // could remove extraHitRPhitolerance configuration parameter
226  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2], extraHitRPhitolerance);
227  const float curvature = predictionRPhi.curvature(ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
228  const float abscurv = std::abs(curvature);
229  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
230  if (theComparitor) {
231  SeedingHitSet tmpTriplet(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
232  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
233  allCells[foundQuadruplets[quadId][2]].getOuterHit());
234 
235  if (!theComparitor->compatible(tmpTriplet)) {
236  continue;
237  }
238  }
239 
240  float chi2 = std::numeric_limits<float>::quiet_NaN();
241  // TODO: Do we have any use case to not use bending correction?
242  if (useBendingCorrection) {
243  // Following PixelFitterByConformalMappingAndLine
244  const float simpleCot = (gps.back().z() - gps.front().z()) / (gps.back().perp() - gps.front().perp());
245  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
246  for (int i = 0; i < 4; ++i) {
247  const GlobalPoint& point = gps[i];
248  const GlobalError& error = ges[i];
249  bc_r[i] = sqrt(sqr(point.x() - region.origin().x()) + sqr(point.y() - region.origin().y()));
251  bc_z[i] = point.z() - region.origin().z();
252  bc_errZ2[i] = (barrels[i]) ? error.czz() : error.rerr(point) * sqr(simpleCot);
253  }
254  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
255  chi2 = rzLine.chi2();
256  } else {
257  RZLine rzLine(gps, ges, barrels);
258  chi2 = rzLine.chi2();
259  }
260  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2) {
261  continue;
262  }
263  // TODO: Do we have any use case to not use circle fit? Maybe
264  // HLT where low-pT inefficiency is not a problem?
265  if (fitFastCircle) {
266  FastCircleFit c(gps, ges);
267  chi2 += c.chi2();
268  if (edm::isNotFinite(chi2))
269  continue;
270  if (fitFastCircleChi2Cut && chi2 > thisMaxChi2)
271  continue;
272  }
273  result[index].emplace_back(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
274  allCells[foundQuadruplets[quadId][1]].getInnerHit(),
275  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
276  allCells[foundQuadruplets[quadId][2]].getOuterHit());
277  }
278  index++;
279  }
280 }
ihd::RegionLayerSets::region
const TrackingRegion & region() const
Definition: IntermediateHitDoublets.h:61
egmGedGsfElectronPFIsolation_cff.vertexIndex
vertexIndex
Definition: egmGedGsfElectronPFIsolation_cff.py:54
DDAxes::y
Handle.h
CAHitQuadrupletGenerator::maxChi2
const QuantityDependsPt maxChi2
Definition: CAHitQuadrupletGenerator.h:126
electrons_cff.bool
bool
Definition: electrons_cff.py:366
mps_fire.i
i
Definition: mps_fire.py:428
edm::ParameterSetDescription::add
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: ParameterSetDescription.h:95
PixelSubdetector::PixelBarrel
Definition: PixelSubdetector.h:11
hfClusterShapes_cfi.hits
hits
Definition: hfClusterShapes_cfi.py:5
CAHitQuadrupletGenerator::minLayers
static constexpr unsigned int minLayers
Definition: CAHitQuadrupletGenerator.h:34
IntermediateHitDoublets::LayerPairHitDoublets
Definition: IntermediateHitDoublets.h:145
SeedingHitSet
Definition: SeedingHitSet.h:6
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
CellularAutomaton::evolve
void evolve(const unsigned int)
Definition: CellularAutomaton.cc:91
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
HLT_FULL_cff.extraHitRPhitolerance
extraHitRPhitolerance
Definition: HLT_FULL_cff.py:9864
CACut::setCutValuesByLayerIds
void setCutValuesByLayerIds(CAGraph &caLayers)
Definition: CACut.h:43
edm
HLT enums.
Definition: AlignableModifier.h:19
edm::ParameterSetDescription
Definition: ParameterSetDescription.h:52
SeedingLayerSetsHits
Definition: SeedingLayerSetsHits.h:18
CAHitQuadrupletGenerator::QuantityDependsPtEval::value
float value(float curvature) const
Definition: CAHitQuadrupletGenerator.h:64
ThirdHitPredictionFromCircle.h
sqr
int sqr(const T &t)
Definition: pfalgo_common_ref.h:9
PixelRecoUtilities::inversePt
T inversePt(T curvature, const edm::EventSetup &iSetup)
Definition: PixelRecoUtilities.h:48
beam_dqm_sourceclient-live_cfg.maxChi2
maxChi2
Definition: beam_dqm_sourceclient-live_cfg.py:144
CellularAutomaton::createAndConnectCells
void createAndConnectCells(const std::vector< const HitDoublets * > &, const TrackingRegion &, const CACut &, const CACut &, const float)
Definition: CellularAutomaton.cc:5
DDAxes::x
findQualityFiles.v
v
Definition: findQualityFiles.py:179
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
HLT_FULL_cff.fitFastCircle
fitFastCircle
Definition: HLT_FULL_cff.py:9866
PixelRecoUtilities::curvature
T curvature(T InversePt, const edm::EventSetup &iSetup)
Definition: PixelRecoUtilities.h:42
spr::find
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
ihd::RegionLayerSets
Definition: IntermediateHitDoublets.h:48
newFWLiteAna.found
found
Definition: newFWLiteAna.py:118
relativeConstraints.error
error
Definition: relativeConstraints.py:53
FastCircleFit
Definition: FastCircleFit.h:27
RZLine::chi2
float chi2() const
Definition: RZLine.h:94
HLT_FULL_cff.fitFastCircleChi2Cut
fitFastCircleChi2Cut
Definition: HLT_FULL_cff.py:9887
CAHitQuadrupletGenerator::fitFastCircleChi2Cut
const bool fitFastCircleChi2Cut
Definition: CAHitQuadrupletGenerator.h:128
CellularAutomaton::getAllCells
std::vector< CACell > & getAllCells()
Definition: CellularAutomaton.h:16
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
CAHitQuadrupletGenerator::hitNtuplets
void hitNtuplets(const IntermediateHitDoublets &regionDoublets, std::vector< OrderedHitSeeds > &result, const edm::EventSetup &es, const SeedingLayerSetsHits &layers)
Definition: CAHitQuadrupletGenerator.cc:160
Point3DBase< float, GlobalTag >
ParameterSetDescription.h
CAHitQuadrupletGenerator::extraHitRPhitolerance
const float extraHitRPhitolerance
Definition: CAHitQuadrupletGenerator.h:124
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
CAHitQuadrupletGenerator::QuantityDependsPt::evaluator
QuantityDependsPtEval evaluator(const edm::EventSetup &es) const
Definition: CAHitQuadrupletGenerator.h:106
Basic2DVector< Scalar >
CAHitQuadrupletGenerator::caThetaCut
CACut caThetaCut
Definition: CAHitQuadrupletGenerator.h:131
PixelPluginsPhase0_cfi.isBarrel
isBarrel
Definition: PixelPluginsPhase0_cfi.py:17
ihd::RegionLayerSets::end
const_iterator end() const
Definition: IntermediateHitDoublets.h:66
edm::ParameterSet
Definition: ParameterSet.h:47
Event.h
ParameterSet
Definition: Functions.h:16
CAHitQuadrupletGenerator::QuantityDependsPtEval
Definition: CAHitQuadrupletGenerator.h:59
HLT_FULL_cff.useBendingCorrection
useBendingCorrection
Definition: HLT_FULL_cff.py:9886
edm::ParameterSet::addParameter
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
CAHitQuadrupletGenerator::useBendingCorrection
const bool useBendingCorrection
Definition: CAHitQuadrupletGenerator.h:129
HLT_FULL_cff.region
region
Definition: HLT_FULL_cff.py:88267
CAHitQuadrupletGenerator.h
CAHitQuadrupletGenerator::fitFastCircle
const bool fitFastCircle
Definition: CAHitQuadrupletGenerator.h:127
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
CAHitQuadrupletGenerator::fillDescriptions
static void fillDescriptions(edm::ParameterSetDescription &desc)
Definition: CAHitQuadrupletGenerator.cc:46
BarrelDetLayer.h
IntermediateHitDoublets
Definition: IntermediateHitDoublets.h:131
GlobalErrorBase< double, ErrorMatrixTag >
CellularAutomaton
Definition: CellularAutomaton.h:12
CAGraph
Definition: CAGraph.h:50
edm::EventSetup
Definition: EventSetup.h:58
get
#define get
CALayerPair
Definition: CAGraph.h:34
CAHitQuadrupletGenerator::initEvent
void initEvent(const edm::Event &ev, const edm::EventSetup &es)
Definition: CAHitQuadrupletGenerator.cc:84
DynArray.h
looper.cfg
cfg
Definition: looper.py:297
CAGraph.h
RZLine
Definition: RZLine.h:12
submitPVResolutionJobs.desc
string desc
Definition: submitPVResolutionJobs.py:251
CAHitQuadrupletGenerator::theComparitor
std::unique_ptr< SeedComparitor > theComparitor
Definition: CAHitQuadrupletGenerator.h:57
std
Definition: JetResolutionObject.h:76
isFinite.h
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
T
long double T
Definition: Basic3DVectorLD.h:48
ihd::RegionLayerSets::begin
const_iterator begin() const
Definition: IntermediateHitDoublets.h:64
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
genVertex_cff.x
x
Definition: genVertex_cff.py:12
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
RZLine::ErrZ2_tag
Definition: RZLine.h:14
CAHitQuadrupletGenerator::caHardPtCut
const float caHardPtCut
Definition: CAHitQuadrupletGenerator.h:133
pixelrecoutilities::LongitudinalBendingCorrection
Definition: LongitudinalBendingCorrection.h:6
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
TrackingRegion
Definition: TrackingRegion.h:41
mps_fire.result
result
Definition: mps_fire.py:311
ConsumesCollector.h
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ThirdHitPredictionFromCircle
Definition: ThirdHitPredictionFromCircle.h:10
c
auto & c
Definition: CAHitNtupletGeneratorKernelsImpl.h:46
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
CellularAutomaton::findNtuplets
void findNtuplets(std::vector< CACell::CAntuplet > &, const unsigned int)
Definition: CellularAutomaton.cc:127
point
*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
edm::Event
Definition: Event.h:73
CAHitQuadrupletGenerator::CAHitQuadrupletGenerator
CAHitQuadrupletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
Definition: CAHitQuadrupletGenerator.h:38
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
hgcalTopologyTester_cfi.layers
layers
Definition: hgcalTopologyTester_cfi.py:8
CellularAutomaton.h
CAHitQuadrupletGenerator::caPhiCut
CACut caPhiCut
Definition: CAHitQuadrupletGenerator.h:132
g
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
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443