CMS 3D CMS Logo

SiStripFineDelayTLA.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 #include <iostream>
4 #include <TMath.h>
6 
24 
25 using namespace std;
27  conf_(conf)
28 {
29 }
30 
32 {
33  // get geometry
35  es.get<TrackerDigiGeometryRecord>().get(estracker);
36  tracker=&(*estracker);
37 }
38 
39 // Virtual destructor needed.
41 {
42 }
43 
44 std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > SiStripFineDelayTLA::findtrackangle(const std::vector<Trajectory>& trajVec)
45 {
46  if (!trajVec.empty()) {
47  return findtrackangle(trajVec.front()); }
48  std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > hitangleassociation;
49  return hitangleassociation;
50 }
51 
52 std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> > SiStripFineDelayTLA::findtrackangle(const Trajectory& traj)
53 {
54  std::vector<std::pair< std::pair<DetId, LocalPoint> ,float> >hitangleassociation;
55  std::vector<TrajectoryMeasurement> TMeas=traj.measurements();
56  std::vector<TrajectoryMeasurement>::iterator itm;
57  for (itm=TMeas.begin();itm!=TMeas.end();itm++){
58  TrajectoryStateOnSurface tsos=itm->updatedState();
59  auto thit=itm->recHit();
60  const SiStripMatchedRecHit2D* matchedhit=dynamic_cast<const SiStripMatchedRecHit2D*>((*thit).hit());
61  const SiStripRecHit2D* hit=dynamic_cast<const SiStripRecHit2D*>((*thit).hit());
62  LocalVector trackdirection=tsos.localDirection();
63  if(matchedhit){//if matched hit...
64  GluedGeomDet * gdet=(GluedGeomDet *)tracker->idToDet(matchedhit->geographicalId());
65  GlobalVector gtrkdir=gdet->toGlobal(trackdirection);
66  // trackdirection on mono det
67  // this the pointer to the mono hit of a matched hit
68  const SiStripRecHit2D monohit=matchedhit->monoHit();
69  const GeomDetUnit * monodet=gdet->monoDet();
70  LocalVector monotkdir=monodet->toLocal(gtrkdir);
71  if(monotkdir.z()!=0){
72  // the local angle (mono)
73  float localpitch = ((StripTopology*)(&monodet->topology()))->localPitch(tsos.localPosition());
74  float thickness = ((((((monohit.geographicalId())>>25)&0x7f)==0xd)||
75  ((((monohit.geographicalId())>>25)&0x7f)==0xe))&&
76  ((((monohit.geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
77  float angle = computeAngleCorr(monotkdir, localpitch, thickness);
78  hitangleassociation.push_back(make_pair(make_pair(monohit.geographicalId(),monohit.localPosition()), angle));
79  // trackdirection on stereo det
80  // this the pointer to the stereo hit of a matched hit
81  const SiStripRecHit2D stereohit=matchedhit->stereoHit();
82  const GeomDetUnit * stereodet=gdet->stereoDet();
83  LocalVector stereotkdir=stereodet->toLocal(gtrkdir);
84  if(stereotkdir.z()!=0){
85  // the local angle (stereo)
86  float localpitch = ((StripTopology*)(&stereodet->topology()))->localPitch(tsos.localPosition());
87  float thickness = ((((((stereohit.geographicalId())>>25)&0x7f)==0xd)||
88  ((((stereohit.geographicalId())>>25)&0x7f)==0xe))&&
89  ((((stereohit.geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
90  float angle = computeAngleCorr(stereotkdir, localpitch, thickness);
91  hitangleassociation.push_back(make_pair(make_pair(stereohit.geographicalId(),stereohit.localPosition()), angle));
92  }
93  }
94  }
95  else if(hit){
97  // hit= pointer to the rechit
98  if(trackdirection.z()!=0){
99  // the local angle (single hit)
100  float localpitch = ((StripTopology*)(&det->topology()))->localPitch(tsos.localPosition());
101  float thickness = ((((((hit->geographicalId())>>25)&0x7f)==0xd)||
102  ((((hit->geographicalId())>>25)&0x7f)==0xe))&&
103  ((((hit->geographicalId())>>5)&0x7)>4)) ? 0.0500 : 0.0320;
104  float angle = computeAngleCorr(trackdirection, localpitch, thickness);
105  hitangleassociation.push_back(make_pair(make_pair(hit->geographicalId(),hit->localPosition()), angle));
106  }
107  }
108  }
109  return hitangleassociation;
110 }
111 
112 double SiStripFineDelayTLA::computeAngleCorr(const LocalVector& v, double pitch, double thickness)
113 {
114  double v_xy = sqrt(v.x()*v.x()+v.y()*v.y());
115  double L = fabs(thickness*v_xy/v.z());
116  double Lmax = fabs(pitch/v.x()*v_xy);
117  if(L<Lmax) {
118  LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax
119  << " Signal contained in strip. Correction is " << v.z()/v.mag();
120  return v.z()/v.mag();
121  } else {
122  LogDebug("SiStripFineDelayTLA ") << L << " vs " << Lmax
123  << " Signal not contained in strip. Correction is " << thickness/pitch*v.x()/v_xy*v.z()/v.mag()
124  << " instead of " << v.z()/v.mag();
125  return thickness/pitch*v.x()/v_xy*v.z()/v.mag();
126  }
127 }
128 
#define LogDebug(id)
SiStripFineDelayTLA(const edm::ParameterSet &conf)
const GeomDetUnit * monoDet() const
Definition: GluedGeomDet.h:20
virtual const Topology & topology() const
Definition: GeomDet.cc:81
LocalVector localDirection() const
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
void init(const edm::Event &e, const edm::EventSetup &c)
T y() const
Definition: PV3DBase.h:63
double computeAngleCorr(const LocalVector &v, double pitch, double thickness)
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
DataContainer const & measurements() const
Definition: Trajectory.h:196
std::vector< std::pair< std::pair< DetId, LocalPoint >,float > > findtrackangle(const std::vector< Trajectory > &traj)
T mag() const
Definition: PV3DBase.h:67
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
SiStripRecHit2D stereoHit() const
const T & get() const
Definition: EventSetup.h:58
SiStripRecHit2D monoHit() const
const TrackerGeomDet * idToDet(DetId) const override
LocalPoint localPosition() const final
DetId geographicalId() const
const TrackerGeometry * tracker
T x() const
Definition: PV3DBase.h:62
const GeomDetUnit * stereoDet() const
Definition: GluedGeomDet.h:21
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11