CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PixelTripletLowPtGenerator.cc
Go to the documentation of this file.
5 
8 
13 
16 
19 
20 #undef Debug
21 
22 using namespace std;
23 using namespace ctfseeding;
24 
25 /*****************************************************************************/
27  HitTripletGeneratorFromPairAndLayers(), // no theMaxElement used in this class
28  theTracker(nullptr),
29  theClusterShapeCacheToken(iC.consumes<SiPixelClusterShapeCache>(cfg.getParameter<edm::InputTag>("clusterShapeCacheSrc")))
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 }
38 
39 /*****************************************************************************/
41 
42 /*****************************************************************************/
43 
44 
45 /*****************************************************************************/
47  (const edm::EventSetup& es)
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 }
63 
64 /*****************************************************************************/
66  (const TrackingRecHit* recHit)
67 {
68  DetId detId = recHit->geographicalId();
69 
70  return
71  theTracker->idToDet(detId)->toGlobal(recHit->localPosition());
72 }
73 
74 /*****************************************************************************/
76  const TrackingRegion& region,
78  const edm::Event & ev,
79  const edm::EventSetup& es,
81  const std::vector<SeedingLayerSetsHits::SeedingLayer>& thirdLayers)
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, ev, 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 }
211 
212 
PixelTripletLowPtGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
T getParameter(std::string const &) const
std::vector< Hit > hits(float phiMin, float phiMax) const
int i
Definition: DBlmapReader.cc:9
tuple cfg
Definition: looper.py:293
T max() const
void getTracker(const edm::EventSetup &es)
virtual Location location() const =0
Which part of the detector (barrel, endcap)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
bool ev
virtual unsigned int size() const
T min() const
#define nullptr
std::unique_ptr< TripletFilter > theFilter
tuple result
Definition: query.py:137
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
SeedingHitSet::ConstRecHitPointer Hit
BaseTrackerRecHit const * Hit
virtual void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es, SeedingLayerSetsHits::SeedingLayerSet pairLayers, const std::vector< SeedingLayerSetsHits::SeedingLayer > &thirdLayers) override
Definition: DetId.h:18
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
GlobalPoint getGlobalPosition(const TrackingRecHit *recHit)
DetId geographicalId() const
bool isCompatibleWithMultipleScattering(GlobalPoint g3, const std::vector< const TrackingRecHit * > &h, std::vector< GlobalVector > &localDirs, const edm::EventSetup &es)
void getRanges(const DetLayer *layer, float phi[], float rz[])
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:23
virtual LocalPoint localPosition() const =0
tuple size
Write out results.