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, and fillDescriptions().

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

Definition at line 30 of file CAHitTripletGenerator.cc.

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

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

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

49  {
50  desc.add<double>("extraHitRPhitolerance", 0.06);
51  desc.add<bool>("useBendingCorrection", false);
52  desc.add<double>("CAThetaCut", 0.00125);
53  desc.add<double>("CAPhiCut", 0.1);
54  desc.add<double>("CAHardPtCut", 0);
55 
56  edm::ParameterSetDescription descMaxChi2;
57  descMaxChi2.add<double>("pt1", 0.8);
58  descMaxChi2.add<double>("pt2", 2);
59  descMaxChi2.add<double>("value1", 50);
60  descMaxChi2.add<double>("value2", 8);
61  descMaxChi2.add<bool>("enabled", true);
62  desc.add<edm::ParameterSetDescription>("maxChi2", descMaxChi2);
63 
64  edm::ParameterSetDescription descComparitor;
65  descComparitor.add<std::string>("ComponentName", "none");
66  descComparitor.setAllowAnything(); // until we have moved SeedComparitor too to EDProducers
67  desc.add<edm::ParameterSetDescription>("SeedComparitorPSet", descComparitor);
68 }
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 175 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, 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().

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