CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
PixelTripletLowPtGenerator Class Reference

#include <PixelTripletLowPtGenerator.h>

Inheritance diagram for PixelTripletLowPtGenerator:
HitTripletGeneratorFromPairAndLayers

Public Member Functions

void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, const SeedingLayerSetsHits::SeedingLayerSet &pairLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers) override
 
void hitTriplets (const TrackingRegion &region, OrderedHitTriplets &result, const edm::EventSetup &es, const HitDoublets &doublets, const RecHitsSortedInPhi **thirdHitMap, const std::vector< const DetLayer * > &thirdLayerDetLayer, const int nThirdLayers) override
 
 PixelTripletLowPtGenerator (const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
 
 ~PixelTripletLowPtGenerator () override
 
- Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
 HitTripletGeneratorFromPairAndLayers (unsigned int maxElement=0)
 
 HitTripletGeneratorFromPairAndLayers (const edm::ParameterSet &pset)
 
void init (std::unique_ptr< HitPairGeneratorFromLayerPair > &&pairs, LayerCacheType *layerCache)
 
const HitPairGeneratorFromLayerPairpairGenerator () const
 
virtual ~HitTripletGeneratorFromPairAndLayers ()
 

Private Member Functions

GlobalPoint getGlobalPosition (const TrackingRecHit *recHit)
 
void getTracker (const edm::EventSetup &es)
 

Private Attributes

std::string builderName
 
bool checkClusterShape
 
bool checkMultipleScattering
 
double maxAngleRatio
 
double nSigMultipleScattering
 
double rzTolerance
 
edm::EDGetTokenT< SiPixelClusterShapeCachetheClusterShapeCacheToken
 
std::unique_ptr< TripletFiltertheFilter
 
const TrackerGeometrytheTracker
 

Additional Inherited Members

- Public Types inherited from HitTripletGeneratorFromPairAndLayers
typedef LayerHitMapCache LayerCacheType
 
- Static Public Member Functions inherited from HitTripletGeneratorFromPairAndLayers
static void fillDescriptions (edm::ParameterSetDescription &desc)
 
- Protected Attributes inherited from HitTripletGeneratorFromPairAndLayers
LayerCacheTypetheLayerCache
 
const unsigned int theMaxElement
 
std::unique_ptr< HitPairGeneratorFromLayerPairthePairGenerator
 

Detailed Description

Definition at line 23 of file PixelTripletLowPtGenerator.h.

Constructor & Destructor Documentation

PixelTripletLowPtGenerator::PixelTripletLowPtGenerator ( const edm::ParameterSet cfg,
edm::ConsumesCollector iC 
)

Definition at line 25 of file PixelTripletLowPtGenerator.cc.

References builderName, checkClusterShape, checkMultipleScattering, edm::ParameterSet::getParameter(), maxAngleRatio, nSigMultipleScattering, and rzTolerance.

26  : HitTripletGeneratorFromPairAndLayers(), // no theMaxElement used in this class
27  theTracker(nullptr),
29  iC.consumes<SiPixelClusterShapeCache>(cfg.getParameter<edm::InputTag>("clusterShapeCacheSrc"))) {
30  checkMultipleScattering = cfg.getParameter<bool>("checkMultipleScattering");
31  nSigMultipleScattering = cfg.getParameter<double>("nSigMultipleScattering");
32  checkClusterShape = cfg.getParameter<bool>("checkClusterShape");
33  rzTolerance = cfg.getParameter<double>("rzTolerance");
34  maxAngleRatio = cfg.getParameter<double>("maxAngleRatio");
35  builderName = cfg.getParameter<string>("TTRHBuilder");
36 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
T getParameter(std::string const &) const
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
const TrackerGeometry * theTracker
PixelTripletLowPtGenerator::~PixelTripletLowPtGenerator ( )
override

Definition at line 39 of file PixelTripletLowPtGenerator.cc.

39 {}

Member Function Documentation

GlobalPoint PixelTripletLowPtGenerator::getGlobalPosition ( const TrackingRecHit recHit)
private

Definition at line 59 of file PixelTripletLowPtGenerator.cc.

References TrackingRecHit::geographicalId(), TrackerGeometry::idToDet(), TrackingRecHit::localPosition(), theTracker, and GeomDet::toGlobal().

Referenced by hitTriplets().

59  {
60  DetId detId = recHit->geographicalId();
61 
62  return theTracker->idToDet(detId)->toGlobal(recHit->localPosition());
63 }
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
virtual LocalPoint localPosition() const =0
Definition: DetId.h:17
const TrackerGeometry * theTracker
const TrackerGeomDet * idToDet(DetId) const override
DetId geographicalId() const
void PixelTripletLowPtGenerator::getTracker ( const edm::EventSetup es)
private

Definition at line 44 of file PixelTripletLowPtGenerator.cc.

References edm::EventSetup::get(), edm::ESHandle< T >::product(), theFilter, theTracker, and PbPb_ZMuSkimMuonDPG_cff::tracker.

Referenced by hitTriplets().

44  {
45  if (theTracker == nullptr) {
46  // Get tracker geometry
48  es.get<TrackerDigiGeometryRecord>().get(tracker);
49 
50  theTracker = tracker.product();
51  }
52 
53  if (!theFilter) {
54  theFilter = std::make_unique<TripletFilter>(es);
55  }
56 }
std::unique_ptr< TripletFilter > theFilter
const TrackerGeometry * theTracker
T get() const
Definition: EventSetup.h:73
T const * product() const
Definition: ESHandle.h:86
void PixelTripletLowPtGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets trs,
const edm::Event ev,
const edm::EventSetup es,
const SeedingLayerSetsHits::SeedingLayerSet pairLayers,
const std::vector< SeedingLayerSetsHits::SeedingLayer > &  thirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 66 of file PixelTripletLowPtGenerator.cc.

References builderName, beam_dqm_sourceclient-live_cfg::cerr, checkClusterShape, checkMultipleScattering, SiPixelPhase1OnlineDQM_cff::clusterShapeCache, edm::EventSetup::get(), edm::Event::getByToken(), getGlobalPosition(), HitInfo::getInfo(), ThirdHitPrediction::getRanges(), getTracker(), RecHitsSortedInPhi::hits(), mps_fire::i, ThirdHitPrediction::isCompatibleWithMultipleScattering(), DetLayer::location(), PixelRecoRange< T >::max(), maxAngleRatio, PixelRecoRange< T >::min(), nSigMultipleScattering, phi, HLT_2018_cff::points, edm::ESHandle< T >::product(), FastTrackerRecHitMaskProducer_cfi::recHits, HLT_2018_cff::region, rzTolerance, findQualityFiles::size, DetLayer::subDetector(), theClusterShapeCacheToken, theFilter, HitTripletGeneratorFromPairAndLayers::theLayerCache, and HitTripletGeneratorFromPairAndLayers::thePairGenerator.

71  {
72  //Retrieve tracker topology from geometry
74  es.get<TrackerTopologyRcd>().get(tTopoHand);
75  const TrackerTopology* tTopo = tTopoHand.product();
76 
78  ev.getByToken(theClusterShapeCacheToken, clusterShapeCache);
79 
80  // Generate pairs
81  OrderedHitPairs pairs;
82  pairs.reserve(30000);
83  thePairGenerator->hitPairs(region, pairs, ev, es, pairLayers);
84 
85  if (pairs.empty())
86  return;
87 
88  int size = thirdLayers.size();
89 
90  // Set aliases
91  const RecHitsSortedInPhi** thirdHitMap = new const RecHitsSortedInPhi*[size];
92  for (int il = 0; il < size; il++)
93  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, es);
94 
95  // Get tracker
96  getTracker(es);
97 
98  // Look at all generated pairs
99  for (OrderedHitPairs::const_iterator ip = pairs.begin(); ip != pairs.end(); ip++) {
100  // Fill rechits and points
101  vector<const TrackingRecHit*> recHits(3);
102  vector<GlobalPoint> points(3);
103 
104  recHits[0] = (*ip).inner()->hit();
105  recHits[1] = (*ip).outer()->hit();
106 
107 #ifdef Debug
108  cerr << " RecHits " + HitInfo::getInfo(*recHits[0]) + HitInfo::getInfo(*recHits[1]) << endl;
109 #endif
110 
111  for (int i = 0; i < 2; i++)
113 
114  // Initialize helix prediction
115  ThirdHitPrediction thePrediction(
117 
118  // Look at all layers
119  for (int il = 0; il < size; il++) {
120  const DetLayer* layer = thirdLayers[il].detLayer();
121 
122 #ifdef Debug
123  cerr << " check layer " << layer->subDetector() << " " << layer->location() << endl;
124 #endif
125 
126  // Get ranges for the third hit
127  float phi[2], rz[2];
128  thePrediction.getRanges(layer, phi, rz);
129 
130  PixelRecoRange<float> phiRange(phi[0], phi[1]);
131  PixelRecoRange<float> rzRange(rz[0] - rzTolerance, rz[1] + rzTolerance);
132 
133  // Get third hit candidates from cache
135  vector<Hit> thirdHits = thirdHitMap[il]->hits(phiRange.min(), phiRange.max());
136  typedef vector<Hit>::const_iterator IH;
137 
138  for (IH th = thirdHits.begin(), eh = thirdHits.end(); th < eh; ++th) {
139  // Fill rechit and point
140  recHits[2] = (*th)->hit();
142 
143 #ifdef Debug
144  cerr << " third hit " + HitInfo::getInfo(*recHits[2]) << endl;
145 #endif
146 
147  // Check if third hit is compatible with multiple scattering
148  vector<GlobalVector> globalDirs;
149  if (thePrediction.isCompatibleWithMultipleScattering(points[2], recHits, globalDirs, es) == false) {
150 #ifdef Debug
151  cerr << " not compatible: multiple scattering" << endl;
152 #endif
154  continue;
155  }
156 
157  // Convert to localDirs
158  /*
159  vector<LocalVector> localDirs;
160  vector<GlobalVector>::const_iterator globalDir = globalDirs.begin();
161  for(vector<const TrackingRecHit *>::const_iterator
162  recHit = recHits.begin();
163  recHit != recHits.end(); recHit++)
164  {
165  localDirs.push_back(theTracker->idToDet(
166  (*recHit)->geographicalId())->toLocal(*globalDir));
167  globalDir++;
168  }
169 */
170 
171  // Check if the cluster shapes are compatible with thrusts
172  if (checkClusterShape) {
173  if (!theFilter->checkTrack(recHits, globalDirs, tTopo, *clusterShapeCache)) {
174 #ifdef Debug
175  cerr << " not compatible: cluster shape" << endl;
176 #endif
177  continue;
178  }
179  }
180 
181  // All checks passed, put triplet back
182  result.push_back(OrderedHitTriplet((*ip).inner(), (*ip).outer(), *th));
183  }
184  }
185  }
186  delete[] thirdHitMap;
187 
188  return;
189 }
size
Write out results.
std::vector< Hit > hits(float phiMin, float phiMax) const
void getTracker(const edm::EventSetup &es)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
virtual Location location() const =0
Which part of the detector (barrel, endcap)
std::unique_ptr< TripletFilter > theFilter
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
T get() const
Definition: EventSetup.h:73
GlobalPoint getGlobalPosition(const TrackingRecHit *recHit)
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:19
T const * product() const
Definition: ESHandle.h:86
void PixelTripletLowPtGenerator::hitTriplets ( const TrackingRegion region,
OrderedHitTriplets result,
const edm::EventSetup es,
const HitDoublets doublets,
const RecHitsSortedInPhi **  thirdHitMap,
const std::vector< const DetLayer * > &  thirdLayerDetLayer,
const int  nThirdLayers 
)
overridevirtual

Implements HitTripletGeneratorFromPairAndLayers.

Definition at line 190 of file PixelTripletLowPtGenerator.cc.

References Exception.

196  {
197  throw cms::Exception("Error") << "PixelTripletLowPtGenerator::hitTriplets is not implemented \n";
198 }

Member Data Documentation

std::string PixelTripletLowPtGenerator::builderName
private

Definition at line 55 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

bool PixelTripletLowPtGenerator::checkClusterShape
private

Definition at line 57 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

bool PixelTripletLowPtGenerator::checkMultipleScattering
private

Definition at line 56 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::maxAngleRatio
private

Definition at line 53 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::nSigMultipleScattering
private

Definition at line 51 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::rzTolerance
private

Definition at line 52 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

edm::EDGetTokenT<SiPixelClusterShapeCache> PixelTripletLowPtGenerator::theClusterShapeCacheToken
private

Definition at line 50 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets().

std::unique_ptr<TripletFilter> PixelTripletLowPtGenerator::theFilter
private

Definition at line 48 of file PixelTripletLowPtGenerator.h.

Referenced by getTracker(), and hitTriplets().

const TrackerGeometry* PixelTripletLowPtGenerator::theTracker
private

Definition at line 47 of file PixelTripletLowPtGenerator.h.

Referenced by getGlobalPosition(), and getTracker().