CMS 3D CMS Logo

CAHitQuadrupletGenerator.cc
Go to the documentation of this file.
1 #include <functional>
2 
12 
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  : 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 }
46 
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 }
84 
86  if (theComparitor)
87  theComparitor->init(ev, es);
89 }
90 namespace {
91  void createGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
92  for (unsigned int i = 0; i < layers.size(); i++) {
93  for (unsigned int j = 0; j < 4; ++j) {
94  auto vertexIndex = 0;
95  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j].name());
96  if (foundVertex == g.theLayers.end()) {
97  g.theLayers.emplace_back(layers[i][j].name(), layers[i][j].detLayer()->seqNum(), layers[i][j].hits().size());
98  vertexIndex = g.theLayers.size() - 1;
99  } else {
100  vertexIndex = foundVertex - g.theLayers.begin();
101  }
102  if (j == 0) {
103  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(), vertexIndex) == g.theRootLayers.end()) {
104  g.theRootLayers.emplace_back(vertexIndex);
105  }
106  }
107  }
108  }
109  }
110  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
111  g.theLayerPairs.clear();
112  for (unsigned int i = 0; i < g.theLayers.size(); i++) {
113  g.theLayers[i].theInnerLayers.clear();
114  g.theLayers[i].theInnerLayerPairs.clear();
115  g.theLayers[i].theOuterLayers.clear();
116  g.theLayers[i].theOuterLayerPairs.clear();
117  for (auto& v : g.theLayers[i].isOuterHitOfCell)
118  v.clear();
119  }
120  }
121  void fillGraph(const SeedingLayerSetsHits& layers,
122  const IntermediateHitDoublets::RegionLayerSets& regionLayerPairs,
123  CAGraph& g,
124  std::vector<const HitDoublets*>& hitDoublets) {
125  for (unsigned int i = 0; i < layers.size(); i++) {
126  for (unsigned int j = 0; j < 4; ++j) {
127  auto vertexIndex = 0;
128  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j].name());
129  if (foundVertex == g.theLayers.end()) {
130  vertexIndex = g.theLayers.size() - 1;
131  } else {
132  vertexIndex = foundVertex - g.theLayers.begin();
133  }
134 
135  if (j > 0) {
136  auto innerVertex = std::find(g.theLayers.begin(), g.theLayers.end(), layers[i][j - 1].name());
137 
138  CALayerPair tmpInnerLayerPair(innerVertex - g.theLayers.begin(), vertexIndex);
139 
140  if (std::find(g.theLayerPairs.begin(), g.theLayerPairs.end(), tmpInnerLayerPair) == g.theLayerPairs.end()) {
141  auto found = std::find_if(regionLayerPairs.begin(),
142  regionLayerPairs.end(),
144  return pair.innerLayerIndex() == layers[i][j - 1].index() &&
145  pair.outerLayerIndex() == layers[i][j].index();
146  });
147  if (found != regionLayerPairs.end()) {
148  hitDoublets.emplace_back(&(found->doublets()));
149  g.theLayerPairs.push_back(tmpInnerLayerPair);
150  g.theLayers[vertexIndex].theInnerLayers.push_back(innerVertex - g.theLayers.begin());
151  innerVertex->theOuterLayers.push_back(vertexIndex);
152  g.theLayers[vertexIndex].theInnerLayerPairs.push_back(g.theLayerPairs.size() - 1);
153  innerVertex->theOuterLayerPairs.push_back(g.theLayerPairs.size() - 1);
154  }
155  }
156  }
157  }
158  }
159  }
160 } // namespace
161 
163  std::vector<OrderedHitSeeds>& result,
164  const SeedingLayerSetsHits& layers) {
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 
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 }
size
Write out results.
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void findNtuplets(std::vector< CACell::CAntuplet > &, const unsigned int)
void setCutValuesByLayerIds(CAGraph &caLayers)
Definition: CACut.h:43
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
std::vector< CACell > & getAllCells()
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
static constexpr unsigned int minLayers
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_iterator begin() const
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > theFieldToken
T sqrt(T t)
Definition: SSEVec.h:19
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:135
void initEvent(const edm::Event &ev, const edm::EventSetup &es)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const MagneticField * theField
const TrackingRegion & region() const
bool getData(T &iHolder) const
Definition: EventSetup.h:122
T inversePt(T curvature, const MagneticField &field)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
Definition: RZLine.h:12
CAHitQuadrupletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
const QuantityDependsPt maxChi2
QuantityDependsPtEval evaluator(const MagneticField &field) const
float chi2() const
Definition: RZLine.h:94
void createAndConnectCells(const std::vector< const HitDoublets *> &, const TrackingRegion &, const CACut &, const CACut &, const float)
HLT enums.
Square< F >::type sqr(const F &f)
Definition: Square.h:14
void hitNtuplets(const IntermediateHitDoublets &regionDoublets, std::vector< OrderedHitSeeds > &result, const SeedingLayerSetsHits &layers)
float x
static void fillDescriptions(edm::ParameterSetDescription &desc)
#define get
long double T
*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)
const_iterator end() const