CMS 3D CMS Logo

RPCRecHitStandardAlgo.cc
Go to the documentation of this file.
1 /*
2  * See header file for a description of this class.
3  *
4  * \author M. Maggi -- INFN
5  */
6 
7 #include "RPCCluster.h"
16 
17 // First Step
19  const RPCCluster& cluster,
22  float& time, float& timeErr) const
23 {
24  // Get Average Strip position
25  const float fstrip = (roll.centreOfStrip(cluster.firstStrip())).x();
26  const float lstrip = (roll.centreOfStrip(cluster.lastStrip())).x();
27  const float centreOfCluster = (fstrip + lstrip)/2;
28  const double y = cluster.hasY() ? cluster.y() : 0;
29  Point = LocalPoint(centreOfCluster, y, 0);
30 
31  if ( !cluster.hasY() ) {
32  error = LocalError(roll.localError((cluster.firstStrip()+cluster.lastStrip())/2.));
33  }
34  else {
35  // Use the default one for local x error
36  float ex2 = roll.localError((cluster.firstStrip()+cluster.lastStrip())/2.).xx();
37  // Maximum estimate of local y error, (distance to the boundary)/sqrt(3)
38  // which gives consistent error to the default one at y=0
39  const float stripLen = roll.specificTopology().stripLength();
40  const float maxDy = stripLen/2 - std::abs(cluster.y());
41 
42  // Apply x-position correction for the endcap
43  if ( roll.id().region() != 0 ) {
44  const auto& topo = dynamic_cast<const TrapezoidalStripTopology&>(roll.topology());
45  const double angle = topo.stripAngle((cluster.firstStrip()+cluster.lastStrip())/2.);
46  const double x = centreOfCluster - y*std::tan(angle);
47  Point = LocalPoint(x, y, 0);
48 
49  // rescale x-error by the change of local pitch
50  const double scale = topo.localPitch(Point)/topo.pitch();
51  ex2 *= scale*scale;
52  }
53 
54  error = LocalError(ex2, 0, maxDy*maxDy/3.);
55  }
56 
57  if ( cluster.hasTime() ) {
58  time = cluster.time();
59  timeErr = cluster.timeRMS();
60  }
61  else {
62  time = 0;
63  timeErr = -1;
64  }
65 
66  return true;
67 }
68 
70  const RPCCluster& cl,
71  const float& angle,
72  const GlobalPoint& globPos,
75  float& time, float& timeErr) const
76 {
77  this->compute(roll,cl,Point,error,time,timeErr);
78  return true;
79 }
80 
float xx() const
Definition: LocalError.h:24
virtual float stripLength() const =0
LocalPoint centreOfStrip(int strip) const
Definition: RPCRoll.cc:52
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
float timeRMS() const
Definition: RPCCluster.cc:29
const StripTopology & specificTopology() const
Definition: RPCRoll.cc:107
bool hasTime() const
Definition: RPCCluster.cc:27
const Topology & topology() const override
Definition: RPCRoll.cc:30
RPCDetId id() const
Definition: RPCRoll.cc:24
LocalError localError(float strip) const
Definition: RPCRoll.cc:65
math::XYZPoint Point
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool hasY() const
Definition: RPCCluster.cc:31
bool compute(const RPCRoll &roll, const RPCCluster &cluster, LocalPoint &point, LocalError &error, float &time, float &timeErr) const override
standard local recHit computation
float y() const
Definition: RPCCluster.cc:32
int lastStrip() const
Definition: RPCCluster.cc:23
int firstStrip() const
Definition: RPCCluster.cc:22
float time() const
Definition: RPCCluster.cc:28
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:63