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
CAHitTripletGenerator Class Reference

#include <CAHitTripletGenerator.h>

Inheritance diagram for CAHitTripletGenerator:
HitTripletGenerator OrderedHitsGenerator

Classes

class  QuantityDependsPt
 
class  QuantityDependsPtEval
 

Public Types

typedef LayerHitMapCache LayerCacheType
 
typedef OrderedHitTriplets ResultType
 

Public Member Functions

 CAHitTripletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC, bool needSeedingLayerSetsHits=true)
 
 CAHitTripletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC, bool needSeedingLayerSetsHits=true)
 
void hitNtuplets (const IntermediateHitDoublets::RegionLayerSets &regionLayerPairs, OrderedHitTriplets &result, const edm::EventSetup &es, const SeedingLayerSetsHits &layers)
 
virtual void hitTriplets (const TrackingRegion &reg, OrderedHitTriplets &triplets, const edm::Event &ev, const edm::EventSetup &es)
 from base class More...
 
void initEvent (const edm::Event &ev, const edm::EventSetup &es)
 
virtual ~CAHitTripletGenerator ()
 
- Public Member Functions inherited from HitTripletGenerator
virtual void clear () final
 
 HitTripletGenerator (unsigned int size=500)
 
 HitTripletGenerator (HitTripletGenerator const &other)
 
virtual void hitTriplets (const TrackingRegion &reg, OrderedHitTriplets &prs, const edm::EventSetup &es)
 
virtual const OrderedHitTripletsrun (const TrackingRegion &region, const edm::Event &ev, const edm::EventSetup &es) final
 
virtual ~HitTripletGenerator ()
 
- 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 = 3
 

Private Member Functions

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

Private Attributes

const float caHardPtCut = 0.f
 
const float caPhiCut = 1.f
 
const float caThetaCut = 0.00125f
 
const float extraHitRPhitolerance
 
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 31 of file CAHitTripletGenerator.h.

Member Typedef Documentation

Definition at line 33 of file CAHitTripletGenerator.h.

Definition at line 36 of file CAHitTripletGenerator.h.

Constructor & Destructor Documentation

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

Definition at line 40 of file CAHitTripletGenerator.h.

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

Definition at line 36 of file CAHitTripletGenerator.cc.

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

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  useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
42  caThetaCut( cfg.getParameter<double>("CAThetaCut")),
43  caPhiCut(cfg.getParameter<double>("CAPhiCut")),
44  caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
45 {
46 
47  if(needSeedingLayerSetsHits)
49 
50  if (cfg.exists("SeedComparitorPSet"))
51  {
52  edm::ParameterSet comparitorPSet = cfg.getParameter < edm::ParameterSet
53  > ("SeedComparitorPSet");
54  std::string comparitorName = comparitorPSet.getParameter < std::string
55  > ("ComponentName");
56  if (comparitorName != "none")
57  {
58  theComparitor.reset(
59  SeedComparitorFactory::get()->create(comparitorName,
60  comparitorPSet, iC));
61  }
62  }
63 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
const QuantityDependsPt maxChi2
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< SeedingLayerSetsHits > theSeedingLayerToken
std::unique_ptr< SeedComparitor > theComparitor
T get(const Candidate &c)
Definition: component.h:55
CAHitTripletGenerator::~CAHitTripletGenerator ( )
virtual

Definition at line 65 of file CAHitTripletGenerator.cc.

66 {
67 }

Member Function Documentation

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

Definition at line 69 of file CAHitTripletGenerator.cc.

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

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

Definition at line 46 of file CAHitTripletGenerator.h.

46 { return "caHitTriplet"; }
void CAHitTripletGenerator::hitNtuplets ( const IntermediateHitDoublets::RegionLayerSets regionLayerPairs,
OrderedHitTriplets result,
const edm::EventSetup es,
const SeedingLayerSetsHits layers 
)

Definition at line 196 of file CAHitTripletGenerator.cc.

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

199  {
200  CAGraph g;
201 
202  std::vector<const HitDoublets *> hitDoublets;
203 
204  auto layerPairEqual = [](const IntermediateHitDoublets::LayerPairHitDoublets& pair,
207  return pair.innerLayerIndex() == inner && pair.outerLayerIndex() == outer;
208  };
209  fillGraph(layers, g, hitDoublets, [&](const SeedingLayerSetsHits::SeedingLayer& inner, const SeedingLayerSetsHits::SeedingLayer& outer, std::vector<const HitDoublets *>& hitDoublets) {
210  using namespace std::placeholders;
211  auto found = std::find_if(regionLayerPairs.begin(), regionLayerPairs.end(), std::bind(layerPairEqual, _1, inner.index(), outer.index()));
212  if(found != regionLayerPairs.end()) {
213  hitDoublets.emplace_back(&(found->doublets()));
214  return true;
215  }
216  return false;
217  });
218 
219  hitTriplets(regionLayerPairs.region(), result, hitDoublets, g, es);
220 }
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
SeedingLayerSetsHits::LayerIndex innerLayerIndex() const
const TrackingRegion & region() const
virtual void hitTriplets(const TrackingRegion &reg, OrderedHitTriplets &triplets, const edm::Event &ev, const edm::EventSetup &es)
from base class
SeedingLayerSetsHits::LayerIndex outerLayerIndex() const
void CAHitTripletGenerator::hitTriplets ( const TrackingRegion reg,
OrderedHitTriplets triplets,
const edm::Event ev,
const edm::EventSetup es 
)
virtual

from base class

Implements HitTripletGenerator.

Definition at line 159 of file CAHitTripletGenerator.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().

162 {
164  ev.getByToken(theSeedingLayerToken, hlayers);
165  const SeedingLayerSetsHits& layers = *hlayers;
166  if (layers.numberOfLayersInSet() != 3)
167  throw cms::Exception("Configuration")
168  << "CAHitTripletGenerator expects SeedingLayerSetsHits::numberOfLayersInSet() to be 3, got "
169  << layers.numberOfLayersInSet();
170 
171  CAGraph g;
172 
173  std::vector<HitDoublets> hitDoublets;
174 
175 
176 
177  HitPairGeneratorFromLayerPair thePairGenerator(0, 1, &theLayerCache);
178  fillGraph(layers, g, hitDoublets, [&](const SeedingLayerSetsHits::SeedingLayer& inner, const SeedingLayerSetsHits::SeedingLayer& outer, std::vector<HitDoublets>& hitDoublets) {
179  hitDoublets.emplace_back(thePairGenerator.doublets(region, ev, es, inner, outer));
180  return true;
181  });
182 
183 
184  if (theComparitor)
185  theComparitor->init(ev, es);
186 
187  std::vector<const HitDoublets *> hitDoubletsPtr;
188  hitDoubletsPtr.reserve(hitDoublets.size());
189  for(const auto& e: hitDoublets)
190  hitDoubletsPtr.emplace_back(&e);
191 
192  hitTriplets(region, result, hitDoubletsPtr, g, es);
194 }
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
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
edm::EDGetTokenT< SeedingLayerSetsHits > theSeedingLayerToken
virtual void hitTriplets(const TrackingRegion &reg, OrderedHitTriplets &triplets, const edm::Event &ev, const edm::EventSetup &es)
from base class
std::unique_ptr< SeedComparitor > theComparitor
void CAHitTripletGenerator::hitTriplets ( const TrackingRegion reg,
OrderedHitTriplets result,
std::vector< const HitDoublets * > &  hitDoublets,
const CAGraph g,
const edm::EventSetup es 
)
private

Definition at line 223 of file CAHitTripletGenerator.cc.

References funct::abs(), caHardPtCut, caPhiCut, caThetaCut, RZLine::chi2(), beam_dqm_sourceclient-live_cfg::chi2, LayerHitMapCache::clear(), ThirdHitPredictionFromCircle::curvature(), PixelRecoUtilities::curvature(), GlobalErrorBase< T, ErrorWeightType >::czz(), relativeConstraints::error, CAHitTripletGenerator::QuantityDependsPt::evaluator(), extraHitRPhitolerance, CellularAutomaton::findTriplets(), 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, theLayerCache, useBendingCorrection, CAHitTripletGenerator::QuantityDependsPtEval::value(), x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

226  {
227  std::vector<CACell::CAntuplet> foundTriplets;
228  CellularAutomaton ca(g);
229 
230  ca.findTriplets(hitDoublets, foundTriplets, region, caThetaCut, caPhiCut,
231  caHardPtCut);
232 
233  unsigned int numberOfFoundTriplets = foundTriplets.size();
234 
235  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
236 
237  // re-used thoughout, need to be vectors because of RZLine interface
238  std::array<float, 3> bc_r;
239  std::array<float, 3> bc_z;
240  std::array<float, 3> bc_errZ2;
241  std::array < GlobalPoint, 3 > gps;
242  std::array < GlobalError, 3 > ges;
243  std::array<bool, 3> barrels;
244 
245  for (unsigned int tripletId = 0; tripletId < numberOfFoundTriplets;
246  ++tripletId)
247  {
248 
249  OrderedHitTriplet tmpTriplet(foundTriplets[tripletId][0]->getInnerHit(),
250  foundTriplets[tripletId][0]->getOuterHit(),
251  foundTriplets[tripletId][1]->getOuterHit());
252 
253  auto isBarrel = [](const unsigned id) -> bool
254  {
255  return id == PixelSubdetector::PixelBarrel;
256  };
257  for (unsigned int i = 0; i < 2; ++i)
258  {
259  auto const& ahit = foundTriplets[tripletId][i]->getInnerHit();
260  gps[i] = ahit->globalPosition();
261  ges[i] = ahit->globalPositionError();
262  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
263  }
264 
265  auto const& ahit = foundTriplets[tripletId][1]->getOuterHit();
266  gps[2] = ahit->globalPosition();
267  ges[2] = ahit->globalPositionError();
268  barrels[2] = isBarrel(ahit->geographicalId().subdetId());
269 
270  PixelRecoLineRZ line(gps[0], gps[2]);
271  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2],
273  const float curvature = predictionRPhi.curvature(
274  ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
275  const float abscurv = std::abs(curvature);
276  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
277  float chi2 = std::numeric_limits<float>::quiet_NaN();
278  // TODO: Do we have any use case to not use bending correction?
280  {
281  // Following PixelFitterByConformalMappingAndLine
282  const float simpleCot = (gps.back().z() - gps.front().z())
283  / (gps.back().perp() - gps.front().perp());
284  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
285  for (int i = 0; i < 3; ++i)
286  {
287  const GlobalPoint & point = gps[i];
288  const GlobalError & error = ges[i];
289  bc_r[i] = sqrt(
290  sqr(point.x() - region.origin().x())
291  + sqr(point.y() - region.origin().y()));
293  es)(bc_r[i]);
294  bc_z[i] = point.z() - region.origin().z();
295  bc_errZ2[i] =
296  (barrels[i]) ?
297  error.czz() :
298  error.rerr(point) * sqr(simpleCot);
299  }
300  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
301  chi2 = rzLine.chi2();
302  }
303  else
304  {
305  RZLine rzLine(gps, ges, barrels);
306  chi2 = rzLine.chi2();
307  }
308 
309  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
310  {
311  continue;
312 
313  }
314 
315  if (theComparitor)
316  {
317  if (!theComparitor->compatible(tmpTriplet))
318  {
319 
320  continue;
321  }
322  }
323 
324  result.push_back(tmpTriplet);
325 
326  }
328 }
int i
Definition: DBlmapReader.cc:9
const QuantityDependsPt maxChi2
bool isBarrel(GeomDetEnumerators::SubDetector m)
T y() const
Definition: PV3DBase.h:63
T inversePt(T curvature, const edm::EventSetup &iSetup)
QuantityDependsPtEval evaluator(const edm::EventSetup &es) const
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
T rerr(const GlobalPoint &aPoint) const
Square< F >::type sqr(const F &f)
Definition: Square.h:13
std::unique_ptr< SeedComparitor > theComparitor
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 CAHitTripletGenerator::initEvent ( const edm::Event ev,
const edm::EventSetup es 
)

Definition at line 90 of file CAHitTripletGenerator.cc.

References theComparitor.

90  {
91  if (theComparitor) theComparitor->init(ev, es);
92 }
std::unique_ptr< SeedComparitor > theComparitor

Member Data Documentation

const float CAHitTripletGenerator::caHardPtCut = 0.f
private

Definition at line 144 of file CAHitTripletGenerator.h.

Referenced by hitTriplets().

const float CAHitTripletGenerator::caPhiCut = 1.f
private

Definition at line 143 of file CAHitTripletGenerator.h.

Referenced by hitTriplets().

const float CAHitTripletGenerator::caThetaCut = 0.00125f
private

Definition at line 142 of file CAHitTripletGenerator.h.

Referenced by hitTriplets().

const float CAHitTripletGenerator::extraHitRPhitolerance
private

Definition at line 137 of file CAHitTripletGenerator.h.

Referenced by hitTriplets().

const QuantityDependsPt CAHitTripletGenerator::maxChi2
private

Definition at line 139 of file CAHitTripletGenerator.h.

Referenced by hitTriplets().

unsigned int CAHitTripletGenerator::minLayers = 3
static

Definition at line 35 of file CAHitTripletGenerator.h.

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

Definition at line 71 of file CAHitTripletGenerator.h.

Referenced by CAHitTripletGenerator(), hitTriplets(), and initEvent().

LayerCacheType CAHitTripletGenerator::theLayerCache
private

Definition at line 69 of file CAHitTripletGenerator.h.

Referenced by hitTriplets().

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

Definition at line 67 of file CAHitTripletGenerator.h.

Referenced by CAHitTripletGenerator(), and hitTriplets().

const bool CAHitTripletGenerator::useBendingCorrection
private

Definition at line 140 of file CAHitTripletGenerator.h.

Referenced by hitTriplets().