CMS 3D CMS Logo

List of all members | Classes | Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Private Attributes
CAHitTripletGenerator Class Reference

#include <CAHitTripletGenerator.h>

Classes

class  QuantityDependsPt
 
class  QuantityDependsPtEval
 

Public Types

typedef LayerHitMapCache LayerCacheType
 
typedef OrderedHitSeeds ResultType
 

Public Member Functions

 CAHitTripletGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
 
 CAHitTripletGenerator (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)
 
 ~CAHitTripletGenerator ()=default
 

Static Public Member Functions

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

Static Public Attributes

static unsigned int minLayers = 3
 

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
 
const bool useBendingCorrection
 

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 
)
inline

Definition at line 39 of file CAHitTripletGenerator.h.

References looper::cfg.

39 : CAHitTripletGenerator(cfg, iC) {}
CAHitTripletGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &&iC)
CAHitTripletGenerator::CAHitTripletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 35 of file CAHitTripletGenerator.cc.

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

36  :
37  extraHitRPhitolerance(cfg.getParameter<double>("extraHitRPhitolerance")), //extra window in ThirdHitPredictionFromCircle range (divide by R to get phi)
38  maxChi2(cfg.getParameter < edm::ParameterSet > ("maxChi2")),
39  useBendingCorrection(cfg.getParameter<bool>("useBendingCorrection")),
40  caThetaCut( cfg.getParameter<double>("CAThetaCut")),
41  caPhiCut(cfg.getParameter<double>("CAPhiCut")),
42  caHardPtCut(cfg.getParameter<double>("CAHardPtCut"))
43 {
44  edm::ParameterSet comparitorPSet = cfg.getParameter < edm::ParameterSet > ("SeedComparitorPSet");
45  std::string comparitorName = comparitorPSet.getParameter < std::string > ("ComponentName");
46  if (comparitorName != "none")
47  {
48  theComparitor.reset(
49  SeedComparitorFactory::get()->create(comparitorName,
50  comparitorPSet, iC));
51  }
52 }
T getParameter(std::string const &) const
const QuantityDependsPt maxChi2
def create(alignables, pedeDump, additionalData, outputFile, config)
std::unique_ptr< SeedComparitor > theComparitor
T get(const Candidate &c)
Definition: component.h:55
CAHitTripletGenerator::~CAHitTripletGenerator ( )
default

Member Function Documentation

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

Definition at line 54 of file CAHitTripletGenerator.cc.

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

54  {
55  desc.add<double>("extraHitRPhitolerance", 0.06);
56  desc.add<bool>("useBendingCorrection", false);
57  desc.add<double>("CAThetaCut", 0.00125);
58  desc.add<double>("CAPhiCut", 0.1);
59  desc.add<double>("CAHardPtCut", 0);
60 
61  edm::ParameterSetDescription descMaxChi2;
62  descMaxChi2.add<double>("pt1", 0.8);
63  descMaxChi2.add<double>("pt2", 2);
64  descMaxChi2.add<double>("value1", 50);
65  descMaxChi2.add<double>("value2", 8);
66  descMaxChi2.add<bool>("enabled", true);
67  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
68 
69  edm::ParameterSetDescription descComparitor;
70  descComparitor.add<std::string>("ComponentName", "none");
71  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
72  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
73 }
ParameterDescriptionBase * add(U const &iLabel, T const &value)
static const char* CAHitTripletGenerator::fillDescriptionsLabel ( )
inlinestatic

Definition at line 45 of file CAHitTripletGenerator.h.

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

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

Definition at line 180 of file CAHitTripletGenerator.cc.

References funct::abs(), caHardPtCut, caPhiCut, caThetaCut, vertices_cff::chi2, RZLine::chi2(), ThirdHitPredictionFromCircle::curvature(), PixelRecoUtilities::curvature(), GlobalErrorBase< T, ErrorWeightType >::czz(), relativeConstraints::error, CAHitTripletGenerator::QuantityDependsPt::evaluator(), extraHitRPhitolerance, CellularAutomaton::findTriplets(), g, CellularAutomaton::getAllCells(), mps_fire::i, triggerObjects_cff::id, diffTreeTool::index, PixelRecoUtilities::inversePt(), gedGsfElectrons_cfi::isBarrel, edm::isNotFinite(), geometryCSVtoXML::line, maxChi2, TrackingRegion::origin(), PixelSubdetector::PixelBarrel, point, EnergyCorrector::pt, GlobalErrorBase< T, ErrorWeightType >::rerr(), funct::sqr(), mathSSE::sqrt(), theComparitor, useBendingCorrection, CAHitTripletGenerator::QuantityDependsPtEval::value(), x, PV3DBase< T, PVType, FrameType >::x(), y, PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

183  {
184  CAGraph g;
185 
186  std::vector<const HitDoublets *> hitDoublets;
187 
188  std::vector<CACell::CAntuplet> foundTriplets;
189 
190  int index =0;
191  for(const auto& regionLayerPairs: regionDoublets) {
192 
193  const TrackingRegion& region = regionLayerPairs.region();
194  hitDoublets.clear();
195  foundTriplets.clear();
196 
197  if (index == 0){
198  createGraphStructure(layers, g);
199  }
200  else{
201  clearGraphStructure(layers, g);
202  }
203  fillGraph(layers, regionLayerPairs, g, hitDoublets);
204  CellularAutomaton ca(g);
205  ca.findTriplets(hitDoublets, foundTriplets, region, caThetaCut, caPhiCut,
206  caHardPtCut);
207 
208  auto & allCells = ca.getAllCells();
209 
210  const QuantityDependsPtEval maxChi2Eval = maxChi2.evaluator(es);
211 
212  // re-used thoughout, need to be vectors because of RZLine interface
213  std::array<float, 3> bc_r;
214  std::array<float, 3> bc_z;
215  std::array<float, 3> bc_errZ2;
216  std::array<GlobalPoint, 3> gps;
217  std::array<GlobalError, 3> ges;
218  std::array<bool, 3> barrels;
219 
220  unsigned int numberOfFoundTriplets = foundTriplets.size();
221  for (unsigned int tripletId = 0; tripletId < numberOfFoundTriplets;
222  ++tripletId)
223  {
224 
225  OrderedHitTriplet tmpTriplet(allCells[foundTriplets[tripletId][0]].getInnerHit(),
226  allCells[foundTriplets[tripletId][0]].getOuterHit(),
227  allCells[foundTriplets[tripletId][1]].getOuterHit());
228 
229  auto isBarrel = [](const unsigned id) -> bool
230  {
231  return id == PixelSubdetector::PixelBarrel;
232  };
233  for (unsigned int i = 0; i < 2; ++i)
234  {
235  auto const& ahit = allCells[foundTriplets[tripletId][i]].getInnerHit();
236  gps[i] = ahit->globalPosition();
237  ges[i] = ahit->globalPositionError();
238  barrels[i] = isBarrel(ahit->geographicalId().subdetId());
239  }
240 
241  auto const& ahit = allCells[foundTriplets[tripletId][1]].getOuterHit();
242  gps[2] = ahit->globalPosition();
243  ges[2] = ahit->globalPositionError();
244  barrels[2] = isBarrel(ahit->geographicalId().subdetId());
245 
246  PixelRecoLineRZ line(gps[0], gps[2]);
247  ThirdHitPredictionFromCircle predictionRPhi(gps[0], gps[2],
249  const float curvature = predictionRPhi.curvature(
250  ThirdHitPredictionFromCircle::Vector2D(gps[1].x(), gps[1].y()));
251  const float abscurv = std::abs(curvature);
252  const float thisMaxChi2 = maxChi2Eval.value(abscurv);
253  float chi2 = std::numeric_limits<float>::quiet_NaN();
254  // TODO: Do we have any use case to not use bending correction?
255 
257  {
258  // Following PixelFitterByConformalMappingAndLine
259  const float simpleCot = (gps.back().z() - gps.front().z())
260  / (gps.back().perp() - gps.front().perp());
261  const float pt = 1.f / PixelRecoUtilities::inversePt(abscurv, es);
262  for (int i = 0; i < 3; ++i)
263  {
264  const GlobalPoint & point = gps[i];
265  const GlobalError & error = ges[i];
266  bc_r[i] = sqrt(
267  sqr(point.x() - region.origin().x())
268  + sqr(point.y() - region.origin().y()));
270  es)(bc_r[i]);
271  bc_z[i] = point.z() - region.origin().z();
272  bc_errZ2[i] =
273  (barrels[i]) ?
274  error.czz() :
275  error.rerr(point) * sqr(simpleCot);
276  }
277  RZLine rzLine(bc_r, bc_z, bc_errZ2, RZLine::ErrZ2_tag());
278  chi2 = rzLine.chi2();
279  }
280  else
281  {
282  RZLine rzLine(gps, ges, barrels);
283  chi2 = rzLine.chi2();
284  }
285 
286  if (edm::isNotFinite(chi2) || chi2 > thisMaxChi2)
287  {
288  continue;
289 
290  }
291 
292  if (theComparitor)
293  {
294  if (!theComparitor->compatible(tmpTriplet))
295  {
296 
297  continue;
298  }
299  }
300  result[index].emplace_back(tmpTriplet);
301 
302  }
303  index++;
304  }
305 }
const QuantityDependsPt maxChi2
GlobalPoint const & origin() const
T y() const
Definition: PV3DBase.h:63
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)
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 
)

Member Data Documentation

const float CAHitTripletGenerator::caHardPtCut = 0.f
private

Definition at line 130 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::caPhiCut = 1.f
private

Definition at line 129 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::caThetaCut = 0.00125f
private

Definition at line 128 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::extraHitRPhitolerance
private

Definition at line 123 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const QuantityDependsPt CAHitTripletGenerator::maxChi2
private

Definition at line 125 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

unsigned int CAHitTripletGenerator::minLayers = 3
static

Definition at line 35 of file CAHitTripletGenerator.h.

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

Definition at line 57 of file CAHitTripletGenerator.h.

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

LayerCacheType CAHitTripletGenerator::theLayerCache
private

Definition at line 55 of file CAHitTripletGenerator.h.

const bool CAHitTripletGenerator::useBendingCorrection
private

Definition at line 126 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().