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 
7 
11 
14 
15 #undef Debug
16 
17 using namespace std;
18 using namespace ctfseeding;
19 
20 /*****************************************************************************/
22  LayerCacheType* layerCache)
23 {
24  thePairGenerator = pairs.clone();
25  theLayerCache = layerCache;
26 
27  checkMultipleScattering = ps.getParameter<bool>("checkMultipleScattering");
28  nSigMultipleScattering = ps.getParameter<double>("nSigMultipleScattering");
29  checkClusterShape = ps.getParameter<bool>("checkClusterShape");
30  rzTolerance = ps.getParameter<double>("rzTolerance");
31  maxAngleRatio = ps.getParameter<double>("maxAngleRatio");
32  builderName = ps.getParameter<string>("TTRHBuilder");
33 }
34 
35 /*****************************************************************************/
37  std::vector<SeedingLayerSetsHits::SeedingLayer> thirdLayers) {
38  thePairGenerator->setSeedingLayers(pairLayers);
39  theLayers = thirdLayers;
40 }
41 
42 /*****************************************************************************/
44  (const edm::EventSetup& es)
45 {
46  if(theTracker == 0)
47  {
48  // Get tracker geometry
50  es.get<TrackerDigiGeometryRecord>().get(tracker);
51 
52  theTracker = tracker.product();
53  }
54 
55  if(theFilter == 0)
56  {
57  theFilter = new TripletFilter(es);
58  }
59 }
60 
61 /*****************************************************************************/
63  (const TrackingRecHit* recHit)
64 {
65  DetId detId = recHit->geographicalId();
66 
67  return
68  theTracker->idToDet(detId)->toGlobal(recHit->localPosition());
69 }
70 
71 /*****************************************************************************/
73  const TrackingRegion& region,
75  const edm::Event & ev,
76  const edm::EventSetup& es)
77 {
78 
79  //Retrieve tracker topology from geometry
81  es.get<IdealGeometryRecord>().get(tTopoHand);
82  const TrackerTopology *tTopo=tTopoHand.product();
83 
84  // Generate pairs
85  OrderedHitPairs pairs; pairs.reserve(30000);
86  thePairGenerator->hitPairs(region,pairs,ev,es);
87 
88  if (pairs.size() == 0) return;
89 
90  int size = theLayers.size();
91 
92  // Set aliases
93  const RecHitsSortedInPhi **thirdHitMap = new const RecHitsSortedInPhi*[size];
94  for(int il=0; il<size; il++)
95  thirdHitMap[il] = &(*theLayerCache)(theLayers[il], region, ev, es);
96 
97  // Get tracker
98  getTracker(es);
99 
100  // Look at all generated pairs
101  for(OrderedHitPairs::const_iterator ip = pairs.begin();
102  ip!= pairs.end(); ip++)
103  {
104  // Fill rechits and points
105  vector<const TrackingRecHit*> recHits(3);
106  vector<GlobalPoint> points(3);
107 
108  recHits[0] = (*ip).inner()->hit();
109  recHits[1] = (*ip).outer()->hit();
110 
111 #ifdef Debug
112  cerr << " RecHits " + HitInfo::getInfo(*recHits[0]) +
113  HitInfo::getInfo(*recHits[1]) << endl;
114 #endif
115 
116  for(int i=0; i<2; i++)
117  points[i] = getGlobalPosition(recHits[i]);
118 
119  // Initialize helix prediction
121  thePrediction(region,
122  points[0],points[1], es,
123  nSigMultipleScattering,maxAngleRatio,builderName);
124 
125  // Look at all layers
126  for(int il=0; il<size; il++)
127  {
128  const DetLayer * layer = theLayers[il].detLayer();
129 
130 #ifdef Debug
131  cerr << " check layer " << layer->subDetector()
132  << " " << layer->location() << endl;
133 #endif
134 
135  // Get ranges for the third hit
136  float phi[2],rz[2];
137  thePrediction.getRanges(layer, phi,rz);
138 
139  PixelRecoRange<float> phiRange(phi[0] , phi[1] );
140  PixelRecoRange<float> rzRange( rz[0] - rzTolerance, rz[1] + rzTolerance);
141 
142  // Get third hit candidates from cache
144  vector<Hit> thirdHits = thirdHitMap[il]->hits(phiRange.min(),phiRange.max());
145  typedef vector<Hit>::const_iterator IH;
146 
147  for (IH th=thirdHits.begin(), eh=thirdHits.end(); th < eh; ++th)
148  {
149  // Fill rechit and point
150  recHits[2] = (*th)->hit();
151  points[2] = getGlobalPosition(recHits[2]);
152 
153 #ifdef Debug
154  cerr << " third hit " + HitInfo::getInfo(*recHits[2]) << endl;
155 #endif
156 
157  // Check if third hit is compatible with multiple scattering
158  vector<GlobalVector> globalDirs;
159  if(thePrediction.isCompatibleWithMultipleScattering
160  (points[2], recHits, globalDirs, es) == false)
161  {
162 #ifdef Debug
163  cerr << " not compatible: multiple scattering" << endl;
164 #endif
165  if(checkMultipleScattering) continue;
166  }
167 
168  // Convert to localDirs
169 /*
170  vector<LocalVector> localDirs;
171  vector<GlobalVector>::const_iterator globalDir = globalDirs.begin();
172  for(vector<const TrackingRecHit *>::const_iterator
173  recHit = recHits.begin();
174  recHit != recHits.end(); recHit++)
175  {
176  localDirs.push_back(theTracker->idToDet(
177  (*recHit)->geographicalId())->toLocal(*globalDir));
178  globalDir++;
179  }
180 */
181 
182  // Check if the cluster shapes are compatible with thrusts
183  if(checkClusterShape)
184  {
185  if(! theFilter->checkTrack(recHits,globalDirs,tTopo))
186  {
187 #ifdef Debug
188  cerr << " not compatible: cluster shape" << endl;
189 #endif
190  continue;
191  }
192  }
193 
194  // All checks passed, put triplet back
195  result.push_back(OrderedHitTriplet((*ip).inner(),(*ip).outer(),*th));
196  }
197  }
198  }
199  delete [] thirdHitMap;
200 
201  return;
202 }
203 
204 
std::vector< Hit > hits(float phiMin, float phiMax) const
int i
Definition: DBlmapReader.cc:9
virtual HitPairGenerator * clone() const =0
T max() const
void getTracker(const edm::EventSetup &es)
virtual Location location() const =0
Which part of the detector (barrel, endcap)
void setSeedingLayers(SeedingLayerSetsHits::SeedingLayerSet pairLayers, std::vector< SeedingLayerSetsHits::SeedingLayer > thirdLayers) override
virtual SubDetector subDetector() const =0
The type of detector (PixelBarrel, PixelEndcap, TIB, TOB, TID, TEC, CSC, DT, RPCBarrel, RPCEndcap)
virtual unsigned int size() const
T min() const
void init(const HitPairGenerator &pairs, LayerCacheType *layerCache) override
tuple result
Definition: query.py:137
Definition: DetId.h:18
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
virtual void hitTriplets(const TrackingRegion &region, OrderedHitTriplets &trs, const edm::Event &ev, const edm::EventSetup &es)
TransientTrackingRecHit::ConstRecHitPointer Hit
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.
Definition: DDAxes.h:10