CMS 3D CMS Logo

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

#include <CAHitQuadrupletGenerator.h>

Classes

class  QuantityDependsPt
 
class  QuantityDependsPtEval
 

Public Types

typedef LayerHitMapCache LayerCacheType
 
typedef OrderedHitSeeds ResultType
 

Public Member Functions

 CAHitQuadrupletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
 CAHitQuadrupletGenerator (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)
 
 ~CAHitQuadrupletGenerator ()=default
 

Static Public Member Functions

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

Static Public Attributes

static unsigned int minLayers = 4
 

Private Attributes

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

Detailed Description

Definition at line 31 of file CAHitQuadrupletGenerator.h.

Member Typedef Documentation

Definition at line 33 of file CAHitQuadrupletGenerator.h.

Definition at line 36 of file CAHitQuadrupletGenerator.h.

Constructor & Destructor Documentation

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

Definition at line 40 of file CAHitQuadrupletGenerator.h.

References looper::cfg.

40 : CAHitQuadrupletGenerator(cfg, iC) {}
CAHitQuadrupletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
CAHitQuadrupletGenerator::CAHitQuadrupletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 38 of file CAHitQuadrupletGenerator.cc.

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

38  :
39 extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
40 maxChi2(cfg.getParameter<edm::ParameterSet>("maxChi2")),
41 fitFastCircle(cfg.getParameter<bool>("fitFastCircle")),
42 fitFastCircleChi2Cut(cfg.getParameter<bool>("fitFastCircleChi2Cut")),
43 useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
44 caThetaCut(cfg.getParameter<double>("CAThetaCut")),
45 caPhiCut(cfg.getParameter<double>("CAPhiCut")),
46 caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
47 {
48  edm::ParameterSet comparitorPSet = cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
49  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
50  if (comparitorName != "none")
51  {
52  theComparitor.reset(SeedComparitorFactory::get()->create(comparitorName, comparitorPSet, iC));
53  }
54 }
T getParameter(std::string const &) const
def create(alignables, pedeDump, additionalData, outputFile, config)
std::unique_ptr< SeedComparitor > theComparitor
const QuantityDependsPt maxChi2
T get(const Candidate &c)
Definition: component.h:55
CAHitQuadrupletGenerator::~CAHitQuadrupletGenerator ( )
default

Member Function Documentation

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

Definition at line 56 of file CAHitQuadrupletGenerator.cc.

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

56  {
57  desc.add<double>("extraHitRPhitolerance", 0.1);
58  desc.add<bool>("fitFastCircle", false);
59  desc.add<bool>("fitFastCircleChi2Cut", false);
60  desc.add<bool>("useBendingCorrection", false);
61  desc.add<double>("CAThetaCut", 0.00125);
62  desc.add<double>("CAPhiCut", 10);
63  desc.add<double>("CAHardPtCut", 0);
64  desc.addOptional<bool>("CAOnlyOneLastHitPerLayerFilter")->setComment("Deprecated and has no effect. To be fully removed later when the parameter is no longer used in HLT configurations.");
65  edm::ParameterSetDescription descMaxChi2;
66  descMaxChi2.add<double>("pt1", 0.2);
67  descMaxChi2.add<double>("pt2", 1.5);
68  descMaxChi2.add<double>("value1", 500);
69  descMaxChi2.add<double>("value2", 50);
70  descMaxChi2.add<bool>("enabled", true);
71  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
72 
73  edm::ParameterSetDescription descComparitor;
74  descComparitor.add<std::string>("ComponentName", "none");
75  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
76  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
77 }
ParameterDescriptionBase * addOptional(U const &iLabel, T const &value)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static const char* CAHitQuadrupletGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 46 of file CAHitQuadrupletGenerator.h.

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

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

Definition at line 180 of file CAHitQuadrupletGenerator.cc.

References funct::abs(), EnergyCorrector::c, caHardPtCut, caPhiCut, caThetaCut, vertices_cff::chi2, FastCircleFit::chi2(), RZLine::chi2(), CellularAutomaton::createAndConnectCells(), PixelRecoUtilities::curvature(), GlobalErrorBase< T, ErrorWeightType >::czz(), relativeConstraints::error, CAHitQuadrupletGenerator::QuantityDependsPt::evaluator(), CellularAutomaton::evolve(), extraHitRPhitolerance, CellularAutomaton::findNtuplets(), fitFastCircle, fitFastCircleChi2Cut, g, CellularAutomaton::getAllCells(), mps_fire::i, triggerObjects_cff::id, diffTreeTool::index, PixelRecoUtilities::inversePt(), gedGsfElectrons_cfi::isBarrel, edm::isNotFinite(), maxChi2, TrackingRegion::origin(), PixelSubdetector::PixelBarrel, point, EnergyCorrector::pt, GlobalErrorBase< T, ErrorWeightType >::rerr(), funct::sqr(), mathSSE::sqrt(), theComparitor, useBendingCorrection, CAHitQuadrupletGenerator::QuantityDependsPtEval::value(), x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

183  {
184  CAGraph g;
185 
186  std::vector<const HitDoublets *> hitDoublets;
187 
188  const int numberOfHitsInNtuplet = 4;
189  std::vector<CACell::CAntuplet> foundQuadruplets;
190 
191 
192 
193  int index =0;
194  for(const auto& regionLayerPairs: regionDoublets) {
195 
196  const TrackingRegion& region = regionLayerPairs.region();
197  hitDoublets.clear();
198  foundQuadruplets.clear();
199  if (index == 0){
200  createGraphStructure(layers, g);
201  }
202  else{
203  clearGraphStructure(layers, g);
204  }
205 
206  fillGraph(layers, regionLayerPairs, g, hitDoublets);
207 
208  CellularAutomaton ca(g);
209 
210  ca.createAndConnectCells(hitDoublets, region, caThetaCut,
212 
213  ca.evolve(numberOfHitsInNtuplet);
214 
215  ca.findNtuplets(foundQuadruplets, numberOfHitsInNtuplet);
216 
217  auto & allCells = ca.getAllCells();
218 
219  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
220 
221  // re-used thoughout
222  std::array<float, 4> bc_r;
223  std::array<float, 4> bc_z;
224  std::array<float, 4> bc_errZ2;
225  std::array<GlobalPoint, 4> gps;
226  std::array<GlobalError, 4> ges;
227  std::array<bool, 4> barrels;
228 
229  unsigned int numberOfFoundQuadruplets = foundQuadruplets.size();
230 
231  // Loop over quadruplets
232  for (unsigned int quadId = 0; quadId < numberOfFoundQuadruplets; ++quadId)
233  {
234 
235  auto isBarrel = [](const unsigned id) -> bool
236  {
237  return id == PixelSubdetector::PixelBarrel;
238  };
239  for(unsigned int i = 0; i< 3; ++i)
240  {
241  auto const& ahit = allCells[foundQuadruplets[quadId][i]].getInnerHit();
242  gps[i] = ahit->globalPosition();
243  ges[i] = ahit->globalPositionError();
244  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
245  }
246 
247  auto const& ahit = allCells[foundQuadruplets[quadId][2]].getOuterHit();
248  gps[3] = ahit->globalPosition();
249  ges[3] = ahit->globalPositionError();
250  barrels[3] = isBarrel(ahit->geographicalId().subdetId());
251  // TODO:
252  // - if we decide to always do the circle fit for 4 hits, we don't
253  // need ThirdHitPredictionFromCircle for the curvature; then we
254  // could remove extraHitRPhitolerance configuration parameter
255  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2], extraHitRPhitolerance);
256  const float curvature = predictionRPhi.curvature(ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
257  const float abscurv = std::abs(curvature);
258  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
259  if (theComparitor)
260  {
261  SeedingHitSet tmpTriplet(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
262  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
263  allCells[foundQuadruplets[quadId][2]].getOuterHit());
264 
265 
266  if (!theComparitor->compatible(tmpTriplet) )
267  {
268  continue;
269  }
270  }
271 
272  float chi2 = std::numeric_limits<float>::quiet_NaN();
273  // TODO: Do we have any use case to not use bending correction?
275  {
276  // Following PixelFitterByConformalMappingAndLine
277  const float simpleCot = ( gps.back().z() - gps.front().z() ) / (gps.back().perp() - gps.front().perp() );
278  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
279  for (int i=0; i < 4; ++i)
280  {
281  const GlobalPoint & point = gps[i];
282  const GlobalError & error = ges[i];
283  bc_r[i] = sqrt( sqr(point.x() - region.origin().x()) + sqr(point.y() - region.origin().y()) );
284  bc_r[i] += pixelrecoutilities::LongitudinalBendingCorrection(pt, es)(bc_r[i]);
285  bc_z[i] = point.z() - region.origin().z();
286  bc_errZ2[i] = (barrels[i]) ? error.czz() : error.rerr(point)*sqr(simpleCot);
287  }
288  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
289  chi2 = rzLine.chi2();
290  }
291  else
292  {
293  RZLine rzLine(gps, ges, barrels);
294  chi2 = rzLine.chi2();
295  }
296  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
297  {
298  continue;
299  }
300  // TODO: Do we have any use case to not use circle fit? Maybe
301  // HLT where low-pT inefficiency is not a problem?
302  if (fitFastCircle)
303  {
304  FastCircleFit c(gps, ges);
305  chi2 += c.chi2();
306  if (edm::isNotFinite(chi2)) continue;
307  if (fitFastCircleChi2Cut && chi2 > thisMaxChi2) continue;
308  }
309  result[index].emplace_back(allCells[foundQuadruplets[quadId][0]].getInnerHit(),
310  allCells[foundQuadruplets[quadId][1]].getInnerHit(),
311  allCells[foundQuadruplets[quadId][2]].getInnerHit(),
312  allCells[foundQuadruplets[quadId][2]].getOuterHit());
313  }
314  index++;
315  }
316 
317 }
GlobalPoint const & origin() const
T y() const
Definition: PV3DBase.h:63
std::unique_ptr< SeedComparitor > theComparitor
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)
bool isNotFinite(T x)
Definition: isFinite.h:10
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
const QuantityDependsPt maxChi2
T rerr(const GlobalPoint &aPoint) const
Square< F >::type sqr(const F &f)
Definition: Square.h:13
QuantityDependsPtEval evaluator(const edm::EventSetup &es) const
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 CAHitQuadrupletGenerator::initEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Member Data Documentation

const float CAHitQuadrupletGenerator::caHardPtCut = 0.f
private

Definition at line 133 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::caPhiCut = 0.1f
private

Definition at line 132 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::caThetaCut = 0.00125f
private

Definition at line 131 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::extraHitRPhitolerance
private

Definition at line 124 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const bool CAHitQuadrupletGenerator::fitFastCircle
private

Definition at line 127 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const bool CAHitQuadrupletGenerator::fitFastCircleChi2Cut
private

Definition at line 128 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const QuantityDependsPt CAHitQuadrupletGenerator::maxChi2
private

Definition at line 126 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

unsigned int CAHitQuadrupletGenerator::minLayers = 4
static

Definition at line 35 of file CAHitQuadrupletGenerator.h.

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

Definition at line 58 of file CAHitQuadrupletGenerator.h.

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

LayerCacheType CAHitQuadrupletGenerator::theLayerCache
private

Definition at line 56 of file CAHitQuadrupletGenerator.h.

const bool CAHitQuadrupletGenerator::useBendingCorrection
private

Definition at line 129 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().