CMS 3D CMS Logo

PixelTripletLowPtGenerator.cc
Go to the documentation of this file.
5 
8 
11 
14 
15 #undef Debug
16 
17 using namespace std;
18 
19 /*****************************************************************************/
21  : HitTripletGeneratorFromPairAndLayers(), // no theMaxElement used in this class
23  m_topoToken(iC.esConsumes<TrackerTopology, TrackerTopologyRcd>()),
24  theTracker(nullptr),
25  theClusterShapeCacheToken(
26  iC.consumes<SiPixelClusterShapeCache>(cfg.getParameter<edm::InputTag>("clusterShapeCacheSrc"))) {
27  checkMultipleScattering = cfg.getParameter<bool>("checkMultipleScattering");
28  nSigMultipleScattering = cfg.getParameter<double>("nSigMultipleScattering");
29  checkClusterShape = cfg.getParameter<bool>("checkClusterShape");
30  rzTolerance = cfg.getParameter<double>("rzTolerance");
31  maxAngleRatio = cfg.getParameter<double>("maxAngleRatio");
32  builderName = cfg.getParameter<string>("TTRHBuilder");
33 }
34 
35 /*****************************************************************************/
37 
38 /*****************************************************************************/
39 
40 /*****************************************************************************/
42  if (theTracker == nullptr) {
43  // Get tracker geometry
45  }
46 
47  if (!theFilter) {
48  theFilter = std::make_unique<TripletFilter>(es);
49  }
50 }
51 
52 /*****************************************************************************/
54  DetId detId = recHit->geographicalId();
55 
56  return theTracker->idToDet(detId)->toGlobal(recHit->localPosition());
57 }
58 
59 /*****************************************************************************/
62  const edm::Event& ev,
63  const edm::EventSetup& es,
64  const SeedingLayerSetsHits::SeedingLayerSet& pairLayers,
65  const std::vector<SeedingLayerSetsHits::SeedingLayer>& thirdLayers) {
66  //Retrieve tracker topology from geometry
67  const TrackerTopology* tTopo = &es.getData(m_topoToken);
68 
71 
72  // Generate pairs
73  OrderedHitPairs pairs;
74  pairs.reserve(30000);
75  thePairGenerator->hitPairs(region, pairs, ev, es, pairLayers);
76 
77  if (pairs.empty())
78  return;
79 
80  int size = thirdLayers.size();
81 
82  // Set aliases
83  const RecHitsSortedInPhi** thirdHitMap = new const RecHitsSortedInPhi*[size];
84  for (int il = 0; il < size; il++)
85  thirdHitMap[il] = &(*theLayerCache)(thirdLayers[il], region, es);
86 
87  // Get tracker
88  getTracker(es);
89 
90  // Look at all generated pairs
91  for (OrderedHitPairs::const_iterator ip = pairs.begin(); ip != pairs.end(); ip++) {
92  // Fill rechits and points
93  vector<const TrackingRecHit*> recHits(3);
94  vector<GlobalPoint> points(3);
95 
96  recHits[0] = (*ip).inner()->hit();
97  recHits[1] = (*ip).outer()->hit();
98 
99 #ifdef Debug
100  cerr << " RecHits " + HitInfo::getInfo(*recHits[0]) + HitInfo::getInfo(*recHits[1]) << endl;
101 #endif
102 
103  for (int i = 0; i < 2; i++)
105 
106  // Initialize helix prediction
107  ThirdHitPrediction thePrediction(
109 
110  // Look at all layers
111  for (int il = 0; il < size; il++) {
112  const DetLayer* layer = thirdLayers[il].detLayer();
113 
114 #ifdef Debug
115  cerr << " check layer " << layer->subDetector() << " " << layer->location() << endl;
116 #endif
117 
118  // Get ranges for the third hit
119  float phi[2], rz[2];
120  thePrediction.getRanges(layer, phi, rz);
121 
122  PixelRecoRange<float> phiRange(phi[0], phi[1]);
123  PixelRecoRange<float> rzRange(rz[0] - rzTolerance, rz[1] + rzTolerance);
124 
125  // Get third hit candidates from cache
127  vector<Hit> thirdHits = thirdHitMap[il]->hits(phiRange.min(), phiRange.max());
128  typedef vector<Hit>::const_iterator IH;
129 
130  for (IH th = thirdHits.begin(), eh = thirdHits.end(); th < eh; ++th) {
131  // Fill rechit and point
132  recHits[2] = (*th)->hit();
134 
135 #ifdef Debug
136  cerr << " third hit " + HitInfo::getInfo(*recHits[2]) << endl;
137 #endif
138 
139  // Check if third hit is compatible with multiple scattering
140  vector<GlobalVector> globalDirs;
141  if (thePrediction.isCompatibleWithMultipleScattering(points[2], recHits, globalDirs, es) == false) {
142 #ifdef Debug
143  cerr << " not compatible: multiple scattering" << endl;
144 #endif
146  continue;
147  }
148 
149  // Convert to localDirs
150  /*
151  vector<LocalVector> localDirs;
152  vector<GlobalVector>::const_iterator globalDir = globalDirs.begin();
153  for(vector<const TrackingRecHit *>::const_iterator
154  recHit = recHits.begin();
155  recHit != recHits.end(); recHit++)
156  {
157  localDirs.push_back(theTracker->idToDet(
158  (*recHit)->geographicalId())->toLocal(*globalDir));
159  globalDir++;
160  }
161 */
162 
163  // Check if the cluster shapes are compatible with thrusts
164  if (checkClusterShape) {
165  if (!theFilter->checkTrack(recHits, globalDirs, tTopo, *clusterShapeCache)) {
166 #ifdef Debug
167  cerr << " not compatible: cluster shape" << endl;
168 #endif
169  continue;
170  }
171  }
172 
173  // All checks passed, put triplet back
174  result.push_back(OrderedHitTriplet((*ip).inner(), (*ip).outer(), *th));
175  }
176  }
177  }
178  delete[] thirdHitMap;
179 
180  return;
181 }
184  const edm::EventSetup& es,
185  const HitDoublets& doublets,
186  const RecHitsSortedInPhi** thirdHitMap,
187  const std::vector<const DetLayer*>& thirdLayerDetLayer,
188  const int nThirdLayers) {
189  throw cms::Exception("Error") << "PixelTripletLowPtGenerator::hitTriplets is not implemented \n";
190 }
SiPixelPhase1OnlineDQM_cff.clusterShapeCache
clusterShapeCache
Definition: SiPixelPhase1OnlineDQM_cff.py:120
TrackerGeometry::idToDet
const TrackerGeomDet * idToDet(DetId) const override
Definition: TrackerGeometry.cc:193
ThirdHitPrediction.h
PixelRecoPointRZ.h
OrderedHitPairs
Definition: OrderedHitPairs.h:8
ThirdHitPrediction::isCompatibleWithMultipleScattering
bool isCompatibleWithMultipleScattering(GlobalPoint g3, const std::vector< const TrackingRecHit * > &h, std::vector< GlobalVector > &localDirs, const edm::EventSetup &es)
Definition: ThirdHitPrediction.cc:345
mps_fire.i
i
Definition: mps_fire.py:428
SiPixelClusterShapeCache
Definition: SiPixelClusterShapeCache.h:43
ESHandle.h
DetLayer
Definition: DetLayer.h:21
HLT_FULL_cff.points
points
Definition: HLT_FULL_cff.py:21449
edm
HLT enums.
Definition: AlignableModifier.h:19
TrackerTopology
Definition: TrackerTopology.h:16
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89281
PixelTripletLowPtGenerator::theFilter
std::unique_ptr< TripletFilter > theFilter
Definition: PixelTripletLowPtGenerator.h:56
HLT_FULL_cff.doublets
doublets
Definition: HLT_FULL_cff.py:9865
PixelTripletLowPtGenerator.h
PixelTripletLowPtGenerator::checkMultipleScattering
bool checkMultipleScattering
Definition: PixelTripletLowPtGenerator.h:64
edm::Handle< SiPixelClusterShapeCache >
rpcPointValidation_cfi.recHit
recHit
Definition: rpcPointValidation_cfi.py:7
PixelRecoRange::min
T min() const
Definition: PixelRecoRange.h:25
DetId
Definition: DetId.h:17
PixelTripletLowPtGenerator::PixelTripletLowPtGenerator
PixelTripletLowPtGenerator(const edm::ParameterSet &cfg, edm::ConsumesCollector &iC)
Definition: PixelTripletLowPtGenerator.cc:20
PixelTripletLowPtGenerator::checkClusterShape
bool checkClusterShape
Definition: PixelTripletLowPtGenerator.h:65
PixelTripletLowPtGenerator::~PixelTripletLowPtGenerator
~PixelTripletLowPtGenerator() override
Definition: PixelTripletLowPtGenerator.cc:36
PixelTripletLowPtGenerator::m_topoToken
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > m_topoToken
Definition: PixelTripletLowPtGenerator.h:50
HitTripletGeneratorFromPairAndLayers::thePairGenerator
std::unique_ptr< HitPairGeneratorFromLayerPair > thePairGenerator
Definition: HitTripletGeneratorFromPairAndLayers.h:55
PixelRecoRange::max
T max() const
Definition: PixelRecoRange.h:26
PixelTripletLowPtGenerator::nSigMultipleScattering
double nSigMultipleScattering
Definition: PixelTripletLowPtGenerator.h:59
TrackerDigiGeometryRecord
Definition: TrackerDigiGeometryRecord.h:15
HitTripletGeneratorFromPairAndLayers
Definition: HitTripletGeneratorFromPairAndLayers.h:25
PixelTripletLowPtGenerator::builderName
std::string builderName
Definition: PixelTripletLowPtGenerator.h:63
RecHitsSortedInPhi::Hit
BaseTrackerRecHit const * Hit
Definition: RecHitsSortedInPhi.h:19
ThirdHitPrediction::getRanges
void getRanges(const DetLayer *layer, float phi[], float rz[])
Definition: ThirdHitPrediction.cc:304
Point3DBase< float, GlobalTag >
PixelTripletLowPtGenerator::rzTolerance
double rzTolerance
Definition: PixelTripletLowPtGenerator.h:60
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
FastTrackerRecHitMaskProducer_cfi.recHits
recHits
Definition: FastTrackerRecHitMaskProducer_cfi.py:8
RecHitsSortedInPhi
Definition: RecHitsSortedInPhi.h:17
PixelRecoRange< float >
GeomDet::toGlobal
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
HitDoublets
Definition: RecHitsSortedInPhi.h:124
edm::ParameterSet
Definition: ParameterSet.h:47
PixelTripletLowPtGenerator::theClusterShapeCacheToken
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
Definition: PixelTripletLowPtGenerator.h:58
PixelTripletLowPtGenerator::getTracker
void getTracker(const edm::EventSetup &es)
Definition: PixelTripletLowPtGenerator.cc:41
Event.h
OrderedHitTriplets
Definition: OrderedHitTriplets.h:9
HLT_FULL_cff.region
region
Definition: HLT_FULL_cff.py:88267
PixelTripletLowPtGenerator::hitTriplets
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
Definition: PixelTripletLowPtGenerator.cc:60
edm::EventSetup
Definition: EventSetup.h:58
TripletFilter.h
edm::EventSetup::getData
bool getData(T &iHolder) const
Definition: EventSetup.h:127
looper.cfg
cfg
Definition: looper.py:297
TrackingRecHit
Definition: TrackingRecHit.h:21
DDAxes::phi
HitTripletGeneratorFromPairAndLayers::theLayerCache
LayerCacheType * theLayerCache
Definition: HitTripletGeneratorFromPairAndLayers.h:56
std
Definition: JetResolutionObject.h:76
RecHitsSortedInPhi::hits
std::vector< Hit > hits(float phiMin, float phiMax) const
Definition: RecHitsSortedInPhi.cc:93
HitPairGeneratorFromLayerPair.h
ev
bool ev
Definition: Hydjet2Hadronizer.cc:95
Exception
Definition: hltDiff.cc:245
SeedingLayerSetsHits::SeedingLayerSet
Definition: SeedingLayerSetsHits.h:65
PixelTripletLowPtGenerator::theTracker
const TrackerGeometry * theTracker
Definition: PixelTripletLowPtGenerator.h:55
OrderedHitTriplet
Definition: OrderedHitTriplet.h:11
PixelTripletLowPtGenerator::m_geomToken
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > m_geomToken
Definition: PixelTripletLowPtGenerator.h:49
TrackingRegion
Definition: TrackingRegion.h:41
mps_fire.result
result
Definition: mps_fire.py:311
ConsumesCollector.h
PixelTripletLowPtGenerator::maxAngleRatio
double maxAngleRatio
Definition: PixelTripletLowPtGenerator.h:61
ThirdHitPrediction
Definition: ThirdHitPrediction.h:33
HitInfo.h
TrackerTopologyRcd
Definition: TrackerTopologyRcd.h:10
DeDxTools::esConsumes
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
edm::Event
Definition: Event.h:73
EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.cerr
cerr
Definition: EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0.py:8
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
Hit
SeedingHitSet::ConstRecHitPointer Hit
Definition: SeedGeneratorFromProtoTracksEDProducer.cc:34
SiPixelClusterShapeCache.h
PixelTripletLowPtGenerator::getGlobalPosition
GlobalPoint getGlobalPosition(const TrackingRecHit *recHit)
Definition: PixelTripletLowPtGenerator.cc:53
findQualityFiles.size
size
Write out results.
Definition: findQualityFiles.py:443
HitInfo::getInfo
static std::string getInfo(const DetId &id, const TrackerTopology *tTopo)
Definition: HitInfo.cc:19
TrackerGeometry
Definition: TrackerGeometry.h:14