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 29 of file CAHitTripletGenerator.h.

Member Typedef Documentation

Definition at line 31 of file CAHitTripletGenerator.h.

Definition at line 34 of file CAHitTripletGenerator.h.

Constructor & Destructor Documentation

CAHitTripletGenerator::CAHitTripletGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector &&  iC 
)
inline

Definition at line 37 of file CAHitTripletGenerator.h.

References looper::cfg.

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

Definition at line 36 of file CAHitTripletGenerator.cc.

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

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

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

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

Definition at line 43 of file CAHitTripletGenerator.h.

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

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

Definition at line 181 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().

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

Referenced by hitNtuplets().

const float CAHitTripletGenerator::caPhiCut = 1.f
private

Definition at line 127 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::caThetaCut = 0.00125f
private

Definition at line 126 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const float CAHitTripletGenerator::extraHitRPhitolerance
private

Definition at line 121 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

const QuantityDependsPt CAHitTripletGenerator::maxChi2
private

Definition at line 123 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().

unsigned int CAHitTripletGenerator::minLayers = 3
static

Definition at line 33 of file CAHitTripletGenerator.h.

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

Definition at line 55 of file CAHitTripletGenerator.h.

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

LayerCacheType CAHitTripletGenerator::theLayerCache
private

Definition at line 53 of file CAHitTripletGenerator.h.

const bool CAHitTripletGenerator::useBendingCorrection
private

Definition at line 124 of file CAHitTripletGenerator.h.

Referenced by hitNtuplets().