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 using namespace ctfseeding;
33 
35 
38  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
39  maxChi2(cfg.getParameter < edm::ParameterSet > ("maxChi2")),
40  useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
41  caThetaCut( cfg.getParameter<double>("CAThetaCut")),
42  caPhiCut(cfg.getParameter<double>("CAPhiCut")),
43  caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
44 {
45  edm::ParameterSet comparitorPSet = cfg.getParameter < edm::ParameterSet > ("SeedComparitorPSet");
46  std::string comparitorName = comparitorPSet.getParameter < std::string > ("ComponentName");
47  if (comparitorName != "none")
48  {
49  theComparitor.reset(
50  SeedComparitorFactory::get()->create(comparitorName,
51  comparitorPSet, iC));
52  }
53 }
54 
56  desc.add<double>("extraHitRPhitolerance", 0.06);
57  desc.add<bool>("useBendingCorrection", false);
58  desc.add<double>("CAThetaCut", 0.00125);
59  desc.add<double>("CAPhiCut", 0.1);
60  desc.add<double>("CAHardPtCut", 0);
61 
62  edm::ParameterSetDescription descMaxChi2;
63  descMaxChi2.add<double>("pt1", 0.8);
64  descMaxChi2.add<double>("pt2", 2);
65  descMaxChi2.add<double>("value1", 50);
66  descMaxChi2.add<double>("value2", 8);
67  descMaxChi2.add<bool>("enabled", true);
68  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
69 
70  edm::ParameterSetDescription descComparitor;
71  descComparitor.add<std::string>("ComponentName", "none");
72  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
73  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
74 }
75 
77  if (theComparitor) theComparitor->init(ev, es);
78 }
79 
80 namespace {
81  void createGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
82  for (unsigned int i = 0; i < layers.size(); i++)
83  {
84  for (unsigned int j = 0; j < 3; ++j)
85  {
86  auto vertexIndex = 0;
87  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
88  layers[i][j].name());
89  if (foundVertex == g.theLayers.end())
90  {
91  g.theLayers.emplace_back(layers[i][j].name(),
92  layers[i][j].hits().size());
93  vertexIndex = g.theLayers.size() - 1;
94  }
95  else
96  {
97  vertexIndex = foundVertex - g.theLayers.begin();
98  }
99  if (j == 0)
100  {
101 
102  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(),
103  vertexIndex) == g.theRootLayers.end())
104  {
105  g.theRootLayers.emplace_back(vertexIndex);
106 
107  }
108 
109  }
110  }
111  }
112  }
113 
114  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
115  g.theLayerPairs.clear();
116  for (unsigned int i = 0; i < g.theLayers.size(); i++ ){
117  g.theLayers[i].theInnerLayers.clear();
118  g.theLayers[i].theInnerLayerPairs.clear();
119  g.theLayers[i].theOuterLayers.clear();
120  g.theLayers[i].theOuterLayerPairs.clear();
121  for (auto & v : g.theLayers[i].isOuterHitOfCell) v.clear();
122  }
123 
124  }
125 
126  void fillGraph(const SeedingLayerSetsHits& layers, const IntermediateHitDoublets::RegionLayerSets& regionLayerPairs,
127  CAGraph& g, std::vector<const HitDoublets *>& hitDoublets) {
128  for (unsigned int i = 0; i < layers.size(); i++)
129  {
130  for (unsigned int j = 0; j < 3; ++j)
131  {
132  auto vertexIndex = 0;
133  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
134  layers[i][j].name());
135 
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 
182  std::vector<OrderedHitSeeds>& result,
183  const edm::EventSetup& es,
184  const SeedingLayerSetsHits& layers) {
185  CAGraph g;
186 
187  std::vector<const HitDoublets *> hitDoublets;
188 
189  std::vector<CACell::CAntuplet> foundTriplets;
190 
191  int index =0;
192  for(const auto& regionLayerPairs: regionDoublets) {
193 
194  const TrackingRegion& region = regionLayerPairs.region();
195  hitDoublets.clear();
196  foundTriplets.clear();
197 
198  if (index == 0){
199  createGraphStructure(layers, g);
200  }
201  else{
202  clearGraphStructure(layers, g);
203  }
204  fillGraph(layers, regionLayerPairs, g, hitDoublets);
205  CellularAutomaton ca(g);
206  ca.findTriplets(hitDoublets, foundTriplets, region, caThetaCut, caPhiCut,
207  caHardPtCut);
208 
209  auto & allCells = ca.getAllCells();
210 
211  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
212 
213  // re-used thoughout, need to be vectors because of RZLine interface
214  std::array<float, 3> bc_r;
215  std::array<float, 3> bc_z;
216  std::array<float, 3> bc_errZ2;
217  std::array<GlobalPoint, 3> gps;
218  std::array<GlobalError, 3> ges;
219  std::array<bool, 3> barrels;
220 
221  unsigned int numberOfFoundTriplets = foundTriplets.size();
222  for (unsigned int tripletId = 0; tripletId < numberOfFoundTriplets;
223  ++tripletId)
224  {
225 
226  OrderedHitTriplet tmpTriplet(allCells[foundTriplets[tripletId][0]].getInnerHit(),
227  allCells[foundTriplets[tripletId][0]].getOuterHit(),
228  allCells[foundTriplets[tripletId][1]].getOuterHit());
229 
230  auto isBarrel = [](const unsigned id) -> bool
231  {
232  return id == PixelSubdetector::PixelBarrel;
233  };
234  for (unsigned int i = 0; i < 2; ++i)
235  {
236  auto const& ahit = allCells[foundTriplets[tripletId][i]].getInnerHit();
237  gps[i] = ahit->globalPosition();
238  ges[i] = ahit->globalPositionError();
239  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
240  }
241 
242  auto const& ahit = allCells[foundTriplets[tripletId][1]].getOuterHit();
243  gps[2] = ahit->globalPosition();
244  ges[2] = ahit->globalPositionError();
245  barrels[2] = isBarrel(ahit->geographicalId().subdetId());
246 
247  PixelRecoLineRZ line(gps[0], gps[2]);
248  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2],
250  const float curvature = predictionRPhi.curvature(
251  ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
252  const float abscurv = std::abs(curvature);
253  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
254  float chi2 = std::numeric_limits<float>::quiet_NaN();
255  // TODO: Do we have any use case to not use bending correction?
256 
258  {
259  // Following PixelFitterByConformalMappingAndLine
260  const float simpleCot = (gps.back().z() - gps.front().z())
261  / (gps.back().perp() - gps.front().perp());
262  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
263  for (int i = 0; i < 3; ++i)
264  {
265  const GlobalPoint & point = gps[i];
266  const GlobalError & error = ges[i];
267  bc_r[i] = sqrt(
268  sqr(point.x() - region.origin().x())
269  + sqr(point.y() - region.origin().y()));
271  es)(bc_r[i]);
272  bc_z[i] = point.z() - region.origin().z();
273  bc_errZ2[i] =
274  (barrels[i]) ?
275  error.czz() :
276  error.rerr(point) * sqr(simpleCot);
277  }
278  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
279  chi2 = rzLine.chi2();
280  }
281  else
282  {
283  RZLine rzLine(gps, ges, barrels);
284  chi2 = rzLine.chi2();
285  }
286 
287  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
288  {
289  continue;
290 
291  }
292 
293  if (theComparitor)
294  {
295  if (!theComparitor->compatible(tmpTriplet))
296  {
297 
298  continue;
299  }
300  }
301  result[index].emplace_back(tmpTriplet);
302 
303  }
304  index++;
305  }
306 }
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)
T x() const
Cartesian x coordinate.
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