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