CMS 3D CMS Logo

VectorHitBuilderAlgorithmBase.cc
Go to the documentation of this file.
5 
9 
11 
13  const edm::ParameterSet& conf,
14  const TrackerGeometry* tkGeomProd,
15  const TrackerTopology* tkTopoProd,
17  : tkGeom_(tkGeomProd),
18  tkTopo_(tkTopoProd),
19  cpe_(cpeProd),
20  nMaxVHforeachStack_(conf.getParameter<int>("maxVectorHitsInAStack")),
21  barrelCut_(conf.getParameter<std::vector<double> >("BarrelCut")),
22  endcapCut_(conf.getParameter<std::vector<double> >("EndcapCut")),
23  cpeTag_(conf.getParameter<edm::ESInputTag>("CPE")) {}
24 
26  const Point3DBase<float, LocalTag>& lPosClu_low,
27  const PixelGeomDetUnit* geomDetUnit_upp,
28  const Point3DBase<float, LocalTag>& lPosClu_upp) const {
29  double parallCorr = 0.0;
30  Global3DPoint origin(0, 0, 0);
31  Global3DPoint gPosClu_low = geomDetUnit_low->surface().toGlobal(lPosClu_low);
32  GlobalVector gV = gPosClu_low - origin;
33  LogTrace("VectorHitsBuilderValidation") << " global vector passing to the origin:" << gV;
34 
35  LocalVector lV = geomDetUnit_low->surface().toLocal(gV);
36  LogTrace("VectorHitsBuilderValidation")
37  << " local vector passing to the origin (in the lower detector system of reference):" << lV;
38  LocalVector lV_norm = lV / lV.z();
39  LogTrace("VectorHitsBuilderValidation")
40  << " normalized local vector passing to the origin (in low the lower detector system of reference):" << lV_norm;
41 
42  Global3DPoint gPosClu_upp = geomDetUnit_upp->surface().toGlobal(lPosClu_upp);
43  Local3DPoint lPosClu_uppInLow = geomDetUnit_low->surface().toLocal(gPosClu_upp);
44  parallCorr = lV_norm.x() * lPosClu_uppInLow.z();
45 
46  return parallCorr;
47 }
48 
50  int nCluster = 0;
51  for (const auto& DSViter : clusters) {
52  // Loop over the clusters in the detector unit
53  for (const auto& clustIt : DSViter) {
54  nCluster++;
55  // get the detector unit's id
56  const GeomDetUnit* geomDetUnit(tkGeom_->idToDetUnit(DSViter.detId()));
57  if (!geomDetUnit)
58  return;
59  printCluster(geomDetUnit, &clustIt);
60  }
61  }
62  LogDebug("VectorHitBuilder") << " Number of input clusters: " << nCluster << std::endl;
63 }
64 
66  const Phase2TrackerCluster1D* clustIt) const {
67  if (!geomDetUnit)
68  return;
69  const PixelGeomDetUnit* pixelGeomDetUnit = dynamic_cast<const PixelGeomDetUnit*>(geomDetUnit);
70  const PixelTopology& topol = pixelGeomDetUnit->specificTopology();
71  if (!pixelGeomDetUnit)
72  return;
73 
74  unsigned int layer = tkTopo_->layer(geomDetUnit->geographicalId());
75  unsigned int module = tkTopo_->module(geomDetUnit->geographicalId());
76  LogTrace("VectorHitBuilder") << "Layer:" << layer << " and DetId: " << geomDetUnit->geographicalId().rawId()
77  << std::endl;
80  LogTrace("VectorHitBuilder") << "Pixel cluster (module:" << module << ") " << std::endl;
82  LogTrace("VectorHitBuilder") << "Strip cluster (module:" << module << ") " << std::endl;
83  else
84  LogTrace("VectorHitBuilder") << "no module?!" << std::endl;
85  LogTrace("VectorHitBuilder") << "with pitch:" << topol.pitch().first << " , " << topol.pitch().second << std::endl;
86  LogTrace("VectorHitBuilder") << " and width:" << pixelGeomDetUnit->surface().bounds().width()
87  << " , lenght:" << pixelGeomDetUnit->surface().bounds().length() << std::endl;
88 
89  auto&& lparams = cpe_->localParameters(*clustIt, *pixelGeomDetUnit);
90  Global3DPoint gparams = pixelGeomDetUnit->surface().toGlobal(lparams.first);
91 
92  LogTrace("VectorHitBuilder") << "\t global pos " << gparams << std::endl;
93  LogTrace("VectorHitBuilder") << "\t local pos " << lparams.first << "with err " << lparams.second << std::endl;
94  LogTrace("VectorHitBuilder") << std::endl;
95 
96  return;
97 }
virtual float length() const =0
const TrackerGeomDet * idToDetUnit(DetId) const override
Return the pointer to the GeomDetUnit corresponding to a given DetId.
VectorHitBuilderAlgorithmBase(const edm::ParameterSet &, const TrackerGeometry *, const TrackerTopology *, const ClusterParameterEstimator< Phase2TrackerCluster1D > *)
T z() const
Definition: PV3DBase.h:61
LocalPoint toLocal(const GlobalPoint &gp) const
unsigned int layer(const DetId &id) const
#define LogTrace(id)
unsigned int module(const DetId &id) const
virtual LocalValues localParameters(const T &, const GeomDetUnit &) const =0
T x() const
Definition: PV3DBase.h:59
ModuleType getDetectorType(DetId) const
const ClusterParameterEstimator< Phase2TrackerCluster1D > * cpe_
void printClusters(const edmNew::DetSetVector< Phase2TrackerCluster1D > &clusters) const
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:64
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
void printCluster(const GeomDet *geomDetUnit, const Phase2TrackerCluster1D *cluster) const
double computeParallaxCorrection(const PixelGeomDetUnit *, const Point3DBase< float, LocalTag > &, const PixelGeomDetUnit *, const Point3DBase< float, LocalTag > &) const
HLT enums.
virtual const PixelTopology & specificTopology() const
Returns a reference to the pixel proxy topology.
virtual std::pair< float, float > pitch() const =0
virtual float width() const =0
#define LogDebug(id)
const Bounds & bounds() const
Definition: Surface.h:87