test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes
CAHitQuadrupletGenerator Class Reference

#include <CAHitQuadrupletGenerator.h>

Inheritance diagram for CAHitQuadrupletGenerator:
HitQuadrupletGenerator OrderedHitsGenerator

Classes

class  QuantityDependsPt
 
class  QuantityDependsPtEval
 

Public Types

typedef LayerHitMapCache LayerCacheType
 
typedef OrderedHitSeeds ResultType
 

Public Member Functions

 CAHitQuadrupletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC, bool needSeedingLayerSetsHits=true)
 
 CAHitQuadrupletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC, bool needSeedingLayerSetsHits=true)
 
void hitNtuplets (const IntermediateHitDoublets::RegionLayerSets &regionLayerPairs, OrderedHitSeeds &result, const edm::EventSetup &es, const SeedingLayerSetsHits &layers)
 
virtual void hitQuadruplets (const TrackingRegion &reg, OrderedHitSeeds &quadruplets, const edm::Event &ev, const edm::EventSetup &es)
 from base class More...
 
void initEvent (const edm::Event &ev, const edm::EventSetup &es)
 
virtual ~CAHitQuadrupletGenerator ()
 
- Public Member Functions inherited from HitQuadrupletGenerator
virtual void clear () final
 
 HitQuadrupletGenerator (unsigned int size=500)
 
virtual const OrderedHitSeedsrun (const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es) final
 
virtual ~HitQuadrupletGenerator ()
 
- Public Member Functions inherited from OrderedHitsGenerator
 OrderedHitsGenerator ()
 
virtual ~OrderedHitsGenerator ()
 

Static Public Member Functions

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

Static Public Attributes

static unsigned int minLayers = 4
 

Private Member Functions

void hitQuadruplets (const TrackingRegion &reg, OrderedHitSeeds &result, std::vector< const HitDoublets * > &hitDoublets, const CAGraph &g, const edm::EventSetup &es)
 

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
 
edm::EDGetTokenT
< SeedingLayerSetsHits
theSeedingLayerToken
 
const bool useBendingCorrection
 

Additional Inherited Members

- Public Attributes inherited from OrderedHitsGenerator
unsigned int theMaxElement
 

Detailed Description

Definition at line 32 of file CAHitQuadrupletGenerator.h.

Member Typedef Documentation

Definition at line 34 of file CAHitQuadrupletGenerator.h.

Definition at line 37 of file CAHitQuadrupletGenerator.h.

Constructor & Destructor Documentation

CAHitQuadrupletGenerator::CAHitQuadrupletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC,
bool  needSeedingLayerSetsHits = true 
)
inline

Definition at line 41 of file CAHitQuadrupletGenerator.h.

41 : CAHitQuadrupletGenerator(cfg, iC, needSeedingLayerSetsHits) {}
CAHitQuadrupletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC, bool needSeedingLayerSetsHits=true)
CAHitQuadrupletGenerator::CAHitQuadrupletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC,
bool  needSeedingLayerSetsHits = true 
)

Definition at line 42 of file CAHitQuadrupletGenerator.cc.

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

42  :
43 extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
44 maxChi2(cfg.getParameter<edm::ParameterSet>("maxChi2")),
45 fitFastCircle(cfg.getParameter<bool>("fitFastCircle")),
46 fitFastCircleChi2Cut(cfg.getParameter<bool>("fitFastCircleChi2Cut")),
47 useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
48 caThetaCut(cfg.getParameter<double>("CAThetaCut")),
49 caPhiCut(cfg.getParameter<double>("CAPhiCut")),
50 caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
51 {
52  if(needSeedingLayerSetsHits)
54 
55  if (cfg.exists("SeedComparitorPSet"))
56  {
57  edm::ParameterSet comparitorPSet =
58  cfg.getParameter<edm::ParameterSet>("SeedComparitorPSet");
59  std::string comparitorName = comparitorPSet.getParameter<std::string>("ComponentName");
60  if (comparitorName != "none")
61  {
62  theComparitor.reset(SeedComparitorFactory::get()->create(comparitorName, comparitorPSet, iC));
63  }
64  }
65 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::unique_ptr< SeedComparitor > theComparitor
const QuantityDependsPt maxChi2
edm::EDGetTokenT< SeedingLayerSetsHits > theSeedingLayerToken
T get(const Candidate &c)
Definition: component.h:55
CAHitQuadrupletGenerator::~CAHitQuadrupletGenerator ( )
virtual

Definition at line 67 of file CAHitQuadrupletGenerator.cc.

68 {
69 }

Member Function Documentation

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

Definition at line 71 of file CAHitQuadrupletGenerator.cc.

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

71  {
72  desc.add<double>("extraHitRPhitolerance", 0.1);
73  desc.add<bool>("fitFastCircle", false);
74  desc.add<bool>("fitFastCircleChi2Cut", false);
75  desc.add<bool>("useBendingCorrection", false);
76  desc.add<double>("CAThetaCut", 0.00125);
77  desc.add<double>("CAPhiCut", 10);
78  desc.add<double>("CAHardPtCut", 0);
79 
80  edm::ParameterSetDescription descMaxChi2;
81  descMaxChi2.add<double>("pt1", 0.2);
82  descMaxChi2.add<double>("pt2", 1.5);
83  descMaxChi2.add<double>("value1", 500);
84  descMaxChi2.add<double>("value2", 50);
85  descMaxChi2.add<bool>("enabled", true);
86  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
87 
88  edm::ParameterSetDescription descComparitor;
89  descComparitor.add<std::string>("ComponentName", "none");
90  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
91  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
92 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static const char* CAHitQuadrupletGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 47 of file CAHitQuadrupletGenerator.h.

47 { return "caHitQuadruplet"; }
void CAHitQuadrupletGenerator::hitNtuplets ( const IntermediateHitDoublets::RegionLayerSets regionLayerPairs,
OrderedHitSeeds result,
const edm::EventSetup es,
const SeedingLayerSetsHits layers 
)

Definition at line 201 of file CAHitQuadrupletGenerator.cc.

References ihd::RegionLayerSets< T >::begin(), ihd::RegionLayerSets< T >::end(), newFWLiteAna::found, g, hitQuadruplets(), SeedingLayerSetsHits::SeedingLayer::index(), SurfaceOrientation::inner, IntermediateHitDoublets::LayerPairHitDoublets::innerLayerIndex(), SurfaceOrientation::outer, IntermediateHitDoublets::LayerPairHitDoublets::outerLayerIndex(), ihd::RegionLayerSets< T >::region(), and mps_fire::result.

204  {
205  CAGraph g;
206 
207  std::vector<const HitDoublets *> hitDoublets;
208 
209  auto layerPairEqual = [](const IntermediateHitDoublets::LayerPairHitDoublets& pair,
212  return pair.innerLayerIndex() == inner && pair.outerLayerIndex() == outer;
213  };
214  fillGraph(layers, g, hitDoublets,
215  [&](const SeedingLayerSetsHits::SeedingLayer& inner,
217  std::vector<const HitDoublets *>& hitDoublets) {
218  using namespace std::placeholders;
219  auto found = std::find_if(regionLayerPairs.begin(), regionLayerPairs.end(),
220  std::bind(layerPairEqual, _1, inner.index(), outer.index()));
221  if(found != regionLayerPairs.end()) {
222  hitDoublets.emplace_back(&(found->doublets()));
223  return true;
224  }
225  return false;
226  });
227 
228  hitQuadruplets(regionLayerPairs.region(), result, hitDoublets, g, es);
229 }
const_iterator begin() const
const_iterator end() const
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
tuple result
Definition: mps_fire.py:84
virtual void hitQuadruplets(const TrackingRegion &reg, OrderedHitSeeds &quadruplets, const edm::Event &ev, const edm::EventSetup &es)
from base class
SeedingLayerSetsHits::LayerIndex innerLayerIndex() const
const TrackingRegion & region() const
SeedingLayerSetsHits::LayerIndex outerLayerIndex() const
void CAHitQuadrupletGenerator::hitQuadruplets ( const TrackingRegion reg,
OrderedHitSeeds quadruplets,
const edm::Event ev,
const edm::EventSetup es 
)
virtual

from base class

Implements HitQuadrupletGenerator.

Definition at line 162 of file CAHitQuadrupletGenerator.cc.

References LayerHitMapCache::clear(), HitPairGeneratorFromLayerPair::doublets(), alignCSCRings::e, g, edm::Event::getByToken(), SurfaceOrientation::inner, SeedingLayerSetsHits::numberOfLayersInSet(), SurfaceOrientation::outer, theComparitor, theLayerCache, and theSeedingLayerToken.

Referenced by hitNtuplets().

165 {
167  ev.getByToken(theSeedingLayerToken, hlayers);
168  const SeedingLayerSetsHits& layers = *hlayers;
169  if (layers.numberOfLayersInSet() != 4)
170  throw cms::Exception("Configuration")
171  << "CAHitQuadrupletsGenerator expects SeedingLayerSetsHits::numberOfLayersInSet() to be 4, got "
172  << layers.numberOfLayersInSet();
173 
174  CAGraph g;
175 
176 
177  std::vector<HitDoublets> hitDoublets;
178 
179 
180  HitPairGeneratorFromLayerPair thePairGenerator(0, 1, &theLayerCache);
181  fillGraph(layers, g, hitDoublets,
184  std::vector<HitDoublets>& hitDoublets) {
185  hitDoublets.emplace_back(thePairGenerator.doublets(region, ev, es, inner, outer));
186  return true;
187  });
188 
189  if (theComparitor)
190  theComparitor->init(ev, es);
191 
192  std::vector<const HitDoublets *> hitDoubletsPtr;
193  hitDoubletsPtr.reserve(hitDoublets.size());
194  for(const auto& e: hitDoublets)
195  hitDoubletsPtr.emplace_back(&e);
196 
197  hitQuadruplets(region, result, hitDoubletsPtr, g, es);
199 }
unsigned short numberOfLayersInSet() const
Get number of layers in each SeedingLayerSets.
std::vector< LayerSetAndLayers > layers(const SeedingLayerSetsHits &sets)
Definition: LayerTriplets.cc:4
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
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
tuple result
Definition: mps_fire.py:84
virtual void hitQuadruplets(const TrackingRegion &reg, OrderedHitSeeds &quadruplets, const edm::Event &ev, const edm::EventSetup &es)
from base class
edm::EDGetTokenT< SeedingLayerSetsHits > theSeedingLayerToken
void CAHitQuadrupletGenerator::hitQuadruplets ( const TrackingRegion reg,
OrderedHitSeeds result,
std::vector< const HitDoublets * > &  hitDoublets,
const CAGraph g,
const edm::EventSetup es 
)
private

Definition at line 231 of file CAHitQuadrupletGenerator.cc.

References funct::abs(), EnergyCorrector::c, caHardPtCut, caPhiCut, caThetaCut, FastCircleFit::chi2(), RZLine::chi2(), beam_dqm_sourceclient-live_cfg::chi2, CellularAutomaton::createAndConnectCells(), ThirdHitPredictionFromCircle::curvature(), PixelRecoUtilities::curvature(), GlobalErrorBase< T, ErrorWeightType >::czz(), relativeConstraints::error, CAHitQuadrupletGenerator::QuantityDependsPt::evaluator(), CellularAutomaton::evolve(), extraHitRPhitolerance, CellularAutomaton::findNtuplets(), fitFastCircle, fitFastCircleChi2Cut, i, PixelRecoUtilities::inversePt(), GeomDetEnumerators::isBarrel(), edm::isNotFinite(), geometryCSVtoXML::line, 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().

234  {
235  const int numberOfHitsInNtuplet = 4;
236  std::vector<CACell::CAntuplet> foundQuadruplets;
237 
238  CellularAutomaton ca(g);
239 
240  ca.createAndConnectCells(hitDoublets, region, caThetaCut,
242 
243  ca.evolve(numberOfHitsInNtuplet);
244 
245  ca.findNtuplets(foundQuadruplets, numberOfHitsInNtuplet);
246 
247 
248  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
249 
250  // re-used thoughout, need to be vectors because of RZLine interface
251  std::array<float, 4> bc_r;
252  std::array<float, 4> bc_z;
253  std::array<float, 4> bc_errZ2;
254  std::array<GlobalPoint, 4> gps;
255  std::array<GlobalError, 4> ges;
256  std::array<bool, 4> barrels;
257 
258  unsigned int numberOfFoundQuadruplets = foundQuadruplets.size();
259 
260  // Loop over quadruplets
261  for (unsigned int quadId = 0; quadId < numberOfFoundQuadruplets; ++quadId)
262  {
263 
264  auto isBarrel = [](const unsigned id) -> bool
265  {
266  return id == PixelSubdetector::PixelBarrel;
267  };
268  for(unsigned int i = 0; i< 3; ++i)
269  {
270  auto const& ahit = foundQuadruplets[quadId][i]->getInnerHit();
271  gps[i] = ahit->globalPosition();
272  ges[i] = ahit->globalPositionError();
273  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
274  }
275 
276  auto const& ahit = foundQuadruplets[quadId][2]->getOuterHit();
277  gps[3] = ahit->globalPosition();
278  ges[3] = ahit->globalPositionError();
279  barrels[3] = isBarrel(ahit->geographicalId().subdetId());
280 
281  // TODO:
282  // - 'line' is not used for anything
283  // - if we decide to always do the circle fit for 4 hits, we don't
284  // need ThirdHitPredictionFromCircle for the curvature; then we
285  // could remove extraHitRPhitolerance configuration parameter
286  PixelRecoLineRZ line(gps[0], gps[2]);
287  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2], extraHitRPhitolerance);
288  const float curvature = predictionRPhi.curvature(ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
289  const float abscurv = std::abs(curvature);
290  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
291 
292  if (theComparitor)
293  {
294  SeedingHitSet tmpTriplet(foundQuadruplets[quadId][0]->getInnerHit(), foundQuadruplets[quadId][2]->getInnerHit(), foundQuadruplets[quadId][2]->getOuterHit());
295 
296 
297  if (!theComparitor->compatible(tmpTriplet) )
298  {
299  continue;
300  }
301  }
302 
303  float chi2 = std::numeric_limits<float>::quiet_NaN();
304  // TODO: Do we have any use case to not use bending correction?
306  {
307  // Following PixelFitterByConformalMappingAndLine
308  const float simpleCot = ( gps.back().z() - gps.front().z() ) / (gps.back().perp() - gps.front().perp() );
309  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
310  for (int i=0; i < 4; ++i)
311  {
312  const GlobalPoint & point = gps[i];
313  const GlobalError & error = ges[i];
314  bc_r[i] = sqrt( sqr(point.x() - region.origin().x()) + sqr(point.y() - region.origin().y()) );
315  bc_r[i] += pixelrecoutilities::LongitudinalBendingCorrection(pt, es)(bc_r[i]);
316  bc_z[i] = point.z() - region.origin().z();
317  bc_errZ2[i] = (barrels[i]) ? error.czz() : error.rerr(point)*sqr(simpleCot);
318  }
319  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
320  chi2 = rzLine.chi2();
321  } else
322  {
323  RZLine rzLine(gps, ges, barrels);
324  chi2 = rzLine.chi2();
325  }
326  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
327  {
328  continue;
329 
330  }
331  // TODO: Do we have any use case to not use circle fit? Maybe
332  // HLT where low-pT inefficiency is not a problem?
333  if (fitFastCircle)
334  {
335  FastCircleFit c(gps, ges);
336  chi2 += c.chi2();
337  if (edm::isNotFinite(chi2))
338  continue;
339  if (fitFastCircleChi2Cut && chi2 > thisMaxChi2)
340  continue;
341  }
342 
343  result.emplace_back(foundQuadruplets[quadId][0]->getInnerHit(), foundQuadruplets[quadId][1]->getInnerHit(), foundQuadruplets[quadId][2]->getInnerHit(), foundQuadruplets[quadId][2]->getOuterHit());
344  }
345 }
int i
Definition: DBlmapReader.cc:9
bool isBarrel(GeomDetEnumerators::SubDetector m)
T y() const
Definition: PV3DBase.h:63
std::unique_ptr< SeedComparitor > theComparitor
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 
)

Definition at line 94 of file CAHitQuadrupletGenerator.cc.

References theComparitor.

94  {
95  if (theComparitor) theComparitor->init(ev, es);
96 }
std::unique_ptr< SeedComparitor > theComparitor

Member Data Documentation

const float CAHitQuadrupletGenerator::caHardPtCut = 0.f
private

Definition at line 148 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

const float CAHitQuadrupletGenerator::caPhiCut = 0.1f
private

Definition at line 147 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

const float CAHitQuadrupletGenerator::caThetaCut = 0.00125f
private

Definition at line 146 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

const float CAHitQuadrupletGenerator::extraHitRPhitolerance
private

Definition at line 139 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

const bool CAHitQuadrupletGenerator::fitFastCircle
private

Definition at line 142 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

const bool CAHitQuadrupletGenerator::fitFastCircleChi2Cut
private

Definition at line 143 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

const QuantityDependsPt CAHitQuadrupletGenerator::maxChi2
private

Definition at line 141 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

unsigned int CAHitQuadrupletGenerator::minLayers = 4
static

Definition at line 36 of file CAHitQuadrupletGenerator.h.

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

Definition at line 73 of file CAHitQuadrupletGenerator.h.

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

LayerCacheType CAHitQuadrupletGenerator::theLayerCache
private

Definition at line 71 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().

edm::EDGetTokenT<SeedingLayerSetsHits> CAHitQuadrupletGenerator::theSeedingLayerToken
private

Definition at line 69 of file CAHitQuadrupletGenerator.h.

Referenced by CAHitQuadrupletGenerator(), and hitQuadruplets().

const bool CAHitQuadrupletGenerator::useBendingCorrection
private

Definition at line 144 of file CAHitQuadrupletGenerator.h.

Referenced by hitQuadruplets().