File indexing completed on 2024-04-06 12:26:14
0001 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
0002 #include "Geometry/CSCGeometry/interface/CSCGeometry.h"
0003 #include "DataFormats/CSCRecHit/interface/CSCSegmentCollection.h"
0004 #include "Geometry/CommonDetUnit/interface/GeomDet.h"
0005 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0006 #include "Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h"
0007 #include "FWCore/Framework/interface/ESHandle.h"
0008 #include "Geometry/RPCGeometry/interface/RPCGeomServ.h"
0009 #include "DataFormats/RPCRecHit/interface/RPCRecHit.h"
0010 #include "DataFormats/RPCRecHit/interface/RPCRecHitCollection.h"
0011 #include "CSCSegtoRPC.h"
0012 #include "CSCStationIndex.h"
0013 #include "CSCObjectMap.h"
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015
0016 CSCSegtoRPC::CSCSegtoRPC(edm::ConsumesCollector iC, const edm::ParameterSet& iConfig)
0017 : rpcGeoToken_(iC.esConsumes()), cscGeoToken_(iC.esConsumes()), cscMapToken_(iC.esConsumes()) {
0018 minBX = iConfig.getParameter<int>("minBX");
0019 maxBX = iConfig.getParameter<int>("maxBX");
0020 }
0021
0022 std::unique_ptr<RPCRecHitCollection> CSCSegtoRPC::thePoints(const CSCSegmentCollection* allCSCSegments,
0023 const edm::EventSetup& iSetup,
0024 bool debug,
0025 double eyr) {
0026 edm::ESHandle<RPCGeometry> rpcGeo = iSetup.getHandle(rpcGeoToken_);
0027 edm::ESHandle<CSCGeometry> cscGeo = iSetup.getHandle(cscGeoToken_);
0028 edm::ESHandle<CSCObjectMap> cscMap = iSetup.getHandle(cscMapToken_);
0029
0030 double MaxD = 80.;
0031
0032 if (debug)
0033 LogDebug("CSCSegtoRPC") << "CSC \t Number of CSC Segments in this event = " << allCSCSegments->size() << std::endl;
0034
0035 auto _ThePoints = std::make_unique<RPCRecHitCollection>();
0036 edm::OwnVector<RPCRecHit> RPCPointVector;
0037
0038 if (allCSCSegments->size() == 0) {
0039 if (debug)
0040 LogDebug("CSCSegtoRPC") << "CSC 0 segments skiping event" << std::endl;
0041 } else {
0042 std::map<CSCDetId, int> CSCSegmentsCounter;
0043 CSCSegmentCollection::const_iterator segment;
0044
0045 for (segment = allCSCSegments->begin(); segment != allCSCSegments->end(); ++segment) {
0046 CSCSegmentsCounter[segment->cscDetId()]++;
0047 }
0048
0049 float myTime = -9999.;
0050 float myTimeErr = -9999.;
0051 int myBx = -99;
0052 const float bunchCrossTimeDiff = 25.;
0053
0054 if (debug)
0055 LogDebug("CSCSegtoRPC") << "CSC \t loop over all the CSCSegments " << std::endl;
0056 for (segment = allCSCSegments->begin(); segment != allCSCSegments->end(); ++segment) {
0057 CSCDetId CSCId = segment->cscDetId();
0058
0059 myTime = segment->time();
0060 myBx = round(myTime / bunchCrossTimeDiff);
0061 if (!(myBx <= maxBX && myBx >= minBX))
0062 continue;
0063
0064 if (debug)
0065 LogDebug("CSCSegtoRPC") << "CSC \t \t This Segment is in Chamber id: " << CSCId << std::endl;
0066 if (debug)
0067 LogDebug("CSCSegtoRPC") << "CSC \t \t Number of segments in this CSC = " << CSCSegmentsCounter[CSCId]
0068 << std::endl;
0069 if (debug)
0070 LogDebug("CSCSegtoRPC")
0071 << "CSC \t \t Is the only one in this CSC? is not ind the ring 1? Are there more than 2 "
0072 "segments in the event?"
0073 << std::endl;
0074
0075 if (CSCSegmentsCounter[CSCId] == 1 && CSCId.ring() != 1 && allCSCSegments->size() >= 2) {
0076 if (debug)
0077 LogDebug("CSCSegtoRPC") << "CSC \t \t yes" << std::endl;
0078 int cscEndCap = CSCId.endcap();
0079 int cscStation = CSCId.station();
0080 int cscRing = CSCId.ring();
0081 int rpcRegion = 1;
0082 if (cscEndCap == 2)
0083 rpcRegion = -1;
0084 int rpcRing = cscRing;
0085 if (cscRing == 4)
0086 rpcRing = 1;
0087 int rpcStation = cscStation;
0088 int rpcSegment = CSCId.chamber();
0089
0090 LocalPoint segmentPosition = segment->localPosition();
0091 LocalVector segmentDirection = segment->localDirection();
0092 float dz = segmentDirection.z();
0093
0094 if (debug)
0095 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Information about the segment"
0096 << "RecHits =" << segment->nRecHits() << "Angle =" << acos(dz) * 180 / 3.1415926
0097 << std::endl;
0098
0099 if (debug)
0100 LogDebug("CSCSegtoRPC")
0101 << "CSC \t \t Is a good Segment? dim = 4, 4 <= nRecHits <= 10 Incident angle int range 45 < "
0102 << acos(dz) * 180 / 3.1415926 << " < 135? " << std::endl;
0103
0104 if ((segment->dimension() == 4) && (segment->nRecHits() <= 10 && segment->nRecHits() >= 4)) {
0105 if (debug)
0106 LogDebug("CSCSegtoRPC") << "CSC \t \t yes" << std::endl;
0107 if (debug)
0108 LogDebug("CSCSegtoRPC") << "CSC \t \t CSC Segment Dimension " << segment->dimension() << std::endl;
0109
0110 float Xo = segmentPosition.x();
0111 float Yo = segmentPosition.y();
0112 float Zo = segmentPosition.z();
0113 float dx = segmentDirection.x();
0114 float dy = segmentDirection.y();
0115 float dz = segmentDirection.z();
0116
0117 if (debug)
0118 LogDebug("CSCSegtoRPC") << "Creating the CSCIndex" << std::endl;
0119 CSCStationIndex theindex(rpcRegion, rpcStation, rpcRing, rpcSegment);
0120 if (debug)
0121 LogDebug("CSCSegtoRPC") << "Getting the Rolls for the given index" << std::endl;
0122 std::set<RPCDetId> rollsForThisCSC = cscMap->getRolls(theindex);
0123 if (debug)
0124 LogDebug("CSCSegtoRPC") << "CSC \t \t Getting chamber from Geometry" << std::endl;
0125 const CSCChamber* TheChamber = cscGeo->chamber(CSCId);
0126 if (debug)
0127 LogDebug("CSCSegtoRPC") << "CSC \t \t Getting ID from Chamber" << std::endl;
0128 const CSCDetId TheId = TheChamber->id();
0129
0130 if (debug)
0131 LogDebug("CSCSegtoRPC") << "CSC \t \t Number of rolls for this CSC = " << rollsForThisCSC.size()
0132 << std::endl;
0133
0134 if (debug)
0135 LogDebug("CSCSegtoRPC") << "CSC \t \t Printing The Id" << TheId << std::endl;
0136
0137 if (rpcRing != 1) {
0138
0139 assert(!rollsForThisCSC.empty());
0140
0141 if (debug)
0142 LogDebug("CSCSegtoRPC") << "CSC \t \t Loop over all the rolls asociated to this CSC" << std::endl;
0143 for (std::set<RPCDetId>::iterator iteraRoll = rollsForThisCSC.begin(); iteraRoll != rollsForThisCSC.end();
0144 iteraRoll++) {
0145 const RPCRoll* rollasociated = rpcGeo->roll(*iteraRoll);
0146 RPCDetId rpcId = rollasociated->id();
0147
0148 if (debug)
0149 LogDebug("CSCSegtoRPC") << "CSC \t \t \t We are in the roll getting the surface" << rpcId << std::endl;
0150 const BoundPlane& RPCSurface = rollasociated->surface();
0151
0152 if (debug)
0153 LogDebug("CSCSegtoRPC") << "CSC \t \t \t RollID: " << rpcId << std::endl;
0154
0155 if (debug)
0156 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Doing the extrapolation to this roll" << std::endl;
0157 if (debug)
0158 LogDebug("CSCSegtoRPC") << "CSC \t \t \t CSC Segment Direction in CSCLocal " << segmentDirection
0159 << std::endl;
0160 if (debug)
0161 LogDebug("CSCSegtoRPC") << "CSC \t \t \t CSC Segment Point in CSCLocal " << segmentPosition
0162 << std::endl;
0163
0164 GlobalPoint CenterPointRollGlobal = RPCSurface.toGlobal(LocalPoint(0, 0, 0));
0165 if (debug)
0166 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Center (0,0,0) of the Roll in Global" << CenterPointRollGlobal
0167 << std::endl;
0168 GlobalPoint CenterPointCSCGlobal = TheChamber->toGlobal(LocalPoint(0, 0, 0));
0169 if (debug)
0170 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Center (0,0,0) of the CSC in Global" << CenterPointCSCGlobal
0171 << std::endl;
0172 GlobalPoint segmentPositionInGlobal =
0173 TheChamber->toGlobal(segmentPosition);
0174 if (debug)
0175 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Segment Position in Global" << segmentPositionInGlobal
0176 << std::endl;
0177 LocalPoint CenterRollinCSCFrame = TheChamber->toLocal(CenterPointRollGlobal);
0178
0179 if (debug) {
0180 float rpcphi = 0;
0181 float cscphi = 0;
0182
0183 (CenterPointRollGlobal.barePhi() < 0) ? rpcphi = 2 * 3.141592 + CenterPointRollGlobal.barePhi()
0184 : rpcphi = CenterPointRollGlobal.barePhi();
0185
0186 (CenterPointCSCGlobal.barePhi() < 0) ? cscphi = 2 * 3.1415926536 + CenterPointCSCGlobal.barePhi()
0187 : cscphi = CenterPointCSCGlobal.barePhi();
0188
0189 float df = fabs(cscphi - rpcphi);
0190 float dr = fabs(CenterPointRollGlobal.perp() - CenterPointCSCGlobal.perp());
0191 float diffz = CenterPointRollGlobal.z() - CenterPointCSCGlobal.z();
0192 float dfg = df * 180. / 3.14159265;
0193
0194 if (debug)
0195 LogDebug("CSCSegtoRPC") << "CSC \t \t \t z of RPC=" << CenterPointRollGlobal.z() << "z of CSC"
0196 << CenterPointCSCGlobal.z() << " dfg=" << dfg << std::endl;
0197
0198 RPCGeomServ rpcsrv(rpcId);
0199
0200 if (dr > 200. || fabs(dz) > 55. || dfg > 1.) {
0201 if (debug)
0202 LogDebug("CSCSegtoRPC")
0203 << "\t \t \t CSC Station= " << CSCId.station() << " Ring= " << CSCId.ring()
0204 << " Chamber= " << CSCId.chamber() << " cscphi=" << cscphi * 180 / 3.14159265
0205 << "\t RPC Station= " << rpcId.station() << " ring= " << rpcId.ring() << " segment =-> "
0206 << rpcsrv.segment() << " rollphi=" << rpcphi * 180 / 3.14159265 << "\t dfg=" << dfg
0207 << " dz=" << diffz << " dr=" << dr << std::endl;
0208 }
0209 }
0210
0211 float D = CenterRollinCSCFrame.z();
0212
0213 float X = Xo + dx * D / dz;
0214 float Y = Yo + dy * D / dz;
0215 float Z = D;
0216
0217 const TrapezoidalStripTopology* top_ =
0218 dynamic_cast<const TrapezoidalStripTopology*>(&(rollasociated->topology()));
0219 LocalPoint xmin = top_->localPosition(0.);
0220 if (debug)
0221 LogDebug("CSCSegtoRPC") << "CSC \t \t \t xmin of this Roll " << xmin << "cm" << std::endl;
0222 LocalPoint xmax = top_->localPosition((float)rollasociated->nstrips());
0223 if (debug)
0224 LogDebug("CSCSegtoRPC") << "CSC \t \t \t xmax of this Roll " << xmax << "cm" << std::endl;
0225 float rsize = fabs(xmax.x() - xmin.x());
0226 if (debug)
0227 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Roll Size " << rsize << "cm" << std::endl;
0228 float stripl = top_->stripLength();
0229 float stripw = top_->pitch();
0230
0231 if (debug)
0232 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Strip Lenght " << stripl << "cm" << std::endl;
0233 if (debug)
0234 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Strip Width " << stripw << "cm" << std::endl;
0235
0236 if (debug)
0237 LogDebug("CSCSegtoRPC") << "CSC \t \t \t X Predicted in CSCLocal= " << X << "cm" << std::endl;
0238 if (debug)
0239 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Y Predicted in CSCLocal= " << Y << "cm" << std::endl;
0240 if (debug)
0241 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Z Predicted in CSCLocal= " << Z << "cm" << std::endl;
0242
0243 float extrapolatedDistance = sqrt((X - Xo) * (X - Xo) + (Y - Yo) * (Y - Yo) + (Z - Zo) * (Z - Zo));
0244
0245 if (debug)
0246 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Is the distance of extrapolation less than MaxD? ="
0247 << extrapolatedDistance << "cm"
0248 << " MaxD=" << MaxD << "cm" << std::endl;
0249
0250 if (extrapolatedDistance <= MaxD) {
0251 if (debug)
0252 LogDebug("CSCSegtoRPC") << "CSC \t \t \t yes" << std::endl;
0253
0254 GlobalPoint GlobalPointExtrapolated = TheChamber->toGlobal(LocalPoint(X, Y, Z));
0255 if (debug)
0256 LogDebug("CSCSegtoRPC")
0257 << "CSC \t \t \t Point ExtraPolated in Global" << GlobalPointExtrapolated << std::endl;
0258
0259 LocalPoint PointExtrapolatedRPCFrame = RPCSurface.toLocal(GlobalPointExtrapolated);
0260
0261 if (debug)
0262 LogDebug("CSCSegtoRPC")
0263 << "CSC \t \t \t Point Extrapolated in RPCLocal" << PointExtrapolatedRPCFrame << std::endl;
0264 if (debug)
0265 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Corner of the Roll = (" << rsize * eyr << "," << stripl * eyr
0266 << ")" << std::endl;
0267 if (debug)
0268 LogDebug("CSCSegtoRPC")
0269 << "CSC \t \t \t Info About the Point Extrapolated in X Abs ("
0270 << fabs(PointExtrapolatedRPCFrame.x()) << "," << fabs(PointExtrapolatedRPCFrame.y()) << ","
0271 << fabs(PointExtrapolatedRPCFrame.z()) << ")" << std::endl;
0272 if (debug)
0273 LogDebug("CSCSegtoRPC") << "CSC \t \t \t dz=" << fabs(PointExtrapolatedRPCFrame.z())
0274 << " dx=" << fabs(PointExtrapolatedRPCFrame.x())
0275 << " dy=" << fabs(PointExtrapolatedRPCFrame.y()) << std::endl;
0276
0277 if (debug)
0278 LogDebug("CSCSegtoRPC") << "CSC \t \t \t Does the extrapolation go inside this roll????" << std::endl;
0279
0280 if (fabs(PointExtrapolatedRPCFrame.z()) < 1. && fabs(PointExtrapolatedRPCFrame.x()) < rsize * eyr &&
0281 fabs(PointExtrapolatedRPCFrame.y()) < stripl * eyr) {
0282 if (debug)
0283 LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t yes" << std::endl;
0284 if (debug)
0285 LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Creating the RecHit" << std::endl;
0286
0287 RPCRecHit RPCPoint(rpcId, myBx, PointExtrapolatedRPCFrame);
0288 RPCPoint.setTimeAndError(myTime, myTimeErr);
0289
0290 if (debug)
0291 LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Clearing the vector" << std::endl;
0292 RPCPointVector.clear();
0293 if (debug)
0294 LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Pushing back" << std::endl;
0295 RPCPointVector.push_back(RPCPoint);
0296 if (debug)
0297 LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Putting the vector" << std::endl;
0298 _ThePoints->put(rpcId, RPCPointVector.begin(), RPCPointVector.end());
0299 }
0300 }
0301 }
0302 }
0303 }
0304 }
0305 }
0306 }
0307 return _ThePoints;
0308 }