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

virtual 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)
 
virtual ~PixelTripletLowPtGenerator ()
 
- 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 26 of file PixelTripletLowPtGenerator.cc.

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

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

Definition at line 40 of file PixelTripletLowPtGenerator.cc.

References getTracker().

40 {}

Member Function Documentation

GlobalPoint PixelTripletLowPtGenerator::getGlobalPosition ( const TrackingRecHit recHit)
private

Definition at line 66 of file PixelTripletLowPtGenerator.cc.

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

Referenced by getTracker(), and hitTriplets().

67 {
68  DetId detId = recHit->geographicalId();
69 
70  return
71  theTracker->idToDet(detId)->toGlobal(recHit->localPosition());
72 }
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
virtual LocalPoint localPosition() const =0
Definition: DetId.h:18
const TrackerGeometry * theTracker
const TrackerGeomDet * idToDet(DetId) const override
DetId geographicalId() const
void PixelTripletLowPtGenerator::getTracker ( const edm::EventSetup es)
private

Definition at line 47 of file PixelTripletLowPtGenerator.cc.

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

Referenced by hitTriplets(), and ~PixelTripletLowPtGenerator().

48 {
49  if(theTracker == 0)
50  {
51  // Get tracker geometry
53  es.get<TrackerDigiGeometryRecord>().get(tracker);
54 
55  theTracker = tracker.product();
56  }
57 
58  if(!theFilter)
59  {
60  theFilter = std::make_unique<TripletFilter>(es);
61  }
62 }
std::unique_ptr< TripletFilter > theFilter
const T & get() const
Definition: EventSetup.h:55
const TrackerGeometry * theTracker
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 75 of file PixelTripletLowPtGenerator.cc.

References builderName, MessageLogger_cfi::cerr, checkClusterShape, checkMultipleScattering, 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, hiPixelPairStep_cff::points, edm::ESHandle< T >::product(), rzTolerance, OrderedHitPairs::size(), findQualityFiles::size, DetLayer::subDetector(), theClusterShapeCacheToken, theFilter, HitTripletGeneratorFromPairAndLayers::theLayerCache, and HitTripletGeneratorFromPairAndLayers::thePairGenerator.

82 {
83 
84  //Retrieve tracker topology from geometry
86  es.get<TrackerTopologyRcd>().get(tTopoHand);
87  const TrackerTopology *tTopo=tTopoHand.product();
88 
89  edm::Handle<SiPixelClusterShapeCache> clusterShapeCache;
90  ev.getByToken(theClusterShapeCacheToken, clusterShapeCache);
91 
92  // Generate pairs
93  OrderedHitPairs pairs; pairs.reserve(30000);
94  thePairGenerator->hitPairs(region,pairs,ev,es, pairLayers);
95 
96  if (pairs.size() == 0) return;
97 
98  int size = thirdLayers.size();
99 
100  // Set aliases
101  const RecHitsSortedInPhi **thirdHitMap = new const RecHitsSortedInPhi*[size];
102  for(int il=0; il<size; il++)
103  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, es);
104 
105  // Get tracker
106  getTracker(es);
107 
108  // Look at all generated pairs
109  for(OrderedHitPairs::const_iterator ip = pairs.begin();
110  ip!= pairs.end(); ip++)
111  {
112  // Fill rechits and points
113  vector<const TrackingRecHit*> recHits(3);
114  vector<GlobalPoint> points(3);
115 
116  recHits[0] = (*ip).inner()->hit();
117  recHits[1] = (*ip).outer()->hit();
118 
119 #ifdef Debug
120  cerr << " RecHits " + HitInfo::getInfo(*recHits[0]) +
121  HitInfo::getInfo(*recHits[1]) << endl;
122 #endif
123 
124  for(int i=0; i<2; i++)
125  points[i] = getGlobalPosition(recHits[i]);
126 
127  // Initialize helix prediction
129  thePrediction(region,
130  points[0],points[1], es,
132 
133  // Look at all layers
134  for(int il=0; il<size; il++)
135  {
136  const DetLayer * layer = thirdLayers[il].detLayer();
137 
138 #ifdef Debug
139  cerr << " check layer " << layer->subDetector()
140  << " " << layer->location() << endl;
141 #endif
142 
143  // Get ranges for the third hit
144  float phi[2],rz[2];
145  thePrediction.getRanges(layer, phi,rz);
146 
147  PixelRecoRange<float> phiRange(phi[0] , phi[1] );
148  PixelRecoRange<float> rzRange( rz[0] - rzTolerance, rz[1] + rzTolerance);
149 
150  // Get third hit candidates from cache
152  vector<Hit> thirdHits = thirdHitMap[il]->hits(phiRange.min(),phiRange.max());
153  typedef vector<Hit>::const_iterator IH;
154 
155  for (IH th=thirdHits.begin(), eh=thirdHits.end(); th < eh; ++th)
156  {
157  // Fill rechit and point
158  recHits[2] = (*th)->hit();
159  points[2] = getGlobalPosition(recHits[2]);
160 
161 #ifdef Debug
162  cerr << " third hit " + HitInfo::getInfo(*recHits[2]) << endl;
163 #endif
164 
165  // Check if third hit is compatible with multiple scattering
166  vector<GlobalVector> globalDirs;
167  if(thePrediction.isCompatibleWithMultipleScattering
168  (points[2], recHits, globalDirs, es) == false)
169  {
170 #ifdef Debug
171  cerr << " not compatible: multiple scattering" << endl;
172 #endif
173  if(checkMultipleScattering) continue;
174  }
175 
176  // Convert to localDirs
177 /*
178  vector<LocalVector> localDirs;
179  vector<GlobalVector>::const_iterator globalDir = globalDirs.begin();
180  for(vector<const TrackingRecHit *>::const_iterator
181  recHit = recHits.begin();
182  recHit != recHits.end(); recHit++)
183  {
184  localDirs.push_back(theTracker->idToDet(
185  (*recHit)->geographicalId())->toLocal(*globalDir));
186  globalDir++;
187  }
188 */
189 
190  // Check if the cluster shapes are compatible with thrusts
192  {
193  if(! theFilter->checkTrack(recHits,globalDirs,tTopo, *clusterShapeCache))
194  {
195 #ifdef Debug
196  cerr << " not compatible: cluster shape" << endl;
197 #endif
198  continue;
199  }
200  }
201 
202  // All checks passed, put triplet back
203  result.push_back(OrderedHitTriplet((*ip).inner(),(*ip).outer(),*th));
204  }
205  }
206  }
207  delete [] thirdHitMap;
208 
209  return;
210 }
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:460
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)
virtual unsigned int size() const
std::unique_ptr< TripletFilter > theFilter
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
const T & get() const
Definition: EventSetup.h:55
GlobalPoint getGlobalPosition(const TrackingRecHit *recHit)
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:23
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 211 of file PixelTripletLowPtGenerator.cc.

References Exception.

219 {
220  throw cms::Exception("Error")<<"PixelTripletLowPtGenerator::hitTriplets is not implemented \n";
221 }

Member Data Documentation

std::string PixelTripletLowPtGenerator::builderName
private

Definition at line 58 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

bool PixelTripletLowPtGenerator::checkClusterShape
private

Definition at line 60 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

bool PixelTripletLowPtGenerator::checkMultipleScattering
private

Definition at line 59 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::maxAngleRatio
private

Definition at line 56 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::nSigMultipleScattering
private

Definition at line 54 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::rzTolerance
private

Definition at line 55 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

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

Definition at line 53 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets().

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

Definition at line 50 of file PixelTripletLowPtGenerator.h.

Referenced by getTracker(), and hitTriplets().

const TrackerGeometry* PixelTripletLowPtGenerator::theTracker
private

Definition at line 49 of file PixelTripletLowPtGenerator.h.

Referenced by getGlobalPosition(), and getTracker().