CMS 3D CMS Logo

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