File indexing completed on 2024-04-06 12:15:19
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #include <memory>
0021
0022
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0025
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030
0031 #include <memory>
0032 #include <DataFormats/MuonDetId/interface/RPCDetId.h>
0033 #include <Geometry/RPCGeometry/interface/RPCGeomServ.h>
0034 #include <Geometry/CommonDetUnit/interface/GeomDet.h>
0035 #include <Geometry/Records/interface/MuonGeometryRecord.h>
0036 #include <Geometry/CommonTopologies/interface/RectangularStripTopology.h>
0037 #include <Geometry/CommonTopologies/interface/TrapezoidalStripTopology.h>
0038 #include <Geometry/CSCGeometry/interface/CSCGeometry.h>
0039 #include <Geometry/RPCGeometry/interface/RPCGeometry.h>
0040 #include <Geometry/CSCGeometry/interface/CSCChamber.h>
0041 #include <DataFormats/MuonDetId/interface/CSCDetId.h>
0042
0043
0044
0045
0046 class CSCStationIndex {
0047 public:
0048 CSCStationIndex() : _region(0), _station(0), _ring(0), _chamber(0) {}
0049 CSCStationIndex(int region, int station, int ring, int chamber)
0050 : _region(region), _station(station), _ring(ring), _chamber(chamber) {}
0051 ~CSCStationIndex() {}
0052 int region() const { return _region; }
0053 int station() const { return _station; }
0054 int ring() const { return _ring; }
0055 int chamber() const { return _chamber; }
0056 bool operator<(const CSCStationIndex& cscind) const {
0057 if (cscind.region() != this->region())
0058 return cscind.region() < this->region();
0059 else if (cscind.station() != this->station())
0060 return cscind.station() < this->station();
0061 else if (cscind.ring() != this->ring())
0062 return cscind.ring() < this->ring();
0063 else if (cscind.chamber() != this->chamber())
0064 return cscind.chamber() < this->chamber();
0065 return false;
0066 }
0067
0068 private:
0069 int _region;
0070 int _station;
0071 int _ring;
0072 int _chamber;
0073 };
0074
0075 class RPCCSC : public edm::one::EDAnalyzer<> {
0076 public:
0077 explicit RPCCSC(const edm::ParameterSet&);
0078 ~RPCCSC() override;
0079
0080 void beginJob() override {}
0081 void analyze(edm::Event const& iEvent, edm::EventSetup const&) override;
0082 void endJob() override {}
0083
0084 private:
0085 const edm::ESGetToken<RPCGeometry, MuonGeometryRecord> tokRPC_;
0086 const edm::ESGetToken<CSCGeometry, MuonGeometryRecord> tokCSC_;
0087 std::map<CSCStationIndex, std::set<RPCDetId> > rollstoreCSC;
0088 };
0089
0090 RPCCSC::RPCCSC(const edm::ParameterSet& )
0091 : tokRPC_{esConsumes<RPCGeometry, MuonGeometryRecord>(edm::ESInputTag{})},
0092 tokCSC_{esConsumes<CSCGeometry, MuonGeometryRecord>(edm::ESInputTag{})} {}
0093
0094 RPCCSC::~RPCCSC() {}
0095
0096
0097 void RPCCSC::analyze(const edm::Event& , const edm::EventSetup& iSetup) {
0098 using namespace std;
0099 const RPCGeometry* rpcGeometry = &iSetup.getData(tokRPC_);
0100 const CSCGeometry* cscGeometry = &iSetup.getData(tokCSC_);
0101
0102 for (TrackingGeometry::DetContainer::const_iterator it = rpcGeometry->dets().begin(); it < rpcGeometry->dets().end();
0103 it++) {
0104 if (dynamic_cast<const RPCChamber*>(*it) != nullptr) {
0105 const RPCChamber* ch = dynamic_cast<const RPCChamber*>(*it);
0106 std::vector<const RPCRoll*> roles = (ch->rolls());
0107 for (std::vector<const RPCRoll*>::const_iterator r = roles.begin(); r != roles.end(); ++r) {
0108 RPCDetId rpcId = (*r)->id();
0109 int region = rpcId.region();
0110
0111 RPCGeomServ rpcsrv(rpcId);
0112 std::string nameRoll = rpcsrv.name();
0113
0114
0115 if (region != 0) {
0116
0117
0118
0119 int region = rpcId.region();
0120 int station = rpcId.station();
0121 int ring = rpcId.ring();
0122 int cscring = ring;
0123 int cscstation = station;
0124 RPCGeomServ rpcsrv(rpcId);
0125 int rpcsegment = rpcsrv.segment();
0126
0127 int cscchamber = rpcsegment;
0128 if ((station == 2 || station == 3) && ring == 3) {
0129 cscring = 2;
0130 }
0131
0132 CSCStationIndex ind(region, cscstation, cscring, cscchamber);
0133 std::set<RPCDetId> myrolls;
0134 if (rollstoreCSC.find(ind) != rollstoreCSC.end()) {
0135 myrolls = rollstoreCSC[ind];
0136 }
0137 myrolls.insert(rpcId);
0138 rollstoreCSC[ind] = myrolls;
0139 }
0140 }
0141 }
0142 }
0143 for (TrackingGeometry::DetContainer::const_iterator it = rpcGeometry->dets().begin(); it < rpcGeometry->dets().end();
0144 it++) {
0145 if (dynamic_cast<const RPCChamber*>(*it) != nullptr) {
0146 const RPCChamber* ch = dynamic_cast<const RPCChamber*>(*it);
0147 std::vector<const RPCRoll*> roles = (ch->rolls());
0148 for (std::vector<const RPCRoll*>::const_iterator r = roles.begin(); r != roles.end(); ++r) {
0149 RPCDetId rpcId = (*r)->id();
0150 int region = rpcId.region();
0151 if (region != 0 && (rpcId.ring() == 2 || rpcId.ring() == 3)) {
0152 int region = rpcId.region();
0153 int station = rpcId.station();
0154 int ring = rpcId.ring();
0155 int cscring = ring;
0156
0157 if ((station == 2 || station == 3) && ring == 3)
0158 cscring = 2;
0159
0160 int cscstation = station;
0161 RPCGeomServ rpcsrv(rpcId);
0162 int rpcsegment = rpcsrv.segment();
0163
0164 int cscchamber = rpcsegment + 1;
0165 if (cscchamber == 37)
0166 cscchamber = 1;
0167 CSCStationIndex ind(region, cscstation, cscring, cscchamber);
0168 std::set<RPCDetId> myrolls;
0169 if (rollstoreCSC.find(ind) != rollstoreCSC.end())
0170 myrolls = rollstoreCSC[ind];
0171 myrolls.insert(rpcId);
0172 rollstoreCSC[ind] = myrolls;
0173
0174 cscchamber = rpcsegment - 1;
0175 if (cscchamber == 0)
0176 cscchamber = 36;
0177 CSCStationIndex indDos(region, cscstation, cscring, cscchamber);
0178 std::set<RPCDetId> myrollsDos;
0179 if (rollstoreCSC.find(indDos) != rollstoreCSC.end())
0180 myrollsDos = rollstoreCSC[indDos];
0181 myrollsDos.insert(rpcId);
0182 rollstoreCSC[indDos] = myrollsDos;
0183 }
0184 }
0185 }
0186 }
0187
0188
0189
0190
0191 const CSCGeometry::ChamberContainer& cscChambers = cscGeometry->chambers();
0192 for (auto cscChamber : cscChambers) {
0193 CSCDetId CSCId = cscChamber->id();
0194
0195 int cscEndCap = CSCId.endcap();
0196 int cscStation = CSCId.station();
0197 int cscRing = CSCId.ring();
0198
0199 int rpcRegion = 1;
0200 if (cscEndCap == 2)
0201 rpcRegion = -1;
0202 int rpcRing = cscRing;
0203 if (cscRing == 4)
0204 rpcRing = 1;
0205 int rpcStation = cscStation;
0206 int rpcSegment = CSCId.chamber();
0207
0208
0209 const CSCChamber* TheChamber = cscGeometry->chamber(CSCId);
0210
0211
0212 std::set<RPCDetId> rollsForThisCSC = rollstoreCSC[CSCStationIndex(rpcRegion, rpcStation, rpcRing, rpcSegment)];
0213 if (CSCId.ring() != 1)
0214 edm::LogVerbatim("RPCGeometry") << "CSC for" << CSCId << " " << rollsForThisCSC.size() << " rolls.";
0215
0216 for (auto iteraRoll : rollsForThisCSC) {
0217 const RPCRoll* rollasociated = rpcGeometry->roll(iteraRoll);
0218 RPCDetId rpcId = rollasociated->id();
0219 RPCGeomServ rpcsrv(rpcId);
0220
0221 const BoundPlane& RPCSurface = rollasociated->surface();
0222
0223 GlobalPoint CenterPointRollGlobal = RPCSurface.toGlobal(LocalPoint(0, 0, 0));
0224 GlobalPoint CenterPointCSCGlobal = TheChamber->toGlobal(LocalPoint(0, 0, 0));
0225
0226
0227
0228 float rpcphi = 0;
0229 float cscphi = 0;
0230
0231 (CenterPointRollGlobal.barePhi() < 0) ? rpcphi = 2 * 3.141592 + CenterPointRollGlobal.barePhi()
0232 : rpcphi = CenterPointRollGlobal.barePhi();
0233
0234 (CenterPointCSCGlobal.barePhi() < 0) ? cscphi = 2 * 3.1415926536 + CenterPointCSCGlobal.barePhi()
0235 : cscphi = CenterPointCSCGlobal.barePhi();
0236
0237 float df = fabs(cscphi - rpcphi);
0238 float dr = fabs(CenterPointRollGlobal.perp() - CenterPointCSCGlobal.perp());
0239 float diffz = CenterPointRollGlobal.z() - CenterPointCSCGlobal.z();
0240 float dfg = df * 180. / 3.14159265;
0241
0242 edm::LogVerbatim("RPCGeometry") << "CSC \t " << rpcsrv.segment() << rpcsrv.name() << " dr=" << dr
0243 << " dz=" << diffz << " dfg=" << dfg;
0244
0245 bool print = false;
0246
0247 if ((dr > 200. || fabs(diffz) > 55. || dfg > 1.) && print) {
0248 edm::LogVerbatim("RPCGeometry") << "\t \t problem CSC Station= " << CSCId.station() << " Ring= " << CSCId.ring()
0249 << " Chamber= " << CSCId.chamber() << " cscphi=" << cscphi * 180 / 3.14159265
0250 << "\t RPC Station= " << rpcId.station() << " ring= " << rpcId.ring()
0251 << " segment =-> " << rpcsrv.segment()
0252 << " rollphi=" << rpcphi * 180 / 3.14159265 << "\t dfg=" << dfg
0253 << " dz=" << diffz << " dr=" << dr;
0254 }
0255 }
0256 }
0257 }
0258
0259
0260 DEFINE_FWK_MODULE(RPCCSC);