File indexing completed on 2024-04-06 12:26:15
0001
0002
0003
0004
0005
0006
0007 #include "RPCCluster.h"
0008 #include "RPCRecHitStandardAlgo.h"
0009 #include "DataFormats/MuonDetId/interface/RPCDetId.h"
0010 #include "Geometry/RPCGeometry/interface/RPCRoll.h"
0011 #include "Geometry/CommonTopologies/interface/StripTopology.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/Utilities/interface/Exception.h"
0015 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
0016
0017
0018 bool RPCRecHitStandardAlgo::compute(const RPCRoll& roll,
0019 const RPCCluster& cluster,
0020 LocalPoint& Point,
0021 LocalError& error,
0022 float& time,
0023 float& timeErr) const {
0024
0025 const float fstrip = (roll.centreOfStrip(cluster.firstStrip())).x();
0026 const float lstrip = (roll.centreOfStrip(cluster.lastStrip())).x();
0027 const float centreOfCluster = (fstrip + lstrip) / 2;
0028 const double y = cluster.hasY() ? cluster.y() : 0;
0029 Point = LocalPoint(centreOfCluster, y, 0);
0030
0031 if (!cluster.hasY()) {
0032 error = LocalError(roll.localError((cluster.firstStrip() + cluster.lastStrip()) / 2.));
0033 } else {
0034
0035 float ex2 = roll.localError((cluster.firstStrip() + cluster.lastStrip()) / 2.).xx();
0036
0037
0038 const float stripLen = roll.specificTopology().stripLength();
0039 const float maxDy = stripLen / 2 - std::abs(cluster.y());
0040
0041
0042 if (roll.id().region() != 0) {
0043 const auto& topo = dynamic_cast<const TrapezoidalStripTopology&>(roll.topology());
0044 const double angle = topo.stripAngle((cluster.firstStrip() + cluster.lastStrip()) / 2.);
0045 const double x = centreOfCluster - y * std::tan(angle);
0046 Point = LocalPoint(x, y, 0);
0047
0048
0049 const double scale = topo.localPitch(Point) / topo.pitch();
0050 ex2 *= scale * scale;
0051 }
0052
0053 error = LocalError(ex2, 0, maxDy * maxDy / 3.);
0054 }
0055
0056 if (cluster.hasTime()) {
0057 time = cluster.time();
0058 timeErr = cluster.timeRMS();
0059 } else {
0060 time = 0;
0061 timeErr = -1;
0062 }
0063
0064 return true;
0065 }
0066
0067 bool RPCRecHitStandardAlgo::compute(const RPCRoll& roll,
0068 const RPCCluster& cl,
0069 const float& angle,
0070 const GlobalPoint& globPos,
0071 LocalPoint& Point,
0072 LocalError& error,
0073 float& time,
0074 float& timeErr) const {
0075 this->compute(roll, cl, Point, error, time, timeErr);
0076 return true;
0077 }