CMS 3D CMS Logo

CAHitTripletGenerator.cc
Go to the documentation of this file.
1 #include <unordered_map>
2 
13 
14 #include "CAGraph.h"
15 #include "CellularAutomaton.h"
16 
17 namespace
18 {
19  template <typename T>
20  T sqr(T x)
21  {
22  return x * x;
23  }
24 }
25 
26 using namespace std;
27 
29 
32  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
33  maxChi2(cfg.getParameter < edm::ParameterSet > ("maxChi2")),
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,
44  comparitorPSet, iC)};
45  }
46 }
47 
49  desc.add<double>("extraHitRPhitolerance", 0.06);
50  desc.add<bool>("useBendingCorrection", false);
51  desc.add<double>("CAThetaCut", 0.00125);
52  desc.add<double>("CAPhiCut", 0.1);
53  desc.add<double>("CAHardPtCut", 0);
54 
55  edm::ParameterSetDescription descMaxChi2;
56  descMaxChi2.add<double>("pt1", 0.8);
57  descMaxChi2.add<double>("pt2", 2);
58  descMaxChi2.add<double>("value1", 50);
59  descMaxChi2.add<double>("value2", 8);
60  descMaxChi2.add<bool>("enabled", true);
61  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
62 
63  edm::ParameterSetDescription descComparitor;
64  descComparitor.add<std::string>("ComponentName", "none");
65  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
66  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
67 }
68 
70  if (theComparitor) theComparitor->init(ev, es);
71 }
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 < 3; ++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(),
85  layers[i][j].hits().size());
86  vertexIndex = g.theLayers.size() - 1;
87  }
88  else
89  {
90  vertexIndex = foundVertex - g.theLayers.begin();
91  }
92  if (j == 0)
93  {
94 
95  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(),
96  vertexIndex) == g.theRootLayers.end())
97  {
98  g.theRootLayers.emplace_back(vertexIndex);
99 
100  }
101 
102  }
103  }
104  }
105  }
106 
107  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
108  g.theLayerPairs.clear();
109  for (unsigned int i = 0; i < g.theLayers.size(); i++ ){
110  g.theLayers[i].theInnerLayers.clear();
111  g.theLayers[i].theInnerLayerPairs.clear();
112  g.theLayers[i].theOuterLayers.clear();
113  g.theLayers[i].theOuterLayerPairs.clear();
114  for (auto & v : g.theLayers[i].isOuterHitOfCell) v.clear();
115  }
116 
117  }
118 
119  void fillGraph(const SeedingLayerSetsHits& layers, const IntermediateHitDoublets::RegionLayerSets& regionLayerPairs,
120  CAGraph& g, std::vector<const HitDoublets *>& hitDoublets) {
121  for (unsigned int i = 0; i < layers.size(); i++)
122  {
123  for (unsigned int j = 0; j < 3; ++j)
124  {
125  auto vertexIndex = 0;
126  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
127  layers[i][j].name());
128 
129  if (foundVertex == g.theLayers.end())
130  {
131  vertexIndex = g.theLayers.size() - 1;
132  }
133  else
134  {
135  vertexIndex = foundVertex - g.theLayers.begin();
136  }
137 
138 
139  if (j > 0)
140  {
141 
142  auto innerVertex = std::find(g.theLayers.begin(),
143  g.theLayers.end(), layers[i][j - 1].name());
144 
145  CALayerPair tmpInnerLayerPair(innerVertex - g.theLayers.begin(),
146  vertexIndex);
147 
148  if (std::find(g.theLayerPairs.begin(), g.theLayerPairs.end(),
149  tmpInnerLayerPair) == g.theLayerPairs.end())
150  {
151  auto found = std::find_if(regionLayerPairs.begin(), regionLayerPairs.end(), [&](const IntermediateHitDoublets::LayerPairHitDoublets& pair) {
152  return pair.innerLayerIndex() == layers[i][j - 1].index() && pair.outerLayerIndex() == layers[i][j].index();
153  });
154  if(found != regionLayerPairs.end()) {
155  hitDoublets.emplace_back(&(found->doublets()));
156  g.theLayerPairs.push_back(tmpInnerLayerPair);
157  g.theLayers[vertexIndex].theInnerLayers.push_back(
158  innerVertex - g.theLayers.begin());
159  innerVertex->theOuterLayers.push_back(vertexIndex);
160  g.theLayers[vertexIndex].theInnerLayerPairs.push_back(
161  g.theLayerPairs.size() - 1);
162  innerVertex->theOuterLayerPairs.push_back(
163  g.theLayerPairs.size() - 1);
164  }
165  }
166 
167  }
168 
169  }
170  }
171  }
172 }
173 
175  std::vector<OrderedHitSeeds>& result,
176  const edm::EventSetup& es,
177  const SeedingLayerSetsHits& layers) {
178  CAGraph g;
179 
180  std::vector<const HitDoublets *> hitDoublets;
181 
182  std::vector<CACell::CAntuplet> foundTriplets;
183 
184  int index =0;
185  for(const auto& regionLayerPairs: regionDoublets) {
186 
187  const TrackingRegion& region = regionLayerPairs.region();
188  hitDoublets.clear();
189  foundTriplets.clear();
190 
191  if (index == 0){
192  createGraphStructure(layers, g);
193  }
194  else{
195  clearGraphStructure(layers, g);
196  }
197  fillGraph(layers, regionLayerPairs, g, hitDoublets);
198  CellularAutomaton ca(g);
199  ca.findTriplets(hitDoublets, foundTriplets, region, caThetaCut, caPhiCut,
200  caHardPtCut);
201 
202  auto & allCells = ca.getAllCells();
203 
204  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
205 
206  // re-used thoughout, need to be vectors because of RZLine interface
207  std::array<float, 3> bc_r;
208  std::array<float, 3> bc_z;
209  std::array<float, 3> bc_errZ2;
210  std::array<GlobalPoint, 3> gps;
211  std::array<GlobalError, 3> ges;
212  std::array<bool, 3> barrels;
213 
214  unsigned int numberOfFoundTriplets = foundTriplets.size();
215  for (unsigned int tripletId = 0; tripletId < numberOfFoundTriplets;
216  ++tripletId)
217  {
218 
219  OrderedHitTriplet tmpTriplet(allCells[foundTriplets[tripletId][0]].getInnerHit(),
220  allCells[foundTriplets[tripletId][0]].getOuterHit(),
221  allCells[foundTriplets[tripletId][1]].getOuterHit());
222 
223  auto isBarrel = [](const unsigned id) -> bool
224  {
225  return id == PixelSubdetector::PixelBarrel;
226  };
227  for (unsigned int i = 0; i < 2; ++i)
228  {
229  auto const& ahit = allCells[foundTriplets[tripletId][i]].getInnerHit();
230  gps[i] = ahit->globalPosition();
231  ges[i] = ahit->globalPositionError();
232  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
233  }
234 
235  auto const& ahit = allCells[foundTriplets[tripletId][1]].getOuterHit();
236  gps[2] = ahit->globalPosition();
237  ges[2] = ahit->globalPositionError();
238  barrels[2] = isBarrel(ahit->geographicalId().subdetId());
239 
240  PixelRecoLineRZ line(gps[0], gps[2]);
241  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2],
243  const float curvature = predictionRPhi.curvature(
244  ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
245  const float abscurv = std::abs(curvature);
246  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
247  float chi2 = std::numeric_limits<float>::quiet_NaN();
248  // TODO: Do we have any use case to not use bending correction?
249 
251  {
252  // Following PixelFitterByConformalMappingAndLine
253  const float simpleCot = (gps.back().z() - gps.front().z())
254  / (gps.back().perp() - gps.front().perp());
255  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
256  for (int i = 0; i < 3; ++i)
257  {
258  const GlobalPoint & point = gps[i];
259  const GlobalError & error = ges[i];
260  bc_r[i] = sqrt(
261  sqr(point.x() - region.origin().x())
262  + sqr(point.y() - region.origin().y()));
264  es)(bc_r[i]);
265  bc_z[i] = point.z() - region.origin().z();
266  bc_errZ2[i] =
267  (barrels[i]) ?
268  error.czz() :
269  error.rerr(point) * sqr(simpleCot);
270  }
271  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
272  chi2 = rzLine.chi2();
273  }
274  else
275  {
276  RZLine rzLine(gps, ges, barrels);
277  chi2 = rzLine.chi2();
278  }
279 
280  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
281  {
282  continue;
283 
284  }
285 
286  if (theComparitor)
287  {
288  if (!theComparitor->compatible(tmpTriplet))
289  {
290 
291  continue;
292  }
293  }
294  result[index].emplace_back(tmpTriplet);
295 
296  }
297  index++;
298  }
299 }
size
Write out results.
T getParameter(std::string const &) const
const QuantityDependsPt maxChi2
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
const_iterator begin() const
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()
void hitNtuplets(const IntermediateHitDoublets &regionDoublets, std::vector< OrderedHitSeeds > &result, const edm::EventSetup &es, const SeedingLayerSetsHits &layers)
bool ev
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
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
CAHitTripletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
Range curvature(double transverseIP) const
T inversePt(T curvature, const edm::EventSetup &iSetup)
QuantityDependsPtEval evaluator(const edm::EventSetup &es) const
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void initEvent(const edm::Event &ev, const edm::EventSetup &es)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
std::vector< CALayerPair > theLayerPairs
Definition: CAGraph.h:54
Definition: RZLine.h:12
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
static void fillDescriptions(edm::ParameterSetDescription &desc)
void findTriplets(const std::vector< const HitDoublets * > &hitDoublets, std::vector< CACell::CAntuplet > &foundTriplets, const TrackingRegion &region, const float thetaCut, const float phiCut, const float hardPtCut)
long double T
T x() const
Definition: PV3DBase.h:62
std::unique_ptr< SeedComparitor > theComparitor
float chi2() const
Definition: RZLine.h:97
unsigned short size() const
Get the number of SeedingLayerSets.
static unsigned int minLayers
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