CMS 3D CMS Logo

CAHitQuadrupletGenerator.cc
Go to the documentation of this file.
2 
4 
9 
11 
12 
13 
14 #include "CellularAutomaton.h"
15 
17 
19 
20 #include <functional>
21 
22 namespace
23 {
24 
25  template <typename T>
26  T sqr(T x)
27  {
28  return x*x;
29  }
30 }
31 
32 
33 
34 using namespace std;
35 
37 
39 extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
40 maxChi2(cfg.getParameter<edm::ParameterSet>("maxChi2")),
41 fitFastCircle(cfg.getParameter<bool>("fitFastCircle")),
42 fitFastCircleChi2Cut(cfg.getParameter<bool>("fitFastCircleChi2Cut")),
43 useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
44 caThetaCut(cfg.getParameter<double>("CAThetaCut")),
45 caPhiCut(cfg.getParameter<double>("CAPhiCut")),
46 caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
47 {
48  edm::ParameterSet comparitorPSet = cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
49  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
50  if (comparitorName != "none")
51  {
52  theComparitor.reset(SeedComparitorFactory::get()->create(comparitorName, comparitorPSet, iC));
53  }
54 }
55 
57  desc.add<double>("extraHitRPhitolerance", 0.1);
58  desc.add<bool>("fitFastCircle", false);
59  desc.add<bool>("fitFastCircleChi2Cut", false);
60  desc.add<bool>("useBendingCorrection", false);
61  desc.add<double>("CAThetaCut", 0.00125);
62  desc.add<double>("CAPhiCut", 10);
63  desc.add<double>("CAHardPtCut", 0);
64  desc.addOptional<bool>("CAOnlyOneLastHitPerLayerFilter")->setComment("Deprecated and has no effect. To be fully removed later when the parameter is no longer used in HLT configurations.");
65  edm::ParameterSetDescription descMaxChi2;
66  descMaxChi2.add<double>("pt1", 0.2);
67  descMaxChi2.add<double>("pt2", 1.5);
68  descMaxChi2.add<double>("value1", 500);
69  descMaxChi2.add<double>("value2", 50);
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 }
78 
80  if (theComparitor) theComparitor->init(ev, es);
81 }
82 namespace {
83  void createGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
84  for (unsigned int i = 0; i < layers.size(); i++)
85  {
86  for (unsigned int j = 0; j < 4; ++j)
87  {
88  auto vertexIndex = 0;
89  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
90  layers[i][j].name());
91  if (foundVertex == g.theLayers.end())
92  {
93  g.theLayers.emplace_back(layers[i][j].name(), layers[i][j].hits().size());
94  vertexIndex = g.theLayers.size() - 1;
95  }
96  else
97  {
98  vertexIndex = foundVertex - g.theLayers.begin();
99  }
100  if (j == 0)
101  {
102 
103  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(),
104  vertexIndex) == g.theRootLayers.end())
105  {
106  g.theRootLayers.emplace_back(vertexIndex);
107 
108  }
109 
110  }
111  }
112  }
113  }
114  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
115  g.theLayerPairs.clear();
116  for (unsigned int i = 0; i < g.theLayers.size(); i++ ){
117  g.theLayers[i].theInnerLayers.clear();
118  g.theLayers[i].theInnerLayerPairs.clear();
119  g.theLayers[i].theOuterLayers.clear();
120  g.theLayers[i].theOuterLayerPairs.clear();
121  for (auto & v : g.theLayers[i].isOuterHitOfCell) v.clear();
122  }
123 
124  }
125  void fillGraph(const SeedingLayerSetsHits& layers, const IntermediateHitDoublets::RegionLayerSets& regionLayerPairs,
126  CAGraph& g, std::vector<const HitDoublets *>& hitDoublets) {
127  for (unsigned int i = 0; i < layers.size(); i++)
128  {
129  for (unsigned int j = 0; j < 4; ++j)
130  {
131  auto vertexIndex = 0;
132  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
133  layers[i][j].name());
134  if (foundVertex == g.theLayers.end())
135  {
136  vertexIndex = g.theLayers.size() - 1;
137  }
138  else
139  {
140  vertexIndex = foundVertex - g.theLayers.begin();
141  }
142 
143 
144  if (j > 0)
145  {
146 
147  auto innerVertex = std::find(g.theLayers.begin(),
148  g.theLayers.end(), layers[i][j - 1].name());
149 
150  CALayerPair tmpInnerLayerPair(innerVertex - g.theLayers.begin(),
151  vertexIndex);
152 
153  if (std::find(g.theLayerPairs.begin(), g.theLayerPairs.end(),
154  tmpInnerLayerPair) == g.theLayerPairs.end())
155  {
156  auto found = std::find_if(regionLayerPairs.begin(), regionLayerPairs.end(), [&](const IntermediateHitDoublets::LayerPairHitDoublets& pair) {
157  return pair.innerLayerIndex() == layers[i][j - 1].index() && pair.outerLayerIndex() == layers[i][j].index();
158  });
159  if(found != regionLayerPairs.end()) {
160  hitDoublets.emplace_back(&(found->doublets()));
161  g.theLayerPairs.push_back(tmpInnerLayerPair);
162  g.theLayers[vertexIndex].theInnerLayers.push_back(
163  innerVertex - g.theLayers.begin());
164  innerVertex->theOuterLayers.push_back(vertexIndex);
165  g.theLayers[vertexIndex].theInnerLayerPairs.push_back(
166  g.theLayerPairs.size() - 1);
167  innerVertex->theOuterLayerPairs.push_back(
168  g.theLayerPairs.size() - 1);
169  }
170  }
171 
172  }
173 
174  }
175  }
176  }
177 }
178 
179 
181  std::vector<OrderedHitSeeds>& result,
182  const edm::EventSetup& es,
183  const SeedingLayerSetsHits& layers) {
184  CAGraph g;
185 
186  std::vector<const HitDoublets *> hitDoublets;
187 
188  const int numberOfHitsInNtuplet = 4;
189  std::vector<CACell::CAntuplet> foundQuadruplets;
190 
191 
192 
193  int index =0;
194  for(const auto& regionLayerPairs: regionDoublets) {
195 
196  const TrackingRegion& region = regionLayerPairs.region();
197  hitDoublets.clear();
198  foundQuadruplets.clear();
199  if (index == 0){
200  createGraphStructure(layers, g);
201  }
202  else{
203  clearGraphStructure(layers, g);
204  }
205 
206  fillGraph(layers, regionLayerPairs, g, hitDoublets);
207 
208  CellularAutomaton ca(g);
209 
210  ca.createAndConnectCells(hitDoublets, region, caThetaCut,
212 
213  ca.evolve(numberOfHitsInNtuplet);
214 
215  ca.findNtuplets(foundQuadruplets, numberOfHitsInNtuplet);
216 
217  auto & allCells = ca.getAllCells();
218 
219  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
220 
221  // re-used thoughout
222  std::array<float, 4> bc_r;
223  std::array<float, 4> bc_z;
224  std::array<float, 4> bc_errZ2;
225  std::array<GlobalPoint, 4> gps;
226  std::array<GlobalError, 4> ges;
227  std::array<bool, 4> barrels;
228 
229  unsigned int numberOfFoundQuadruplets = foundQuadruplets.size();
230 
231  // Loop over quadruplets
232  for (unsigned int quadId = 0; quadId < numberOfFoundQuadruplets; ++quadId)
233  {
234 
235  auto isBarrel = [](const unsigned id) -> bool
236  {
237  return id == PixelSubdetector::PixelBarrel;
238  };
239  for(unsigned int i = 0; i< 3; ++i)
240  {
241  auto const& ahit = allCells[foundQuadruplets[quadId][i]].getInnerHit();
242  gps[i] = ahit->globalPosition();
243  ges[i] = ahit->globalPositionError();
244  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
245  }
246 
247  auto const& ahit = allCells[foundQuadruplets[quadId][2]].getOuterHit();
248  gps[3] = ahit->globalPosition();
249  ges[3] = ahit->globalPositionError();
250  barrels[3] = isBarrel(ahit->geographicalId().subdetId());
251  // TODO:
252  // - if we decide to always do the circle fit for 4 hits, we don't
253  // need ThirdHitPredictionFromCircle for the curvature; then we
254  // could remove extraHitRPhitolerance configuration parameter
255  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2], extraHitRPhitolerance);
256  const float curvature = predictionRPhi.curvature(ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
257  const float abscurv = std::abs(curvature);
258  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
259  if (theComparitor)
260  {
261  SeedingHitSet tmpTriplet(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
262  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
263  allCells[foundQuadruplets[quadId][2]].getOuterHit());
264 
265 
266  if (!theComparitor->compatible(tmpTriplet) )
267  {
268  continue;
269  }
270  }
271 
272  float chi2 = std::numeric_limits<float>::quiet_NaN();
273  // TODO: Do we have any use case to not use bending correction?
275  {
276  // Following PixelFitterByConformalMappingAndLine
277  const float simpleCot = ( gps.back().z() - gps.front().z() ) / (gps.back().perp() - gps.front().perp() );
278  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
279  for (int i=0; i < 4; ++i)
280  {
281  const GlobalPoint & point = gps[i];
282  const GlobalError & error = ges[i];
283  bc_r[i] = sqrt( sqr(point.x() - region.origin().x()) + sqr(point.y() - region.origin().y()) );
284  bc_r[i] += pixelrecoutilities::LongitudinalBendingCorrection(pt, es)(bc_r[i]);
285  bc_z[i] = point.z() - region.origin().z();
286  bc_errZ2[i] = (barrels[i]) ? error.czz() : error.rerr(point)*sqr(simpleCot);
287  }
288  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
289  chi2 = rzLine.chi2();
290  }
291  else
292  {
293  RZLine rzLine(gps, ges, barrels);
294  chi2 = rzLine.chi2();
295  }
296  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
297  {
298  continue;
299  }
300  // TODO: Do we have any use case to not use circle fit? Maybe
301  // HLT where low-pT inefficiency is not a problem?
302  if (fitFastCircle)
303  {
304  FastCircleFit c(gps, ges);
305  chi2 += c.chi2();
306  if (edm::isNotFinite(chi2)) continue;
307  if (fitFastCircleChi2Cut && chi2 > thisMaxChi2) continue;
308  }
309  result[index].emplace_back(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
310  allCells[foundQuadruplets[quadId][1]].getInnerHit(),
311  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
312  allCells[foundQuadruplets[quadId][2]].getOuterHit());
313  }
314  index++;
315  }
316 
317 }
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)
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
const_iterator begin() const
void createAndConnectCells(const std::vector< const HitDoublets * > &, const TrackingRegion &, const float, const float, const float)
GlobalPoint const & origin() const
def create(alignables, pedeDump, additionalData, outputFile, config)
const_iterator end() const
T y() const
Definition: PV3DBase.h:63
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:20
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
#define constexpr
T inversePt(T curvature, const edm::EventSetup &iSetup)
bool isNotFinite(T x)
Definition: isFinite.h:10
T curvature(T InversePt, const edm::EventSetup &iSetup)
T sqrt(T t)
Definition: SSEVec.h:18
std::vector< CALayer > theLayers
Definition: CAGraph.h:63
T z() const
Definition: PV3DBase.h:64
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:64
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:13
std::vector< int > theRootLayers
Definition: CAGraph.h:65
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:62
float chi2() const
Definition: RZLine.h:97
unsigned short size() const
Get the number of SeedingLayerSets.
T get(const Candidate &c)
Definition: component.h:55
*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)