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 
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  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 }
43 
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 }
69 
71  if (theComparitor)
72  theComparitor->init(ev, es);
73 }
74 namespace {
75  void createGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
76  for (unsigned int i = 0; i < layers.size(); i++) {
77  for (unsigned int j = 0; j < 4; ++j) {
78  auto vertexIndex = 0;
79  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j].name());
80  if (foundVertex == g.theLayers.end()) {
81  g.theLayers.emplace_back(layers[i][j].name(), layers[i][j].hits().size());
82  vertexIndex = g.theLayers.size() - 1;
83  } else {
84  vertexIndex = foundVertex - g.theLayers.begin();
85  }
86  if (j == 0) {
87  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(), vertexIndex) == g.theRootLayers.end()) {
88  g.theRootLayers.emplace_back(vertexIndex);
89  }
90  }
91  }
92  }
93  }
94  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
95  g.theLayerPairs.clear();
96  for (unsigned int i = 0; i < g.theLayers.size(); i++) {
97  g.theLayers[i].theInnerLayers.clear();
98  g.theLayers[i].theInnerLayerPairs.clear();
99  g.theLayers[i].theOuterLayers.clear();
100  g.theLayers[i].theOuterLayerPairs.clear();
101  for (auto& v : g.theLayers[i].isOuterHitOfCell)
102  v.clear();
103  }
104  }
105  void fillGraph(const SeedingLayerSetsHits& layers,
106  const IntermediateHitDoublets::RegionLayerSets& regionLayerPairs,
107  CAGraph& g,
108  std::vector<const HitDoublets*>& hitDoublets) {
109  for (unsigned int i = 0; i < layers.size(); i++) {
110  for (unsigned int j = 0; j < 4; ++j) {
111  auto vertexIndex = 0;
112  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j].name());
113  if (foundVertex == g.theLayers.end()) {
114  vertexIndex = g.theLayers.size() - 1;
115  } else {
116  vertexIndex = foundVertex - g.theLayers.begin();
117  }
118 
119  if (j > 0) {
120  auto innerVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j - 1].name());
121 
122  CALayerPair tmpInnerLayerPair(innerVertex - g.theLayers.begin(), vertexIndex);
123 
124  if (std::find(g.theLayerPairs.begin(), g.theLayerPairs.end(), tmpInnerLayerPair) == g.theLayerPairs.end()) {
125  auto found = std::find_if(regionLayerPairs.begin(),
126  regionLayerPairs.end(),
128  return pair.innerLayerIndex() == layers[i][j - 1].index() &&
129  pair.outerLayerIndex() == layers[i][j].index();
130  });
131  if (found != regionLayerPairs.end()) {
132  hitDoublets.emplace_back(&(found->doublets()));
133  g.theLayerPairs.push_back(tmpInnerLayerPair);
134  g.theLayers[vertexIndex].theInnerLayers.push_back(innerVertex - g.theLayers.begin());
135  innerVertex->theOuterLayers.push_back(vertexIndex);
136  g.theLayers[vertexIndex].theInnerLayerPairs.push_back(g.theLayerPairs.size() - 1);
137  innerVertex->theOuterLayerPairs.push_back(g.theLayerPairs.size() - 1);
138  }
139  }
140  }
141  }
142  }
143  }
144 } // namespace
145 
147  std::vector<OrderedHitSeeds>& result,
148  const edm::EventSetup& es,
149  const SeedingLayerSetsHits& layers) {
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 }
size
Write out results.
T getParameter(std::string const &) const
void findNtuplets(std::vector< CACell::CAntuplet > &, const unsigned int)
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
const_iterator begin() const
void createAndConnectCells(const std::vector< const HitDoublets * > &, const TrackingRegion &, const float, const float, const float)
GlobalPoint const & origin() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
const_iterator end() const
T y() const
Definition: PV3DBase.h:60
std::vector< CACell > & getAllCells()
bool ev
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
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
std::vector< CALayer > theLayers
Definition: CAGraph.h:45
T z() const
Definition: PV3DBase.h:61
void initEvent(const edm::Event &ev, const edm::EventSetup &es)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< CALayerPair > theLayerPairs
Definition: CAGraph.h:46
Definition: RZLine.h:12
CAHitQuadrupletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
const QuantityDependsPt maxChi2
float chi2() const
Definition: FastCircleFit.h:66
T rerr(const GlobalPoint &aPoint) const
HLT enums.
Square< F >::type sqr(const F &f)
Definition: Square.h:14
std::vector< int > theRootLayers
Definition: CAGraph.h:47
void hitNtuplets(const IntermediateHitDoublets &regionDoublets, std::vector< OrderedHitSeeds > &result, const edm::EventSetup &es, const SeedingLayerSetsHits &layers)
static void fillDescriptions(edm::ParameterSetDescription &desc)
QuantityDependsPtEval evaluator(const edm::EventSetup &es) const
long double T
T x() const
Definition: PV3DBase.h:59
float chi2() const
Definition: RZLine.h:94
unsigned short size() const
Get the number of SeedingLayerSets.
#define constexpr
*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 evolve(const unsigned int)