CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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
HitPairGeneratorFromLayerPair
pairGenerator () const
 
virtual ~HitTripletGeneratorFromPairAndLayers ()
 

Private Member Functions

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

Private Attributes

bool checkClusterShape
 
bool checkMultipleScattering
 
edm::ESGetToken
< TrackerGeometry,
TrackerDigiGeometryRecord
m_geomToken
 
edm::ESGetToken< MagneticField,
IdealMagneticFieldRecord
m_magfieldToken
 
edm::ESGetToken
< MultipleScatteringParametrisationMaker,
TrackerMultipleScatteringRecord
m_msmakerToken
 
edm::ESGetToken
< TrackerTopology,
TrackerTopologyRcd
m_topoToken
 
edm::ESGetToken
< TransientTrackingRecHitBuilder,
TransientRecHitRecord
m_ttrhBuilderToken
 
double maxAngleRatio
 
double nSigMultipleScattering
 
double rzTolerance
 
edm::EDGetTokenT
< SiPixelClusterShapeCache
theClusterShapeCacheToken
 
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
< HitPairGeneratorFromLayerPair
thePairGenerator
 

Detailed Description

Definition at line 33 of file PixelTripletLowPtGenerator.h.

Constructor & Destructor Documentation

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

Definition at line 24 of file PixelTripletLowPtGenerator.cc.

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

25  : HitTripletGeneratorFromPairAndLayers(), // no theMaxElement used in this class
29  m_ttrhBuilderToken(iC.esConsumes(edm::ESInputTag("", cfg.getParameter<string>("TTRHBuilder")))),
31  theTracker(nullptr),
33  iC.consumes<SiPixelClusterShapeCache>(cfg.getParameter<edm::InputTag>("clusterShapeCacheSrc"))) {
34  checkMultipleScattering = cfg.getParameter<bool>("checkMultipleScattering");
35  nSigMultipleScattering = cfg.getParameter<double>("nSigMultipleScattering");
36  checkClusterShape = cfg.getParameter<bool>("checkClusterShape");
37  rzTolerance = cfg.getParameter<double>("rzTolerance");
38  maxAngleRatio = cfg.getParameter<double>("maxAngleRatio");
39 }
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > m_geomToken
edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > m_magfieldToken
edm::EDGetTokenT< SiPixelClusterShapeCache > theClusterShapeCacheToken
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const TrackerGeometry * theTracker
edm::ESGetToken< MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord > m_msmakerToken
edm::ESGetToken< TrackerTopology, TrackerTopologyRcd > m_topoToken
edm::ESGetToken< TransientTrackingRecHitBuilder, TransientRecHitRecord > m_ttrhBuilderToken
PixelTripletLowPtGenerator::~PixelTripletLowPtGenerator ( )
override

Definition at line 42 of file PixelTripletLowPtGenerator.cc.

42 {}

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
const TrackerGeomDet * idToDet(DetId) const override
Definition: DetId.h:17
const TrackerGeometry * theTracker
DetId geographicalId() const
virtual LocalPoint localPosition() const =0
void PixelTripletLowPtGenerator::getTracker ( const edm::EventSetup es)
private

Definition at line 47 of file PixelTripletLowPtGenerator.cc.

References edm::EventSetup::getData(), m_geomToken, theFilter, and theTracker.

Referenced by hitTriplets().

47  {
48  if (theTracker == nullptr) {
49  // Get tracker geometry
51  }
52 
53  if (!theFilter) {
54  theFilter = std::make_unique<TripletFilter>(es);
55  }
56 }
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > m_geomToken
std::unique_ptr< TripletFilter > theFilter
bool getData(T &iHolder) const
Definition: EventSetup.h:128
const TrackerGeometry * theTracker
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 EcnaPython_AdcPeg12_S1_10_R170298_1_0_150_Dee0::cerr, checkClusterShape, checkMultipleScattering, edm::Event::getByToken(), edm::EventSetup::getData(), getGlobalPosition(), HitInfo::getInfo(), ThirdHitPrediction::getRanges(), getTracker(), RecHitsSortedInPhi::hits(), mps_fire::i, ThirdHitPrediction::isCompatibleWithMultipleScattering(), phase1PixelTopology::layer, DetLayer::location(), m_magfieldToken, m_msmakerToken, m_topoToken, m_ttrhBuilderToken, volumeBasedMagneticField_160812_cfi::magfield, PixelRecoRange< T >::max(), maxAngleRatio, PixelRecoRange< T >::min(), nSigMultipleScattering, phi, HLT_FULL_cff::recHits, HLT_FULL_cff::region, rzTolerance, findQualityFiles::size, DetLayer::subDetector(), theClusterShapeCacheToken, theFilter, HitTripletGeneratorFromPairAndLayers::theLayerCache, and HitTripletGeneratorFromPairAndLayers::thePairGenerator.

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

References Exception.

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

Member Data Documentation

bool PixelTripletLowPtGenerator::checkClusterShape
private

Definition at line 72 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

bool PixelTripletLowPtGenerator::checkMultipleScattering
private

Definition at line 71 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> PixelTripletLowPtGenerator::m_geomToken
private

Definition at line 54 of file PixelTripletLowPtGenerator.h.

Referenced by getTracker().

edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> PixelTripletLowPtGenerator::m_magfieldToken
private

Definition at line 56 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets().

edm::ESGetToken<MultipleScatteringParametrisationMaker, TrackerMultipleScatteringRecord> PixelTripletLowPtGenerator::m_msmakerToken
private

Definition at line 58 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets().

edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> PixelTripletLowPtGenerator::m_topoToken
private

Definition at line 55 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets().

edm::ESGetToken<TransientTrackingRecHitBuilder, TransientRecHitRecord> PixelTripletLowPtGenerator::m_ttrhBuilderToken
private

Definition at line 57 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets().

double PixelTripletLowPtGenerator::maxAngleRatio
private

Definition at line 69 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::nSigMultipleScattering
private

Definition at line 67 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

double PixelTripletLowPtGenerator::rzTolerance
private

Definition at line 68 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets(), and PixelTripletLowPtGenerator().

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

Definition at line 66 of file PixelTripletLowPtGenerator.h.

Referenced by hitTriplets().

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

Definition at line 64 of file PixelTripletLowPtGenerator.h.

Referenced by getTracker(), and hitTriplets().

const TrackerGeometry* PixelTripletLowPtGenerator::theTracker
private

Definition at line 63 of file PixelTripletLowPtGenerator.h.

Referenced by getGlobalPosition(), and getTracker().