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.reset(
44  SeedComparitorFactory::get()->create(comparitorName,
45  comparitorPSet, iC));
46  }
47 }
48 
50  desc.add<double>("extraHitRPhitolerance", 0.06);
51  desc.add<bool>("useBendingCorrection", false);
52  desc.add<double>("CAThetaCut", 0.00125);
53  desc.add<double>("CAPhiCut", 0.1);
54  desc.add<double>("CAHardPtCut", 0);
55 
56  edm::ParameterSetDescription descMaxChi2;
57  descMaxChi2.add<double>("pt1", 0.8);
58  descMaxChi2.add<double>("pt2", 2);
59  descMaxChi2.add<double>("value1", 50);
60  descMaxChi2.add<double>("value2", 8);
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 
74 namespace {
75  void createGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
76  for (unsigned int i = 0; i < layers.size(); i++)
77  {
78  for (unsigned int j = 0; j < 3; ++j)
79  {
80  auto vertexIndex = 0;
81  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
82  layers[i][j].name());
83  if (foundVertex == g.theLayers.end())
84  {
85  g.theLayers.emplace_back(layers[i][j].name(),
86  layers[i][j].hits().size());
87  vertexIndex = g.theLayers.size() - 1;
88  }
89  else
90  {
91  vertexIndex = foundVertex - g.theLayers.begin();
92  }
93  if (j == 0)
94  {
95 
96  if (std::find(g.theRootLayers.begin(), g.theRootLayers.end(),
97  vertexIndex) == g.theRootLayers.end())
98  {
99  g.theRootLayers.emplace_back(vertexIndex);
100 
101  }
102 
103  }
104  }
105  }
106  }
107 
108  void clearGraphStructure(const SeedingLayerSetsHits& layers, CAGraph& g) {
109  g.theLayerPairs.clear();
110  for (unsigned int i = 0; i < g.theLayers.size(); i++ ){
111  g.theLayers[i].theInnerLayers.clear();
112  g.theLayers[i].theInnerLayerPairs.clear();
113  g.theLayers[i].theOuterLayers.clear();
114  g.theLayers[i].theOuterLayerPairs.clear();
115  for (auto & v : g.theLayers[i].isOuterHitOfCell) v.clear();
116  }
117 
118  }
119 
120  void fillGraph(const SeedingLayerSetsHits& layers, const IntermediateHitDoublets::RegionLayerSets& regionLayerPairs,
121  CAGraph& g, std::vector<const HitDoublets *>& hitDoublets) {
122  for (unsigned int i = 0; i < layers.size(); i++)
123  {
124  for (unsigned int j = 0; j < 3; ++j)
125  {
126  auto vertexIndex = 0;
127  auto foundVertex = std::find(g.theLayers.begin(), g.theLayers.end(),
128  layers[i][j].name());
129 
130  if (foundVertex == g.theLayers.end())
131  {
132  vertexIndex = g.theLayers.size() - 1;
133  }
134  else
135  {
136  vertexIndex = foundVertex - g.theLayers.begin();
137  }
138 
139 
140  if (j > 0)
141  {
142 
143  auto innerVertex = std::find(g.theLayers.begin(),
144  g.theLayers.end(), layers[i][j - 1].name());
145 
146  CALayerPair tmpInnerLayerPair(innerVertex - g.theLayers.begin(),
147  vertexIndex);
148 
149  if (std::find(g.theLayerPairs.begin(), g.theLayerPairs.end(),
150  tmpInnerLayerPair) == g.theLayerPairs.end())
151  {
152  auto found = std::find_if(regionLayerPairs.begin(), regionLayerPairs.end(), [&](const IntermediateHitDoublets::LayerPairHitDoublets& pair) {
153  return pair.innerLayerIndex() == layers[i][j - 1].index() && pair.outerLayerIndex() == layers[i][j].index();
154  });
155  if(found != regionLayerPairs.end()) {
156  hitDoublets.emplace_back(&(found->doublets()));
157  g.theLayerPairs.push_back(tmpInnerLayerPair);
158  g.theLayers[vertexIndex].theInnerLayers.push_back(
159  innerVertex - g.theLayers.begin());
160  innerVertex->theOuterLayers.push_back(vertexIndex);
161  g.theLayers[vertexIndex].theInnerLayerPairs.push_back(
162  g.theLayerPairs.size() - 1);
163  innerVertex->theOuterLayerPairs.push_back(
164  g.theLayerPairs.size() - 1);
165  }
166  }
167 
168  }
169 
170  }
171  }
172  }
173 }
174 
176  std::vector<OrderedHitSeeds>& result,
177  const edm::EventSetup& es,
178  const SeedingLayerSetsHits& layers) {
179  CAGraph g;
180 
181  std::vector<const HitDoublets *> hitDoublets;
182 
183  std::vector<CACell::CAntuplet> foundTriplets;
184 
185  int index =0;
186  for(const auto& regionLayerPairs: regionDoublets) {
187 
188  const TrackingRegion& region = regionLayerPairs.region();
189  hitDoublets.clear();
190  foundTriplets.clear();
191 
192  if (index == 0){
193  createGraphStructure(layers, g);
194  }
195  else{
196  clearGraphStructure(layers, g);
197  }
198  fillGraph(layers, regionLayerPairs, g, hitDoublets);
199  CellularAutomaton ca(g);
200  ca.findTriplets(hitDoublets, foundTriplets, region, caThetaCut, caPhiCut,
201  caHardPtCut);
202 
203  auto & allCells = ca.getAllCells();
204 
205  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
206 
207  // re-used thoughout, need to be vectors because of RZLine interface
208  std::array<float, 3> bc_r;
209  std::array<float, 3> bc_z;
210  std::array<float, 3> bc_errZ2;
211  std::array<GlobalPoint, 3> gps;
212  std::array<GlobalError, 3> ges;
213  std::array<bool, 3> barrels;
214 
215  unsigned int numberOfFoundTriplets = foundTriplets.size();
216  for (unsigned int tripletId = 0; tripletId < numberOfFoundTriplets;
217  ++tripletId)
218  {
219 
220  OrderedHitTriplet tmpTriplet(allCells[foundTriplets[tripletId][0]].getInnerHit(),
221  allCells[foundTriplets[tripletId][0]].getOuterHit(),
222  allCells[foundTriplets[tripletId][1]].getOuterHit());
223 
224  auto isBarrel = [](const unsigned id) -> bool
225  {
226  return id == PixelSubdetector::PixelBarrel;
227  };
228  for (unsigned int i = 0; i < 2; ++i)
229  {
230  auto const& ahit = allCells[foundTriplets[tripletId][i]].getInnerHit();
231  gps[i] = ahit->globalPosition();
232  ges[i] = ahit->globalPositionError();
233  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
234  }
235 
236  auto const& ahit = allCells[foundTriplets[tripletId][1]].getOuterHit();
237  gps[2] = ahit->globalPosition();
238  ges[2] = ahit->globalPositionError();
239  barrels[2] = isBarrel(ahit->geographicalId().subdetId());
240 
241  PixelRecoLineRZ line(gps[0], gps[2]);
242  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2],
244  const float curvature = predictionRPhi.curvature(
245  ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
246  const float abscurv = std::abs(curvature);
247  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
248  float chi2 = std::numeric_limits<float>::quiet_NaN();
249  // TODO: Do we have any use case to not use bending correction?
250 
252  {
253  // Following PixelFitterByConformalMappingAndLine
254  const float simpleCot = (gps.back().z() - gps.front().z())
255  / (gps.back().perp() - gps.front().perp());
256  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
257  for (int i = 0; i < 3; ++i)
258  {
259  const GlobalPoint & point = gps[i];
260  const GlobalError & error = ges[i];
261  bc_r[i] = sqrt(
262  sqr(point.x() - region.origin().x())
263  + sqr(point.y() - region.origin().y()));
265  es)(bc_r[i]);
266  bc_z[i] = point.z() - region.origin().z();
267  bc_errZ2[i] =
268  (barrels[i]) ?
269  error.czz() :
270  error.rerr(point) * sqr(simpleCot);
271  }
272  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
273  chi2 = rzLine.chi2();
274  }
275  else
276  {
277  RZLine rzLine(gps, ges, barrels);
278  chi2 = rzLine.chi2();
279  }
280 
281  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
282  {
283  continue;
284 
285  }
286 
287  if (theComparitor)
288  {
289  if (!theComparitor->compatible(tmpTriplet))
290  {
291 
292  continue;
293  }
294  }
295  result[index].emplace_back(tmpTriplet);
296 
297  }
298  index++;
299  }
300 }
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: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
*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