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 41 of file CAHitQuadrupletGenerator.cc.

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

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

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

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

References funct::abs(), EnergyCorrector::c, caHardPtCut, caPhiCut, caThetaCut, HiEvtPlane_cfi::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, hcalTTPDigis_cfi::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().

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