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, and fillDescriptions().

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

Definition at line 29 of file CAHitQuadrupletGenerator.cc.

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

29  :
30 extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
31 maxChi2(cfg.getParameter<edm::ParameterSet>("maxChi2")),
32 fitFastCircle(cfg.getParameter<bool>("fitFastCircle")),
33 fitFastCircleChi2Cut(cfg.getParameter<bool>("fitFastCircleChi2Cut")),
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, comparitorPSet, iC)};
44  }
45 }
T getParameter(std::string const &) const
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 47 of file CAHitQuadrupletGenerator.cc.

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

47  {
48  desc.add<double>("extraHitRPhitolerance", 0.1);
49  desc.add<bool>("fitFastCircle", false);
50  desc.add<bool>("fitFastCircleChi2Cut", false);
51  desc.add<bool>("useBendingCorrection", false);
52  desc.add<double>("CAThetaCut", 0.00125);
53  desc.add<double>("CAPhiCut", 10);
54  desc.add<double>("CAHardPtCut", 0);
55  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.");
56  edm::ParameterSetDescription descMaxChi2;
57  descMaxChi2.add<double>("pt1", 0.2);
58  descMaxChi2.add<double>("pt2", 1.5);
59  descMaxChi2.add<double>("value1", 500);
60  descMaxChi2.add<double>("value2", 50);
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 }
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 171 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, PixelRecoUtilities::inversePt(), GeomDetEnumerators::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().

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