CMS 3D CMS Logo

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