CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes
CAHitTripletGenerator Class Reference

#include <CAHitTripletGenerator.h>

Classes

class  QuantityDependsPt
 
class  QuantityDependsPtEval
 

Public Types

typedef LayerHitMapCache LayerCacheType
 
typedef OrderedHitSeeds ResultType
 

Public Member Functions

 CAHitTripletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
 CAHitTripletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
void hitNtuplets (const IntermediateHitDoublets &regionDoublets, std::vector< OrderedHitSeeds > &result, const edm::EventSetup &es, const SeedingLayerSetsHits &layers)
 
void initEvent (const edm::Event &ev, const edm::EventSetup &es)
 
 ~CAHitTripletGenerator ()=default
 

Static Public Member Functions

static void fillDescriptions (edm::ParameterSetDescription &desc)
 
static const char * fillDescriptionsLabel ()
 

Static Public Attributes

static unsigned int minLayers = 3
 

Private Attributes

const float caHardPtCut = 0.f
 
const float caPhiCut = 1.f
 
const float caThetaCut = 0.00125f
 
const float extraHitRPhitolerance
 
const QuantityDependsPt maxChi2
 
std::unique_ptr< SeedComparitortheComparitor
 
LayerCacheType theLayerCache
 
const bool useBendingCorrection
 

Detailed Description

Definition at line 29 of file CAHitTripletGenerator.h.

Member Typedef Documentation

Definition at line 31 of file CAHitTripletGenerator.h.

Definition at line 34 of file CAHitTripletGenerator.h.

Constructor & Destructor Documentation

CAHitTripletGenerator::CAHitTripletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 37 of file CAHitTripletGenerator.h.

References looper::cfg, and fillDescriptions().

37 : CAHitTripletGenerator(cfg, iC) {}
CAHitTripletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
CAHitTripletGenerator::CAHitTripletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 30 of file CAHitTripletGenerator.cc.

References reco::get(), edm::ParameterSet::getParameter(), AlCaHLTBitMon_QueryRunRegistry::string, and theComparitor.

31  :
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 }
T getParameter(std::string const &) const
const QuantityDependsPt maxChi2
std::unique_ptr< SeedComparitor > theComparitor
T get(const Candidate &c)
Definition: component.h:55
CAHitTripletGenerator::~CAHitTripletGenerator ( )
default

Member Function Documentation

void CAHitTripletGenerator::fillDescriptions ( edm::ParameterSetDescription desc)
static

Definition at line 48 of file CAHitTripletGenerator.cc.

References edm::ParameterSetDescription::add(), and AlCaHLTBitMon_QueryRunRegistry::string.

48  {
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 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static const char* CAHitTripletGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 43 of file CAHitTripletGenerator.h.

References ev, LayerTriplets::layers(), and mps_fire::result.

43 { return "caHitTriplet"; }
void CAHitTripletGenerator::hitNtuplets ( const IntermediateHitDoublets regionDoublets,
std::vector< OrderedHitSeeds > &  result,
const edm::EventSetup es,
const SeedingLayerSetsHits layers 
)

Definition at line 174 of file CAHitTripletGenerator.cc.

References funct::abs(), caHardPtCut, caPhiCut, caThetaCut, vertices_cff::chi2, RZLine::chi2(), ThirdHitPredictionFromCircle::curvature(), PixelRecoUtilities::curvature(), GlobalErrorBase< T, ErrorWeightType >::czz(), relativeConstraints::error, CAHitTripletGenerator::QuantityDependsPt::evaluator(), extraHitRPhitolerance, CellularAutomaton::findTriplets(), g, CellularAutomaton::getAllCells(), mps_fire::i, triggerObjects_cff::id, PixelRecoUtilities::inversePt(), GeomDetEnumerators::isBarrel(), edm::isNotFinite(), mps_splice::line, maxChi2, TrackingRegion::origin(), PixelSubdetector::PixelBarrel, point, EnergyCorrector::pt, GlobalErrorBase< T, ErrorWeightType >::rerr(), funct::sqr(), mathSSE::sqrt(), theComparitor, useBendingCorrection, CAHitTripletGenerator::QuantityDependsPtEval::value(), x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

177  {
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 }
const QuantityDependsPt maxChi2
GlobalPoint const & origin() const
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
bool isBarrel(GeomDetEnumerators::SubDetector m)
T y() const
Definition: PV3DBase.h:63
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
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
T z() const
Definition: PV3DBase.h:64
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Definition: RZLine.h:12
T rerr(const GlobalPoint &aPoint) const
Square< F >::type sqr(const F &f)
Definition: Square.h:13
std::unique_ptr< SeedComparitor > theComparitor
T x() const
Definition: PV3DBase.h:62
*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 CAHitTripletGenerator::initEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Member Data Documentation

const float CAHitTripletGenerator::caHardPtCut = 0.f
private

Definition at line 128 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::caPhiCut = 1.f
private

Definition at line 127 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::caThetaCut = 0.00125f
private

Definition at line 126 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::extraHitRPhitolerance
private

Definition at line 121 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const QuantityDependsPt CAHitTripletGenerator::maxChi2
private

Definition at line 123 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

unsigned int CAHitTripletGenerator::minLayers = 3
static

Definition at line 33 of file CAHitTripletGenerator.h.

std::unique_ptr<SeedComparitor> CAHitTripletGenerator::theComparitor
private

Definition at line 55 of file CAHitTripletGenerator.h.

Referenced by CAHitTripletGenerator(), hitNtuplets(), and initEvent().

LayerCacheType CAHitTripletGenerator::theLayerCache
private

Definition at line 53 of file CAHitTripletGenerator.h.

const bool CAHitTripletGenerator::useBendingCorrection
private

Definition at line 124 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().