File indexing completed on 2024-04-06 12:07:46
0001 #include "FWCore/Framework/interface/Event.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004
0005 #include "Geometry/RPCGeometry/interface/RPCRoll.h"
0006 #include "Geometry/RPCGeometry/interface/RPCGeometry.h"
0007 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0008
0009 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0010 #include "DQMServices/Core/interface/MonitorElement.h"
0011
0012 #include "DataFormats/L1TMuon/interface/CPPFDigi.h"
0013 #include <string>
0014
0015 class L1TdeStage2CPPF : public DQMEDAnalyzer {
0016 public:
0017 L1TdeStage2CPPF(const edm::ParameterSet& ps);
0018 ~L1TdeStage2CPPF() override;
0019
0020 protected:
0021 void bookHistograms(DQMStore::IBooker&, const edm::Run&, const edm::EventSetup&) override;
0022 void analyze(const edm::Event&, const edm::EventSetup&) override;
0023
0024 private:
0025 int occupancy_value(int region_, int station_, int ring_);
0026 int bx_value(int region_, int emtfsector_);
0027 int GetSubsector(int emtfsector_, int lsubsector_);
0028
0029 edm::EDGetTokenT<l1t::CPPFDigiCollection> dataToken;
0030 edm::EDGetTokenT<l1t::CPPFDigiCollection> emulToken;
0031 std::string monitorDir;
0032 bool verbose;
0033
0034 MonitorElement* h2_Matching_SameKey_OnPhi_phi_Ce_phi_Cu_bx;
0035 MonitorElement* h2_Matching_SameKey_OnPhi_theta_Ce_theta_Cu_bx;
0036 MonitorElement* h2_Matching_SameKey_OnPhi_zone_Ce_zone_Cu_bx;
0037 MonitorElement* h2_Matching_SameKey_OnPhi_ID_Ce_ID_Cu_bx;
0038 MonitorElement* h2_Matching_SameKey_OnPhi_ID_Ce_roll_Ce_bx;
0039
0040 MonitorElement* h2_Matching_SameKey_OffPhi_phi_Ce_phi_Cu_bx;
0041 MonitorElement* h2_Matching_SameKey_OffPhi_theta_Ce_theta_Cu_bx;
0042 MonitorElement* h2_Matching_SameKey_OffPhi_zone_Ce_zone_Cu_bx;
0043 MonitorElement* h2_Matching_SameKey_OffPhi_ID_Ce_ID_Cu_bx;
0044 MonitorElement* h2_Matching_SameKey_OffPhi_ID_Ce_roll_Ce_bx;
0045
0046 MonitorElement* h2_Matching_SameKey_OnTheta_phi_Ce_phi_Cu_bx;
0047 MonitorElement* h2_Matching_SameKey_OnTheta_theta_Ce_theta_Cu_bx;
0048 MonitorElement* h2_Matching_SameKey_OnTheta_zone_Ce_zone_Cu_bx;
0049
0050 MonitorElement* h2_Matching_SameKey_OffTheta_phi_Ce_phi_Cu_bx;
0051 MonitorElement* h2_Matching_SameKey_OffTheta_theta_Ce_theta_Cu_bx;
0052 MonitorElement* h2_Matching_SameKey_OffTheta_zone_Ce_zone_Cu_bx;
0053
0054 MonitorElement* h1_Matching_SameKey_bx_Summary;
0055 };
0056
0057 L1TdeStage2CPPF::L1TdeStage2CPPF(const edm::ParameterSet& ps)
0058 : dataToken(consumes<l1t::CPPFDigiCollection>(ps.getParameter<edm::InputTag>("dataSource"))),
0059 emulToken(consumes<l1t::CPPFDigiCollection>(ps.getParameter<edm::InputTag>("emulSource"))),
0060 monitorDir(ps.getUntrackedParameter<std::string>("monitorDir", "")),
0061 verbose(ps.getUntrackedParameter<bool>("verbose", false)) {}
0062
0063 L1TdeStage2CPPF::~L1TdeStage2CPPF() {}
0064
0065 void L1TdeStage2CPPF::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run&, const edm::EventSetup&) {
0066 ibooker.setCurrentFolder(monitorDir);
0067
0068 h2_Matching_SameKey_OnPhi_phi_Ce_phi_Cu_bx =
0069 ibooker.book2D("h2_Matching_SameKey_OnPhi_phi_Ce_phi_Cu_bx",
0070 "Matching && Same SubSector && OnPhi ; Emulator #phi ; Unpacker #phi",
0071 62,
0072 0.0,
0073 1240.,
0074 62,
0075 0.0,
0076 1240.);
0077 h2_Matching_SameKey_OnPhi_theta_Ce_theta_Cu_bx =
0078 ibooker.book2D("h2_Matching_SameKey_OnPhi_theta_Ce_theta_Cu_bx",
0079 "Matching && Same SubSector && bx==0 && OnPhi ; Emulator #theta ; Unpacker #theta ",
0080 32,
0081 0,
0082 32.,
0083 32,
0084 0,
0085 32.);
0086 h2_Matching_SameKey_OnPhi_zone_Ce_zone_Cu_bx =
0087 ibooker.book2D("h2_Matching_SameKey_OnPhi_zone_Ce_zone_Cu_bx",
0088 "Matching && Same SubSector && bx==0 && OnPhi ;Emulator Zone ;Unpacker Zone ",
0089 15,
0090 0,
0091 15,
0092 15,
0093 0,
0094 15);
0095 h2_Matching_SameKey_OnPhi_ID_Ce_ID_Cu_bx =
0096 ibooker.book2D("h2_Matching_SameKey_OnPhi_ID_Ce_ID_Cu_bx",
0097 "Matching && Same SubSector && bx==0 && OnPhi ; Emulator Chamber ID;Unpacker Chamber ID ",
0098 38,
0099 0,
0100 38,
0101 38,
0102 0,
0103 38);
0104 h2_Matching_SameKey_OnPhi_ID_Ce_roll_Ce_bx =
0105 ibooker.book2D("h2_Matching_SameKey_OnPhi_ID_Ce_roll_Ce_bx",
0106 "Matching && Same SubSector && bx==0 && OnPhi ;Emulator Chamber ID ;Emulator Roll ",
0107 38,
0108 0,
0109 38,
0110 4,
0111 0,
0112 4);
0113
0114 h2_Matching_SameKey_OffPhi_phi_Ce_phi_Cu_bx =
0115 ibooker.book2D("h2_Matching_SameKey_OffPhi_phi_Ce_phi_Cu_bx",
0116 "Matching && Same SubSector && OffPhi ; Emulator #phi ; Unpacker #phi",
0117 62,
0118 0.0,
0119 1240.,
0120 62,
0121 0.0,
0122 1240.);
0123 h2_Matching_SameKey_OffPhi_theta_Ce_theta_Cu_bx =
0124 ibooker.book2D("h2_Matching_SameKey_OffPhi_theta_Ce_theta_Cu_bx",
0125 "Matching && Same SubSector && bx==0 && OffPhi ; Emulator #theta ; Unpacker #theta ",
0126 32,
0127 0,
0128 32.,
0129 32,
0130 0,
0131 32.);
0132 h2_Matching_SameKey_OffPhi_zone_Ce_zone_Cu_bx =
0133 ibooker.book2D("h2_Matching_SameKey_OffPhi_zone_Ce_zone_Cu_bx",
0134 "Matching && Same SubSector && bx==0 && OffPhi ;Emulator Zone ;Unpacker Zone ",
0135 15,
0136 0,
0137 15,
0138 15,
0139 0,
0140 15);
0141 h2_Matching_SameKey_OffPhi_ID_Ce_ID_Cu_bx =
0142 ibooker.book2D("h2_Matching_SameKey_OffPhi_ID_Ce_ID_Cu_bx",
0143 "Matching && Same SubSector && bx==0 && OffPhi ; Emulator Chamber ID;Unpacker Chamber ID ",
0144 38,
0145 0,
0146 38,
0147 38,
0148 0,
0149 38);
0150 h2_Matching_SameKey_OffPhi_ID_Ce_roll_Ce_bx =
0151 ibooker.book2D("h2_Matching_SameKey_OffPhi_ID_Ce_roll_Ce_bx",
0152 "Matching && Same SubSector && bx==0 && OffPhi ;Emulator Chamber ID ;Emulator Roll ",
0153 38,
0154 0,
0155 38,
0156 4,
0157 0,
0158 4);
0159
0160 h2_Matching_SameKey_OnTheta_phi_Ce_phi_Cu_bx =
0161 ibooker.book2D("h2_Matching_SameKey_OnTheta_phi_Ce_phi_Cu_bx",
0162 "Matching && Same SubSector && OnTheta ; Emulator #phi ; Unpacker #phi",
0163 62,
0164 0.0,
0165 1240.,
0166 62,
0167 0.0,
0168 1240.);
0169 h2_Matching_SameKey_OnTheta_theta_Ce_theta_Cu_bx =
0170 ibooker.book2D("h2_Matching_SameKey_OnTheta_theta_Ce_theta_Cu_bx",
0171 "Matching && Same SubSector && bx==0 && OnTheta ; Emulator #theta ; Unpacker #theta ",
0172 32,
0173 0,
0174 32.,
0175 32,
0176 0,
0177 32.);
0178 h2_Matching_SameKey_OnTheta_zone_Ce_zone_Cu_bx =
0179 ibooker.book2D("h2_Matching_SameKey_OnTheta_zone_Ce_zone_Cu_bx",
0180 "Matching && Same SubSector && bx==0 && OnTheta ;Emulator Zone ;Unpacker Zone ",
0181 15,
0182 0,
0183 15,
0184 15,
0185 0,
0186 15);
0187
0188 h2_Matching_SameKey_OffTheta_phi_Ce_phi_Cu_bx =
0189 ibooker.book2D("h2_Matching_SameKey_OffTheta_phi_Ce_phi_Cu_bx",
0190 "Matching && Same SubSector && OffTheta ; Emulator #phi ; Unpacker #phi",
0191 62,
0192 0.0,
0193 1240.,
0194 62,
0195 0.0,
0196 1240.);
0197 h2_Matching_SameKey_OffTheta_theta_Ce_theta_Cu_bx =
0198 ibooker.book2D("h2_Matching_SameKey_OffTheta_theta_Ce_theta_Cu_bx",
0199 "Matching && Same SubSector && bx==0 && OffTheta ; Emulator #theta ; Unpacker #theta ",
0200 32,
0201 0,
0202 32.,
0203 32,
0204 0,
0205 32.);
0206 h2_Matching_SameKey_OffTheta_zone_Ce_zone_Cu_bx =
0207 ibooker.book2D("h2_Matching_SameKey_OffTheta_zone_Ce_zone_Cu_bx",
0208 "Matching && Same SubSector && bx==0 && OffTheta ;Emulator Zone ;Unpacker Zone ",
0209 15,
0210 0,
0211 15,
0212 15,
0213 0,
0214 15);
0215
0216 h1_Matching_SameKey_bx_Summary =
0217 ibooker.book1D("h1_Matching_SameKey_bx_Summary",
0218 "cppf data-emul mismatch fraction summary; ; Fraction events with mismatch",
0219 2,
0220 1,
0221 3);
0222 h1_Matching_SameKey_bx_Summary->setBinLabel(1, "off/on-phi");
0223 h1_Matching_SameKey_bx_Summary->setBinLabel(2, "off/on-theta");
0224 h1_Matching_SameKey_bx_Summary->setAxisRange(0, 0.01, 2);
0225 }
0226
0227 void L1TdeStage2CPPF::analyze(const edm::Event& e, const edm::EventSetup& c) {
0228 if (verbose)
0229 edm::LogInfo("L1TdeStage2CPPF") << "L1TdeStage2CPPF: analyze...";
0230
0231 edm::Handle<l1t::CPPFDigiCollection> dataCPPFs;
0232 e.getByToken(dataToken, dataCPPFs);
0233
0234 edm::Handle<l1t::CPPFDigiCollection> emulCPPFs;
0235 e.getByToken(emulToken, emulCPPFs);
0236
0237 std::unordered_map<int, int> _nHit_Ce;
0238 std::unordered_map<int, int> _nHit_Cu;
0239 std::unordered_map<int, std::vector<int>> _phi_Ce;
0240 std::unordered_map<int, std::vector<int>> _phi_Cu;
0241 std::unordered_map<int, std::vector<int>> _phi_glob_Ce;
0242 std::unordered_map<int, std::vector<int>> _phi_glob_Cu;
0243 std::unordered_map<int, std::vector<int>> _theta_Ce;
0244 std::unordered_map<int, std::vector<int>> _theta_Cu;
0245 std::unordered_map<int, std::vector<int>> _theta_glob_Ce;
0246 std::unordered_map<int, std::vector<int>> _theta_glob_Cu;
0247 std::unordered_map<int, std::vector<int>> _roll_Ce;
0248 std::unordered_map<int, std::vector<int>> _roll_Cu;
0249 std::unordered_map<int, std::vector<int>> _zone_Ce;
0250 std::unordered_map<int, std::vector<int>> _zone_Cu;
0251 std::unordered_map<int, std::vector<int>> _ID_Ce;
0252 std::unordered_map<int, std::vector<int>> _ID_Cu;
0253 std::unordered_map<int, std::vector<int>> _emtfSubsector_Ce;
0254 std::unordered_map<int, std::vector<int>> _emtfSubsector_Cu;
0255 std::unordered_map<int, std::vector<int>> _emtfSector_Ce;
0256 std::unordered_map<int, std::vector<int>> _emtfSector_Cu;
0257 std::unordered_map<int, std::vector<int>> _bx_Ce;
0258 std::unordered_map<int, std::vector<int>> _bx_Cu;
0259 std::unordered_map<int, std::vector<int>> _cluster_size_Ce;
0260 std::unordered_map<int, std::vector<int>> _cluster_size_Cu;
0261
0262 for (auto& cppf_digis : *emulCPPFs) {
0263 RPCDetId rpcIdCe = (int)cppf_digis.rpcId();
0264 int regionCe = (int)rpcIdCe.region();
0265 int stationCe = (int)rpcIdCe.station();
0266 int sectorCe = (int)rpcIdCe.sector();
0267 int subsectorCe = (int)rpcIdCe.subsector();
0268 int ringCe = (int)rpcIdCe.ring();
0269 int rollCe = (int)(rpcIdCe.roll());
0270 int phiIntCe = (int)cppf_digis.phi_int();
0271 int thetaIntCe = (int)cppf_digis.theta_int();
0272 int phiGlobalCe = (int)cppf_digis.phi_glob();
0273 int thetaGlobalCe = (int)cppf_digis.theta_glob();
0274 int cluster_sizeCe = (int)cppf_digis.cluster_size();
0275 int bxCe = cppf_digis.bx();
0276 int emtfSectorCe = (int)cppf_digis.emtf_sector();
0277 int emtfSubsectorCe = GetSubsector(emtfSectorCe, subsectorCe);
0278 int fillOccupancyCe = occupancy_value(regionCe, stationCe, ringCe);
0279
0280 int nsubCe = 6;
0281 (ringCe == 1 && stationCe > 1) ? nsubCe = 3 : nsubCe = 6;
0282 int chamberIDCe = subsectorCe + nsubCe * (sectorCe - 1);
0283
0284 std::ostringstream oss;
0285 oss << regionCe << stationCe << ringCe << sectorCe << subsectorCe << emtfSectorCe << emtfSubsectorCe;
0286 std::istringstream iss(oss.str());
0287 int unique_id;
0288 iss >> unique_id;
0289
0290 if (_nHit_Ce.find(unique_id) == _nHit_Ce.end()) {
0291 _nHit_Ce.insert({unique_id, 1});
0292 _phi_Ce[unique_id].push_back(phiIntCe);
0293 _phi_glob_Ce[unique_id].push_back(phiGlobalCe);
0294 _theta_Ce[unique_id].push_back(thetaIntCe);
0295 _theta_glob_Ce[unique_id].push_back(thetaGlobalCe);
0296 _roll_Ce[unique_id].push_back(rollCe);
0297 _ID_Ce[unique_id].push_back(chamberIDCe);
0298 _zone_Ce[unique_id].push_back(fillOccupancyCe);
0299 _emtfSubsector_Ce[unique_id].push_back(emtfSubsectorCe);
0300 _emtfSector_Ce[unique_id].push_back(emtfSectorCe);
0301 _bx_Ce[unique_id].push_back(bxCe);
0302 _cluster_size_Ce[unique_id].push_back(cluster_sizeCe);
0303 } else {
0304 _nHit_Ce.at(unique_id) += 1;
0305 _phi_Ce[unique_id].push_back(phiIntCe);
0306 _phi_glob_Ce[unique_id].push_back(phiGlobalCe);
0307 _theta_Ce[unique_id].push_back(thetaIntCe);
0308 _theta_glob_Ce[unique_id].push_back(thetaGlobalCe);
0309 _roll_Ce[unique_id].push_back(rollCe);
0310 _ID_Ce[unique_id].push_back(chamberIDCe);
0311 _zone_Ce[unique_id].push_back(fillOccupancyCe);
0312 _emtfSubsector_Ce[unique_id].push_back(emtfSubsectorCe);
0313 _emtfSector_Ce[unique_id].push_back(emtfSectorCe);
0314 _bx_Ce[unique_id].push_back(bxCe);
0315 _cluster_size_Ce[unique_id].push_back(cluster_sizeCe);
0316 }
0317
0318 }
0319
0320 for (auto& cppf_digis2 : *dataCPPFs) {
0321 RPCDetId rpcIdCu = cppf_digis2.rpcId();
0322 int regionCu = (int)rpcIdCu.region();
0323 int stationCu = (int)rpcIdCu.station();
0324 int sectorCu = (int)rpcIdCu.sector();
0325 int subsectorCu = (int)rpcIdCu.subsector();
0326 int ringCu = (int)rpcIdCu.ring();
0327 int rollCu = (int)(rpcIdCu.roll());
0328 int phiIntCu = (int)cppf_digis2.phi_int();
0329 int thetaIntCu = (int)cppf_digis2.theta_int();
0330 int phiGlobalCu = (int)cppf_digis2.phi_glob();
0331 int thetaGlobalCu = (int)cppf_digis2.theta_glob();
0332 int cluster_sizeCu = (int)cppf_digis2.cluster_size();
0333 int bxCu = (int)cppf_digis2.bx();
0334 int emtfSectorCu = (int)cppf_digis2.emtf_sector();
0335 int emtfSubsectorCu = GetSubsector(emtfSectorCu, subsectorCu);
0336 int fillOccupancyCu = occupancy_value(regionCu, stationCu, ringCu);
0337
0338 int nsubCu = 6;
0339 (ringCu == 1 && stationCu > 1) ? nsubCu = 3 : nsubCu = 6;
0340 int chamberIDCu = subsectorCu + nsubCu * (sectorCu - 1);
0341
0342 std::ostringstream oss2;
0343 oss2 << regionCu << stationCu << ringCu << sectorCu << subsectorCu << emtfSectorCu << emtfSubsectorCu;
0344 std::istringstream iss2(oss2.str());
0345 int unique_id;
0346 iss2 >> unique_id;
0347
0348 if (_nHit_Cu.find(unique_id) == _nHit_Cu.end()) {
0349 _nHit_Cu.insert({unique_id, 1});
0350 _phi_Cu[unique_id].push_back(phiIntCu);
0351 _theta_Cu[unique_id].push_back(thetaIntCu);
0352 _phi_glob_Cu[unique_id].push_back(phiGlobalCu);
0353 _theta_glob_Cu[unique_id].push_back(thetaGlobalCu);
0354 _ID_Cu[unique_id].push_back(chamberIDCu);
0355 _zone_Cu[unique_id].push_back(fillOccupancyCu);
0356 _roll_Cu[unique_id].push_back(rollCu);
0357 _emtfSubsector_Cu[unique_id].push_back(emtfSubsectorCu);
0358 _emtfSector_Cu[unique_id].push_back(emtfSectorCu);
0359 _bx_Cu[unique_id].push_back(bxCu);
0360 _cluster_size_Cu[unique_id].push_back(cluster_sizeCu);
0361 } else {
0362 _nHit_Cu.at(unique_id) += 1;
0363 _phi_Cu[unique_id].push_back(phiIntCu);
0364 _theta_Cu[unique_id].push_back(thetaIntCu);
0365 _phi_glob_Cu[unique_id].push_back(phiGlobalCu);
0366 _theta_glob_Cu[unique_id].push_back(thetaGlobalCu);
0367 _ID_Cu[unique_id].push_back(chamberIDCu);
0368 _zone_Cu[unique_id].push_back(fillOccupancyCu);
0369 _roll_Cu[unique_id].push_back(rollCu);
0370 _emtfSubsector_Cu[unique_id].push_back(emtfSubsectorCu);
0371 _emtfSector_Cu[unique_id].push_back(emtfSectorCu);
0372 _bx_Cu[unique_id].push_back(bxCu);
0373 _cluster_size_Cu[unique_id].push_back(cluster_sizeCu);
0374 }
0375 }
0376
0377 for (auto const& Ce : _nHit_Ce) {
0378 int key_Ce = Ce.first;
0379 int nHit_Ce = Ce.second;
0380
0381 for (auto const& Cu : _nHit_Cu) {
0382 int key_Cu = Cu.first;
0383 int nHit_Cu = Cu.second;
0384
0385 if (key_Ce != key_Cu)
0386 continue;
0387 if (nHit_Ce != nHit_Cu)
0388 continue;
0389
0390 for (int vecSize = 0; vecSize < nHit_Cu; ++vecSize) {
0391 if (_bx_Cu.at(key_Cu)[vecSize] != _bx_Ce.at(key_Ce)[vecSize])
0392 continue;
0393
0394 bool OnPhi_Matching = false;
0395 int index_Ce = vecSize;
0396 int index_Cu = vecSize;
0397 for (int i = 0; i < nHit_Ce; ++i) {
0398 if (_phi_Ce.at(key_Ce)[i] == _phi_Cu.at(key_Cu)[vecSize]) {
0399 OnPhi_Matching = true;
0400 index_Cu = vecSize;
0401 index_Ce = i;
0402 }
0403 }
0404 if (OnPhi_Matching) {
0405 h2_Matching_SameKey_OnPhi_phi_Ce_phi_Cu_bx->Fill(_phi_Ce.at(key_Ce)[index_Ce], _phi_Cu.at(key_Cu)[index_Cu]);
0406 h2_Matching_SameKey_OnPhi_theta_Ce_theta_Cu_bx->Fill(_theta_Ce.at(key_Ce)[index_Ce],
0407 _theta_Cu.at(key_Cu)[index_Cu]);
0408 h2_Matching_SameKey_OnPhi_zone_Ce_zone_Cu_bx->Fill(_zone_Ce.at(key_Ce)[index_Ce],
0409 _zone_Cu.at(key_Cu)[index_Cu]);
0410 h2_Matching_SameKey_OnPhi_ID_Ce_ID_Cu_bx->Fill(_ID_Ce.at(key_Ce)[index_Ce], _ID_Cu.at(key_Cu)[index_Cu]);
0411 h2_Matching_SameKey_OnPhi_ID_Ce_roll_Ce_bx->Fill(_ID_Ce.at(key_Ce)[index_Ce], _roll_Ce.at(key_Ce)[index_Ce]);
0412 } else {
0413 h2_Matching_SameKey_OffPhi_phi_Ce_phi_Cu_bx->Fill(_phi_Ce.at(key_Ce)[index_Ce], _phi_Cu.at(key_Cu)[index_Cu]);
0414 h2_Matching_SameKey_OffPhi_theta_Ce_theta_Cu_bx->Fill(_theta_Ce.at(key_Ce)[index_Ce],
0415 _theta_Cu.at(key_Cu)[index_Cu]);
0416 h2_Matching_SameKey_OffPhi_zone_Ce_zone_Cu_bx->Fill(_zone_Ce.at(key_Ce)[index_Ce],
0417 _zone_Cu.at(key_Cu)[index_Cu]);
0418 h2_Matching_SameKey_OffPhi_ID_Ce_ID_Cu_bx->Fill(_ID_Ce.at(key_Ce)[index_Ce], _ID_Cu.at(key_Cu)[index_Cu]);
0419 h2_Matching_SameKey_OffPhi_ID_Ce_roll_Ce_bx->Fill(_ID_Ce.at(key_Ce)[index_Ce], _roll_Ce.at(key_Ce)[index_Ce]);
0420 }
0421
0422 bool OnTheta_Matching = false;
0423 for (int i = 0; i < nHit_Ce; ++i) {
0424 if (_theta_Ce.at(key_Ce)[i] == _theta_Cu.at(key_Cu)[index_Cu]) {
0425 OnTheta_Matching = true;
0426 index_Cu = vecSize;
0427 index_Ce = i;
0428 }
0429 }
0430 if (OnTheta_Matching) {
0431 h2_Matching_SameKey_OnTheta_phi_Ce_phi_Cu_bx->Fill(_phi_Ce.at(key_Ce)[index_Ce],
0432 _phi_Cu.at(key_Cu)[index_Cu]);
0433 h2_Matching_SameKey_OnTheta_theta_Ce_theta_Cu_bx->Fill(_theta_Ce.at(key_Ce)[index_Ce],
0434 _theta_Cu.at(key_Cu)[index_Cu]);
0435 h2_Matching_SameKey_OnTheta_zone_Ce_zone_Cu_bx->Fill(_zone_Ce.at(key_Ce)[index_Ce],
0436 _zone_Cu.at(key_Cu)[index_Cu]);
0437 } else {
0438 h2_Matching_SameKey_OffTheta_phi_Ce_phi_Cu_bx->Fill(_phi_Ce.at(key_Ce)[index_Ce],
0439 _phi_Cu.at(key_Cu)[index_Cu]);
0440 h2_Matching_SameKey_OffTheta_theta_Ce_theta_Cu_bx->Fill(_theta_Ce.at(key_Ce)[index_Ce],
0441 _theta_Cu.at(key_Cu)[index_Cu]);
0442 h2_Matching_SameKey_OffTheta_zone_Ce_zone_Cu_bx->Fill(_zone_Ce.at(key_Ce)[index_Ce],
0443 _zone_Cu.at(key_Cu)[index_Cu]);
0444 }
0445 }
0446
0447 double off_phi, on_phi, off_theta, on_theta;
0448 on_phi = h2_Matching_SameKey_OnPhi_phi_Ce_phi_Cu_bx->getEntries();
0449 off_phi = h2_Matching_SameKey_OffPhi_phi_Ce_phi_Cu_bx->getEntries();
0450 on_theta = h2_Matching_SameKey_OnTheta_theta_Ce_theta_Cu_bx->getEntries();
0451 off_theta = h2_Matching_SameKey_OffTheta_theta_Ce_theta_Cu_bx->getEntries();
0452 if (on_phi == 0)
0453 on_phi = 0.0001;
0454 if (on_theta == 0)
0455 on_theta = 0.0001;
0456 h1_Matching_SameKey_bx_Summary->setBinContent(1, off_phi / on_phi);
0457 h1_Matching_SameKey_bx_Summary->setBinContent(2, off_theta / on_theta);
0458 }
0459 }
0460 }
0461
0462 int L1TdeStage2CPPF::GetSubsector(int emtfsector_, int lsubsector_) {
0463 const int nsectors = 6;
0464 int gsubsector = 0;
0465 if ((emtfsector_ != -99) and (lsubsector_ != 0)) {
0466 gsubsector = (emtfsector_ - 1) * nsectors + lsubsector_;
0467 }
0468 return gsubsector;
0469 }
0470
0471 int L1TdeStage2CPPF::occupancy_value(int region_, int station_, int ring_) {
0472 int fill_val = 0;
0473 if (region_ == -1) {
0474 if ((station_ == 4) && (ring_ == 3))
0475 fill_val = 1;
0476 else if ((station_ == 4) && (ring_ == 2))
0477 fill_val = 2;
0478 else if ((station_ == 3) && (ring_ == 3))
0479 fill_val = 3;
0480 else if ((station_ == 3) && (ring_ == 2))
0481 fill_val = 4;
0482 else if ((station_ == 2) && (ring_ == 2))
0483 fill_val = 5;
0484 else if ((station_ == 1) && (ring_ == 2))
0485 fill_val = 6;
0486
0487 } else if (region_ == +1) {
0488 if ((station_ == 1) && (ring_ == 2))
0489 fill_val = 7;
0490 else if ((station_ == 2) && (ring_ == 2))
0491 fill_val = 8;
0492 else if ((station_ == 3) && (ring_ == 2))
0493 fill_val = 9;
0494 else if ((station_ == 3) && (ring_ == 3))
0495 fill_val = 10;
0496 else if ((station_ == 4) && (ring_ == 2))
0497 fill_val = 11;
0498 else if ((station_ == 4) && (ring_ == 3))
0499 fill_val = 12;
0500 }
0501 return fill_val;
0502 }
0503
0504 int L1TdeStage2CPPF::bx_value(int region_, int emtfsector_) {
0505 int fill_val = 0;
0506
0507 if (region_ == -1) {
0508 fill_val = 7 - emtfsector_;
0509 } else if (region_ == +1) {
0510 fill_val = 6 + emtfsector_;
0511 }
0512 return fill_val;
0513 }
0514 #include "FWCore/Framework/interface/MakerMacros.h"
0515 DEFINE_FWK_MODULE(L1TdeStage2CPPF);