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,
23  float& timeErr) const {
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  } else {
34  // Use the default one for local x error
35  float ex2 = roll.localError((cluster.firstStrip() + cluster.lastStrip()) / 2.).xx();
36  // Maximum estimate of local y error, (distance to the boundary)/sqrt(3)
37  // which gives consistent error to the default one at y=0
38  const float stripLen = roll.specificTopology().stripLength();
39  const float maxDy = stripLen / 2 - std::abs(cluster.y());
40 
41  // Apply x-position correction for the endcap
42  if (roll.id().region() != 0) {
43  const auto& topo = dynamic_cast<const TrapezoidalStripTopology&>(roll.topology());
44  const double angle = topo.stripAngle((cluster.firstStrip() + cluster.lastStrip()) / 2.);
45  const double x = centreOfCluster - y * std::tan(angle);
46  Point = LocalPoint(x, y, 0);
47 
48  // rescale x-error by the change of local pitch
49  const double scale = topo.localPitch(Point) / topo.pitch();
50  ex2 *= scale * scale;
51  }
52 
53  error = LocalError(ex2, 0, maxDy * maxDy / 3.);
54  }
55 
56  if (cluster.hasTime()) {
57  time = cluster.time();
58  timeErr = cluster.timeRMS();
59  } else {
60  time = 0;
61  timeErr = -1;
62  }
63 
64  return true;
65 }
66 
68  const RPCCluster& cl,
69  const float& angle,
70  const GlobalPoint& globPos,
73  float& time,
74  float& timeErr) const {
75  this->compute(roll, cl, Point, error, time, timeErr);
76  return true;
77 }
LocalError localError(float strip) const
Definition: RPCRoll.cc:33
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
float timeRMS() const
Definition: RPCCluster.cc:23
bool hasTime() const
Definition: RPCCluster.cc:21
int firstStrip() const
Definition: RPCCluster.cc:16
virtual float stripLength() const =0
bool hasY() const
Definition: RPCCluster.cc:27
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const StripTopology & specificTopology() const
Definition: RPCRoll.cc:49
float y() const
Definition: RPCCluster.cc:28
float time() const
Definition: RPCCluster.cc:22
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:53
bool compute(const RPCRoll &roll, const RPCCluster &cluster, LocalPoint &point, LocalError &error, float &time, float &timeErr) const override
standard local recHit computation
Structure Point Contains parameters of Gaussian fits to DMRs.
Definition: DMRtrends.cc:57
const Topology & topology() const override
Definition: RPCRoll.cc:18
RPCDetId id() const
Definition: RPCRoll.cc:16
int lastStrip() const
Definition: RPCCluster.cc:17
float xx() const
Definition: LocalError.h:22
LocalPoint centreOfStrip(int strip) const
Definition: RPCRoll.cc:26
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11