CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
SiStripFineDelayTLA Class Reference

#include <SiStripFineDelayTLA.h>

Public Member Functions

std::vector< std::pair< std::pair< DetId, LocalPoint >, float > > findtrackangle (const std::vector< Trajectory > &traj)
 
std::vector< std::pair< std::pair< DetId, LocalPoint >, float > > findtrackangle (const Trajectory &traj)
 
void init (const edm::Event &e, const edm::EventSetup &c)
 
 SiStripFineDelayTLA (const edm::ParameterSet &conf, edm::ConsumesCollector iC)
 
virtual ~SiStripFineDelayTLA ()
 

Private Member Functions

double computeAngleCorr (const LocalVector &v, double pitch, double thickness)
 

Private Attributes

edm::ParameterSet conf_
 
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecordtkGeomToken_
 
const TrackerGeometrytracker
 

Detailed Description

Definition at line 21 of file SiStripFineDelayTLA.h.

Constructor & Destructor Documentation

◆ SiStripFineDelayTLA()

SiStripFineDelayTLA::SiStripFineDelayTLA ( const edm::ParameterSet conf,
edm::ConsumesCollector  iC 
)
explicit

Definition at line 24 of file SiStripFineDelayTLA.cc.

◆ ~SiStripFineDelayTLA()

SiStripFineDelayTLA::~SiStripFineDelayTLA ( )
virtual

Definition at line 30 of file SiStripFineDelayTLA.cc.

30 {}

Member Function Documentation

◆ computeAngleCorr()

double SiStripFineDelayTLA::computeAngleCorr ( const LocalVector v,
double  pitch,
double  thickness 
)
private

Definition at line 107 of file SiStripFineDelayTLA.cc.

References dttmaxenums::L, LogDebug, mathSSE::sqrt(), ppsPixelTopologyESSourceRun2_cfi::thickness, and findQualityFiles::v.

Referenced by findtrackangle().

107  {
108  double v_xy = sqrt(v.x() * v.x() + v.y() * v.y());
109  double L = fabs(thickness * v_xy / v.z());
110  double Lmax = fabs(pitch / v.x() * v_xy);
111  if (L < Lmax) {
112  LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax << " Signal contained in strip. Correction is "
113  << v.z() / v.mag();
114  return v.z() / v.mag();
115  } else {
116  LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax << " Signal not contained in strip. Correction is "
117  << thickness / pitch * v.x() / v_xy * v.z() / v.mag() << " instead of "
118  << v.z() / v.mag();
119  return thickness / pitch * v.x() / v_xy * v.z() / v.mag();
120  }
121 }
T sqrt(T t)
Definition: SSEVec.h:23
#define LogDebug(id)

◆ findtrackangle() [1/2]

std::vector< std::pair< std::pair< DetId, LocalPoint >, float > > SiStripFineDelayTLA::findtrackangle ( const std::vector< Trajectory > &  traj)

Definition at line 32 of file SiStripFineDelayTLA.cc.

Referenced by SiStripFineDelayHit::detId().

33  {
34  if (!trajVec.empty()) {
35  return findtrackangle(trajVec.front());
36  }
37  std::vector<std::pair<std::pair<DetId, LocalPoint>, float> > hitangleassociation;
38  return hitangleassociation;
39 }
std::vector< std::pair< std::pair< DetId, LocalPoint >, float > > findtrackangle(const std::vector< Trajectory > &traj)

◆ findtrackangle() [2/2]

std::vector< std::pair< std::pair< DetId, LocalPoint >, float > > SiStripFineDelayTLA::findtrackangle ( const Trajectory traj)

Definition at line 41 of file SiStripFineDelayTLA.cc.

References angle(), computeAngleCorr(), TrackingRecHit::geographicalId(), TrackerGeometry::idToDet(), TrajectoryStateOnSurface::localDirection(), BaseTrackerRecHit::localPosition(), TrajectoryStateOnSurface::localPosition(), Trajectory::measurements(), GluedGeomDet::monoDet(), SiStripMatchedRecHit2D::monoHit(), GluedGeomDet::stereoDet(), SiStripMatchedRecHit2D::stereoHit(), ppsPixelTopologyESSourceRun2_cfi::thickness, GeomDet::toGlobal(), GeomDet::toLocal(), GeomDet::topology(), tracker, and PV3DBase< T, PVType, FrameType >::z().

42  {
43  std::vector<std::pair<std::pair<DetId, LocalPoint>, float> > hitangleassociation;
44  std::vector<TrajectoryMeasurement> TMeas = traj.measurements();
45  std::vector<TrajectoryMeasurement>::iterator itm;
46  for (itm = TMeas.begin(); itm != TMeas.end(); itm++) {
47  TrajectoryStateOnSurface tsos = itm->updatedState();
48  auto thit = itm->recHit();
49  const SiStripMatchedRecHit2D* matchedhit = dynamic_cast<const SiStripMatchedRecHit2D*>((*thit).hit());
50  const SiStripRecHit2D* hit = dynamic_cast<const SiStripRecHit2D*>((*thit).hit());
51  LocalVector trackdirection = tsos.localDirection();
52  if (matchedhit) { //if matched hit...
53  const GluedGeomDet* gdet = static_cast<const GluedGeomDet*>(tracker->idToDet(matchedhit->geographicalId()));
54  GlobalVector gtrkdir = gdet->toGlobal(trackdirection);
55  // trackdirection on mono det
56  // this the pointer to the mono hit of a matched hit
57  const SiStripRecHit2D monohit = matchedhit->monoHit();
58  const GeomDetUnit* monodet = gdet->monoDet();
59  LocalVector monotkdir = monodet->toLocal(gtrkdir);
60  if (monotkdir.z() != 0) {
61  // the local angle (mono)
62  float localpitch = static_cast<const StripTopology&>(monodet->topology()).localPitch(tsos.localPosition());
63  float thickness = ((((((monohit.geographicalId()) >> 25) & 0x7f) == 0xd) ||
64  ((((monohit.geographicalId()) >> 25) & 0x7f) == 0xe)) &&
65  ((((monohit.geographicalId()) >> 5) & 0x7) > 4))
66  ? 0.0500
67  : 0.0320;
68  float angle = computeAngleCorr(monotkdir, localpitch, thickness);
69  hitangleassociation.push_back(make_pair(make_pair(monohit.geographicalId(), monohit.localPosition()), angle));
70  // trackdirection on stereo det
71  // this the pointer to the stereo hit of a matched hit
72  const SiStripRecHit2D stereohit = matchedhit->stereoHit();
73  const GeomDetUnit* stereodet = gdet->stereoDet();
74  LocalVector stereotkdir = stereodet->toLocal(gtrkdir);
75  if (stereotkdir.z() != 0) {
76  // the local angle (stereo)
77  float localpitch = static_cast<const StripTopology&>(stereodet->topology()).localPitch(tsos.localPosition());
78  float thickness = ((((((stereohit.geographicalId()) >> 25) & 0x7f) == 0xd) ||
79  ((((stereohit.geographicalId()) >> 25) & 0x7f) == 0xe)) &&
80  ((((stereohit.geographicalId()) >> 5) & 0x7) > 4))
81  ? 0.0500
82  : 0.0320;
83  float angle = computeAngleCorr(stereotkdir, localpitch, thickness);
84  hitangleassociation.push_back(
85  make_pair(make_pair(stereohit.geographicalId(), stereohit.localPosition()), angle));
86  }
87  }
88  } else if (hit) {
89  const GeomDetUnit* det = tracker->idToDet(hit->geographicalId());
90  // hit= pointer to the rechit
91  if (trackdirection.z() != 0) {
92  // the local angle (single hit)
93  float localpitch = static_cast<const StripTopology&>(det->topology()).localPitch(tsos.localPosition());
94  float thickness =
95  ((((((hit->geographicalId()) >> 25) & 0x7f) == 0xd) || ((((hit->geographicalId()) >> 25) & 0x7f) == 0xe)) &&
96  ((((hit->geographicalId()) >> 5) & 0x7) > 4))
97  ? 0.0500
98  : 0.0320;
99  float angle = computeAngleCorr(trackdirection, localpitch, thickness);
100  hitangleassociation.push_back(make_pair(make_pair(hit->geographicalId(), hit->localPosition()), angle));
101  }
102  }
103  }
104  return hitangleassociation;
105 }
SiStripRecHit2D stereoHit() const
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
T z() const
Definition: PV3DBase.h:61
virtual const Topology & topology() const
Definition: GeomDet.cc:67
double computeAngleCorr(const LocalVector &v, double pitch, double thickness)
DataContainer const & measurements() const
Definition: Trajectory.h:178
LocalVector localDirection() const
const TrackerGeomDet * idToDet(DetId) const override
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
const GeomDetUnit * monoDet() const
Definition: GluedGeomDet.h:19
DetId geographicalId() const
const GeomDetUnit * stereoDet() const
Definition: GluedGeomDet.h:20
LocalPoint localPosition() const override
SiStripRecHit2D monoHit() const
const TrackerGeometry * tracker
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11

◆ init()

void SiStripFineDelayTLA::init ( const edm::Event e,
const edm::EventSetup c 
)

Definition at line 27 of file SiStripFineDelayTLA.cc.

References edm::EventSetup::getData(), tkGeomToken_, and tracker.

Referenced by SiStripFineDelayHit::produce().

27 { tracker = &es.getData(tkGeomToken_); }
edm::ESGetToken< TrackerGeometry, TrackerDigiGeometryRecord > tkGeomToken_
const TrackerGeometry * tracker

Member Data Documentation

◆ conf_

edm::ParameterSet SiStripFineDelayTLA::conf_
private

Definition at line 34 of file SiStripFineDelayTLA.h.

◆ tkGeomToken_

edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> SiStripFineDelayTLA::tkGeomToken_
private

Definition at line 36 of file SiStripFineDelayTLA.h.

Referenced by init().

◆ tracker

const TrackerGeometry* SiStripFineDelayTLA::tracker
private

Definition at line 35 of file SiStripFineDelayTLA.h.

Referenced by findtrackangle(), and init().