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 29 of file CAHitQuadrupletGenerator.h.

Member Typedef Documentation

Definition at line 31 of file CAHitQuadrupletGenerator.h.

Definition at line 34 of file CAHitQuadrupletGenerator.h.

Constructor & Destructor Documentation

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

Definition at line 38 of file CAHitQuadrupletGenerator.h.

References looper::cfg.

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

Definition at line 40 of file CAHitQuadrupletGenerator.cc.

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

40  :
41 extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
42 maxChi2(cfg.getParameter<edm::ParameterSet>("maxChi2")),
43 fitFastCircle(cfg.getParameter<bool>("fitFastCircle")),
44 fitFastCircleChi2Cut(cfg.getParameter<bool>("fitFastCircleChi2Cut")),
45 useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
46 caThetaCut(cfg.getParameter<double>("CAThetaCut")),
47 caPhiCut(cfg.getParameter<double>("CAPhiCut")),
48 caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
49 {
50  edm::ParameterSet comparitorPSet = cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
51  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
52  if (comparitorName != "none")
53  {
54  theComparitor.reset(SeedComparitorFactory::get()->create(comparitorName, comparitorPSet, iC));
55  }
56 }
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 58 of file CAHitQuadrupletGenerator.cc.

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

58  {
59  desc.add<double>("extraHitRPhitolerance", 0.1);
60  desc.add<bool>("fitFastCircle", false);
61  desc.add<bool>("fitFastCircleChi2Cut", false);
62  desc.add<bool>("useBendingCorrection", false);
63  desc.add<double>("CAThetaCut", 0.00125);
64  desc.add<double>("CAPhiCut", 10);
65  desc.add<double>("CAHardPtCut", 0);
66  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.");
67  edm::ParameterSetDescription descMaxChi2;
68  descMaxChi2.add<double>("pt1", 0.2);
69  descMaxChi2.add<double>("pt2", 1.5);
70  descMaxChi2.add<double>("value1", 500);
71  descMaxChi2.add<double>("value2", 50);
72  descMaxChi2.add<bool>("enabled", true);
73  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
74 
75  edm::ParameterSetDescription descComparitor;
76  descComparitor.add<std::string>("ComponentName", "none");
77  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
78  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
79 }
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 44 of file CAHitQuadrupletGenerator.h.

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

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

Definition at line 182 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().

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

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::caPhiCut = 0.1f
private

Definition at line 130 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::caThetaCut = 0.00125f
private

Definition at line 129 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const float CAHitQuadrupletGenerator::extraHitRPhitolerance
private

Definition at line 122 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const bool CAHitQuadrupletGenerator::fitFastCircle
private

Definition at line 125 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const bool CAHitQuadrupletGenerator::fitFastCircleChi2Cut
private

Definition at line 126 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

const QuantityDependsPt CAHitQuadrupletGenerator::maxChi2
private

Definition at line 124 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().

unsigned int CAHitQuadrupletGenerator::minLayers = 4
static

Definition at line 33 of file CAHitQuadrupletGenerator.h.

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

Definition at line 56 of file CAHitQuadrupletGenerator.h.

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

LayerCacheType CAHitQuadrupletGenerator::theLayerCache
private

Definition at line 54 of file CAHitQuadrupletGenerator.h.

const bool CAHitQuadrupletGenerator::useBendingCorrection
private

Definition at line 127 of file CAHitQuadrupletGenerator.h.

Referenced by hitNtuplets().