CMS 3D CMS Logo

CAHitTripletGenerator.cc
Go to the documentation of this file.
1 #include <unordered_map>
2 
4 
6 
12 
14 
15 #include "CellularAutomaton.h"
16 
18 
20 
21 namespace
22 {
23 
24 template<typename T>
25 T sqr(T x)
26 {
27  return x * x;
28 }
29 }
30 
31 using namespace std;
32 
34 
37  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
38  maxChi2(cfg.getParameter < edm::ParameterSet > ("maxChi2")),
39  useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
40  caThetaCut( cfg.getParameter<double>("CAThetaCut")),
41  caPhiCut(cfg.getParameter<double>("CAPhiCut")),
42  caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
43 {
44  edm::ParameterSet comparitorPSet = cfg.getParameter < edm::ParameterSet > ("SeedComparitorPSet");
45  std::string comparitorName = comparitorPSet.getParameter < std::string > ("ComponentName");
46  if (comparitorName != "none")
47  {
48  theComparitor.reset(
49  SeedComparitorFactory::get()->create(comparitorName,
50  comparitorPSet, iC));
51  }
52 }
53 
55  desc.add<double>("extraHitRPhitolerance", 0.06);
56  desc.add<bool>("useBendingCorrection", false);
57  desc.add<double>("CAThetaCut", 0.00125);
58  desc.add<double>("CAPhiCut", 0.1);
59  desc.add<double>("CAHardPtCut", 0);
60 
61  edm::ParameterSetDescription descMaxChi2;
62  descMaxChi2.add<double>("pt1", 0.8);
63  descMaxChi2.add<double>("pt2", 2);
64  descMaxChi2.add<double>("value1", 50);
65  descMaxChi2.add<double>("value2", 8);
66  descMaxChi2.add<bool>("enabled", true);
67  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
68 
69  edm::ParameterSetDescription descComparitor;
70  descComparitor.add<std::string>("ComponentName", "none");
71  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
72  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
73 }
74 
76  if (theComparitor) theComparitor->init(ev, es);
77 }
78 
79 namespace {
80  void createGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
81  for (unsigned int i = 0; i < layers.size(); i++)
82  {
83  for (unsigned int j = 0; j < 3; ++j)
84  {
85  auto vertexIndex = 0;
86  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
87  layers[i][j].name());
88  if (foundVertex == g.theLayers.end())
89  {
90  g.theLayers.emplace_back(layers[i][j].name(),
91  layers[i][j].hits().size());
92  vertexIndex = g.theLayers.size() - 1;
93  }
94  else
95  {
96  vertexIndex = foundVertex - g.theLayers.begin();
97  }
98  if (j == 0)
99  {
100 
101  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(),
102  vertexIndex) == g.theRootLayers.end())
103  {
104  g.theRootLayers.emplace_back(vertexIndex);
105 
106  }
107 
108  }
109  }
110  }
111  }
112 
113  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
114  g.theLayerPairs.clear();
115  for (unsigned int i = 0; i < g.theLayers.size(); i++ ){
116  g.theLayers[i].theInnerLayers.clear();
117  g.theLayers[i].theInnerLayerPairs.clear();
118  g.theLayers[i].theOuterLayers.clear();
119  g.theLayers[i].theOuterLayerPairs.clear();
120  for (auto & v : g.theLayers[i].isOuterHitOfCell) v.clear();
121  }
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 < 3; ++j)
130  {
131  auto vertexIndex = 0;
132  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
133  layers[i][j].name());
134 
135  if (foundVertex == g.theLayers.end())
136  {
137  vertexIndex = g.theLayers.size() - 1;
138  }
139  else
140  {
141  vertexIndex = foundVertex - g.theLayers.begin();
142  }
143 
144 
145  if (j > 0)
146  {
147 
148  auto innerVertex = std::find(g.theLayers.begin(),
149  g.theLayers.end(), layers[i][j - 1].name());
150 
151  CALayerPair tmpInnerLayerPair(innerVertex - g.theLayers.begin(),
152  vertexIndex);
153 
154  if (std::find(g.theLayerPairs.begin(), g.theLayerPairs.end(),
155  tmpInnerLayerPair) == g.theLayerPairs.end())
156  {
157  auto found = std::find_if(regionLayerPairs.begin(), regionLayerPairs.end(), [&](const IntermediateHitDoublets::LayerPairHitDoublets& pair) {
158  return pair.innerLayerIndex() == layers[i][j - 1].index() && pair.outerLayerIndex() == layers[i][j].index();
159  });
160  if(found != regionLayerPairs.end()) {
161  hitDoublets.emplace_back(&(found->doublets()));
162  g.theLayerPairs.push_back(tmpInnerLayerPair);
163  g.theLayers[vertexIndex].theInnerLayers.push_back(
164  innerVertex - g.theLayers.begin());
165  innerVertex->theOuterLayers.push_back(vertexIndex);
166  g.theLayers[vertexIndex].theInnerLayerPairs.push_back(
167  g.theLayerPairs.size() - 1);
168  innerVertex->theOuterLayerPairs.push_back(
169  g.theLayerPairs.size() - 1);
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  std::vector<CACell::CAntuplet> foundTriplets;
189 
190  int index =0;
191  for(const auto& regionLayerPairs: regionDoublets) {
192 
193  const TrackingRegion& region = regionLayerPairs.region();
194  hitDoublets.clear();
195  foundTriplets.clear();
196 
197  if (index == 0){
198  createGraphStructure(layers, g);
199  }
200  else{
201  clearGraphStructure(layers, g);
202  }
203  fillGraph(layers, regionLayerPairs, g, hitDoublets);
204  CellularAutomaton ca(g);
205  ca.findTriplets(hitDoublets, foundTriplets, region, caThetaCut, caPhiCut,
206  caHardPtCut);
207 
208  auto & allCells = ca.getAllCells();
209 
210  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
211 
212  // re-used thoughout, need to be vectors because of RZLine interface
213  std::array<float, 3> bc_r;
214  std::array<float, 3> bc_z;
215  std::array<float, 3> bc_errZ2;
216  std::array<GlobalPoint, 3> gps;
217  std::array<GlobalError, 3> ges;
218  std::array<bool, 3> barrels;
219 
220  unsigned int numberOfFoundTriplets = foundTriplets.size();
221  for (unsigned int tripletId = 0; tripletId < numberOfFoundTriplets;
222  ++tripletId)
223  {
224 
225  OrderedHitTriplet tmpTriplet(allCells[foundTriplets[tripletId][0]].getInnerHit(),
226  allCells[foundTriplets[tripletId][0]].getOuterHit(),
227  allCells[foundTriplets[tripletId][1]].getOuterHit());
228 
229  auto isBarrel = [](const unsigned id) -> bool
230  {
231  return id == PixelSubdetector::PixelBarrel;
232  };
233  for (unsigned int i = 0; i < 2; ++i)
234  {
235  auto const& ahit = allCells[foundTriplets[tripletId][i]].getInnerHit();
236  gps[i] = ahit->globalPosition();
237  ges[i] = ahit->globalPositionError();
238  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
239  }
240 
241  auto const& ahit = allCells[foundTriplets[tripletId][1]].getOuterHit();
242  gps[2] = ahit->globalPosition();
243  ges[2] = ahit->globalPositionError();
244  barrels[2] = isBarrel(ahit->geographicalId().subdetId());
245 
246  PixelRecoLineRZ line(gps[0], gps[2]);
247  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2],
249  const float curvature = predictionRPhi.curvature(
250  ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
251  const float abscurv = std::abs(curvature);
252  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
253  float chi2 = std::numeric_limits<float>::quiet_NaN();
254  // TODO: Do we have any use case to not use bending correction?
255 
257  {
258  // Following PixelFitterByConformalMappingAndLine
259  const float simpleCot = (gps.back().z() - gps.front().z())
260  / (gps.back().perp() - gps.front().perp());
261  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
262  for (int i = 0; i < 3; ++i)
263  {
264  const GlobalPoint & point = gps[i];
265  const GlobalError & error = ges[i];
266  bc_r[i] = sqrt(
267  sqr(point.x() - region.origin().x())
268  + sqr(point.y() - region.origin().y()));
270  es)(bc_r[i]);
271  bc_z[i] = point.z() - region.origin().z();
272  bc_errZ2[i] =
273  (barrels[i]) ?
274  error.czz() :
275  error.rerr(point) * sqr(simpleCot);
276  }
277  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
278  chi2 = rzLine.chi2();
279  }
280  else
281  {
282  RZLine rzLine(gps, ges, barrels);
283  chi2 = rzLine.chi2();
284  }
285 
286  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
287  {
288  continue;
289 
290  }
291 
292  if (theComparitor)
293  {
294  if (!theComparitor->compatible(tmpTriplet))
295  {
296 
297  continue;
298  }
299  }
300  result[index].emplace_back(tmpTriplet);
301 
302  }
303  index++;
304  }
305 }
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
def create(alignables, pedeDump, additionalData, outputFile, config)
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
#define constexpr
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
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
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:64
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:65
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
*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