Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:55

0001 #include <string>
0002 #include <vector>
0003 #include <iostream>
0004 #include <map>
0005 
0006 #include "DQM/L1TMonitor/interface/L1TStage2EMTF.h"
0007 
0008 L1TStage2EMTF::L1TStage2EMTF(const edm::ParameterSet& ps)
0009     : daqToken(consumes<l1t::EMTFDaqOutCollection>(ps.getParameter<edm::InputTag>("emtfSource"))),
0010       hitToken(consumes<l1t::EMTFHitCollection>(ps.getParameter<edm::InputTag>("emtfSource"))),
0011       trackToken(consumes<l1t::EMTFTrackCollection>(ps.getParameter<edm::InputTag>("emtfSource"))),
0012       muonToken(consumes<l1t::RegionalMuonCandBxCollection>(ps.getParameter<edm::InputTag>("emtfSource"))),
0013       monitorDir(ps.getUntrackedParameter<std::string>("monitorDir", "")),
0014       verbose(ps.getUntrackedParameter<bool>("verbose", false)) {}
0015 
0016 void L1TStage2EMTF::bookHistograms(DQMStore::IBooker& ibooker, const edm::Run&, const edm::EventSetup&) {
0017   // Monitor Dir
0018   ibooker.setCurrentFolder(monitorDir);
0019 
0020   const std::array<std::string, 6> binNamesErrors{
0021       {"Corruptions", "Synch. Err.", "Synch. Mod.", "BX Mismatch", "Time Misalign", "FMM != Ready"}};
0022 
0023   // DAQ Output Monitor Elements
0024   emtfErrors = ibooker.book1D("emtfErrors", "EMTF Errors", 6, 0, 6);
0025   emtfErrors->setAxisTitle("Error Type (Corruptions Not Implemented)", 1);
0026   emtfErrors->setAxisTitle("Number of Errors", 2);
0027   for (unsigned int bin = 0; bin < binNamesErrors.size(); ++bin) {
0028     emtfErrors->setBinLabel(bin + 1, binNamesErrors[bin], 1);
0029   }
0030 
0031   // CSC LCT Monitor Elements
0032   int nChambs, nWires, nStrips;  // Number of chambers, wiregroups, and halfstrips in each station/ring pair
0033   std::string name, label;
0034   const std::array<std::string, 10> suffix_name{{"42", "41", "32", "31", "22", "21", "13", "12", "11b", "11a"}};
0035   const std::array<std::string, 10> suffix_label{
0036       {"4/2", "4/1", "3/2", "3/1", " 2/2", "2/1", "1/3", "1/2", "1/1b", "1/1a"}};
0037   const std::array<std::string, 12> binNames{
0038       {"ME-N", "ME-4", "ME-3", "ME-2", "ME-1b", "ME-1a", "ME+1a", "ME+1b", "ME+2", "ME+3", "ME+4", "ME+N"}};
0039 
0040   cscLCTBX = ibooker.book2D("cscLCTBX", "CSC LCT BX", 7, -3, 4, 20, 0, 20);
0041   cscLCTBX->setAxisTitle("BX", 1);
0042   for (int xbin = 1, xbin_label = -3; xbin <= 7; ++xbin, ++xbin_label) {
0043     cscLCTBX->setBinLabel(xbin, std::to_string(xbin_label), 1);
0044   }
0045   for (int ybin = 1; ybin <= 10; ++ybin) {
0046     cscLCTBX->setBinLabel(ybin, "ME-" + suffix_label[ybin - 1], 2);
0047     cscLCTBX->setBinLabel(21 - ybin, "ME+" + suffix_label[ybin - 1], 2);
0048   }
0049 
0050   cscLCTOccupancy = ibooker.book2D("cscLCTOccupancy", "CSC Chamber Occupancy", 54, 1, 55, 12, -6, 6);
0051   cscLCTOccupancy->setAxisTitle("Sector (CSCID 1-9 Unlabelled)", 1);
0052   for (int xbin = 1; xbin < 7; ++xbin) {
0053     cscLCTOccupancy->setBinLabel(xbin * 9 - 8, std::to_string(xbin), 1);
0054   }
0055   for (unsigned int ybin = 0; ybin < binNames.size(); ++ybin) {
0056     cscLCTOccupancy->setBinLabel(ybin + 1, binNames[ybin], 2);
0057   }
0058 
0059   //cscOccupancy designed to match the cscDQM plot
0060   cscDQMOccupancy = ibooker.book2D("cscDQMOccupancy", "CSC Chamber Occupancy", 42, 1, 43, 20, 0, 20);
0061   cscDQMOccupancy->setAxisTitle("10#circ Chamber (Ni = Neighbor of Sector i)", 1);
0062   int count = 0;
0063   for (int xbin = 1; xbin < 43; ++xbin) {
0064     cscDQMOccupancy->setBinLabel(xbin, std::to_string(xbin - count), 1);
0065     if (xbin == 2 || xbin == 9 || xbin == 16 || xbin == 23 || xbin == 30 || xbin == 37) {
0066       ++xbin;
0067       ++count;
0068       cscDQMOccupancy->setBinLabel(xbin, "N" + std::to_string(count), 1);
0069     }
0070   }
0071   for (int ybin = 1; ybin <= 10; ++ybin) {
0072     cscDQMOccupancy->setBinLabel(ybin, "ME-" + suffix_label[ybin - 1], 2);
0073     cscDQMOccupancy->setBinLabel(21 - ybin, "ME+" + suffix_label[ybin - 1], 2);
0074   }
0075   cscDQMOccupancy->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0076 
0077   mpcLinkErrors = ibooker.book2D("mpcLinkErrors", "MPC Link Errors", 54, 1, 55, 12, -6, 6);
0078   mpcLinkErrors->setAxisTitle("Sector (CSCID 1-9 Unlabelled)", 1);
0079   for (int xbin = 1; xbin < 7; ++xbin) {
0080     mpcLinkErrors->setBinLabel(xbin * 9 - 8, std::to_string(xbin), 1);
0081   }
0082   for (unsigned int ybin = 0; ybin < binNames.size(); ++ybin) {
0083     mpcLinkErrors->setBinLabel(ybin + 1, binNames[ybin], 2);
0084   }
0085 
0086   mpcLinkGood = ibooker.book2D("mpcLinkGood", "MPC Good Links", 54, 1, 55, 12, -6, 6);
0087   mpcLinkGood->setAxisTitle("Sector (CSCID 1-9 Unlabelled)", 1);
0088   for (int xbin = 1; xbin < 7; ++xbin) {
0089     mpcLinkGood->setBinLabel(xbin * 9 - 8, std::to_string(xbin), 1);
0090   }
0091   for (unsigned int ybin = 0; ybin < binNames.size(); ++ybin) {
0092     mpcLinkGood->setBinLabel(ybin + 1, binNames[ybin], 2);
0093   }
0094 
0095   // RPC Monitor Elements
0096   const std::array<std::string, 6> rpc_name{{"43", "42", "33", "32", "22", "12"}};
0097   const std::array<std::string, 6> rpc_label{{"4/3", "4/2", "3/3", "3/2", "2/2", "1/2"}};
0098 
0099   rpcHitBX = ibooker.book2D("rpcHitBX", "RPC Hit BX", 7, -3, 4, 12, 0, 12);
0100   rpcHitBX->setAxisTitle("BX", 1);
0101   for (int xbin = 1, xbin_label = -3; xbin <= 7; ++xbin, ++xbin_label) {
0102     rpcHitBX->setBinLabel(xbin, std::to_string(xbin_label), 1);
0103   }
0104   for (int ybin = 1; ybin <= 6; ++ybin) {
0105     rpcHitBX->setBinLabel(ybin, "RE-" + rpc_label[ybin - 1], 2);
0106     rpcHitBX->setBinLabel(13 - ybin, "RE+" + rpc_label[ybin - 1], 2);
0107   }
0108 
0109   rpcHitOccupancy = ibooker.book2D("rpcHitOccupancy", "RPC Chamber Occupancy", 42, 1, 43, 12, 0, 12);
0110   rpcHitOccupancy->setAxisTitle("Sector (Ni = Neighbor of Sector i)", 1);
0111   for (int bin = 1; bin <= 6; ++bin) {
0112     rpcHitOccupancy->setBinLabel(bin * 7 - 6, std::to_string(bin), 1);
0113     rpcHitOccupancy->setBinLabel(bin * 7, "N" + std::to_string(bin), 1);
0114     rpcHitOccupancy->setBinLabel(bin, "RE-" + rpc_label[bin - 1], 2);
0115     rpcHitOccupancy->setBinLabel(13 - bin, "RE+" + rpc_label[bin - 1], 2);
0116   }
0117   rpcHitOccupancy->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0118 
0119   // GEM Monitor Elements
0120   // Add GEM Oct 27 2020
0121   hitTypeBX = ibooker.book2D("hitTypeBX", "Hit Type BX", 4, 0.5, 4.5, 7, -3, 4);
0122   hitTypeBX->setBinLabel(1, "CSC", 1);
0123   hitTypeBX->setBinLabel(2, "RPC", 1);
0124   hitTypeBX->setBinLabel(3, "GEM", 1);
0125   hitTypeBX->setBinLabel(4, "Tot", 1);
0126   for (int ybin = 1; ybin < 8; ybin++)
0127     hitTypeBX->setBinLabel(ybin, std::to_string(ybin - 4), 2);
0128 
0129   gemHitBX = ibooker.book2D("gemHitBX", "GEM Hit BX", 7, -3, 4, 2, 0, 2);
0130   gemHitBX->setAxisTitle("BX", 1);
0131   for (int xbin = 1, xbin_label = -3; xbin <= 7; ++xbin, ++xbin_label) {
0132     gemHitBX->setBinLabel(xbin, std::to_string(xbin_label), 1);
0133   }
0134   gemHitBX->setBinLabel(1, "GE-1/1", 2);
0135   gemHitBX->setBinLabel(2, "GE+1/1", 2);
0136 
0137   gemHitOccupancy = ibooker.book2D("gemHitOccupancy", "GEM Chamber Occupancy", 42, 1, 43, 2, 0, 2);
0138   gemHitOccupancy->setAxisTitle("10#circ Chambers (Ni = Neighbor of Sector i)", 1);
0139   count = 0;
0140   for (int xbin = 1; xbin < 43; ++xbin) {
0141     gemHitOccupancy->setBinLabel(xbin, std::to_string(xbin - count), 1);
0142     if (xbin == 2 || xbin == 9 || xbin == 16 || xbin == 23 || xbin == 30 || xbin == 37) {
0143       ++xbin;
0144       ++count;
0145       gemHitOccupancy->setBinLabel(xbin, "N" + std::to_string(count), 1);
0146     }
0147   }
0148 
0149   gemHitOccupancy->setBinLabel(1, "GE-1/1", 2);
0150   gemHitOccupancy->setBinLabel(2, "GE+1/1", 2);
0151   gemHitOccupancy->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0152 
0153   // Track Monitor Elements
0154   emtfnTracks = ibooker.book1D("emtfnTracks", "Number of EMTF Tracks per Event", 11, 0, 11);
0155   for (int xbin = 1; xbin <= 10; ++xbin) {
0156     emtfnTracks->setBinLabel(xbin, std::to_string(xbin - 1), 1);
0157   }
0158   emtfnTracks->setBinLabel(11, "Overflow", 1);
0159 
0160   emtfTracknHits = ibooker.book1D("emtfTracknHits", "Number of Hits per EMTF Track", 5, 0, 5);
0161   for (int xbin = 1; xbin <= 5; ++xbin) {
0162     emtfTracknHits->setBinLabel(xbin, std::to_string(xbin - 1), 1);
0163   }
0164 
0165   emtfTrackBX = ibooker.book2D("emtfTrackBX", "EMTF Track Bunch Crossing", 12, -6, 6, 7, -3, 4);
0166   emtfTrackBX->setAxisTitle("Sector (Endcap)", 1);
0167   for (int xbin = 0; xbin < 6; ++xbin) {
0168     emtfTrackBX->setBinLabel(xbin + 1, std::to_string(6 - xbin) + " (-)", 1);
0169     emtfTrackBX->setBinLabel(12 - xbin, std::to_string(6 - xbin) + " (+)", 1);
0170   }
0171   emtfTrackBX->setAxisTitle("Track BX", 2);
0172   for (int ybin = 1, i = -3; ybin <= 7; ++ybin, ++i) {
0173     emtfTrackBX->setBinLabel(ybin, std::to_string(i), 2);
0174   }
0175 
0176   emtfTrackPt = ibooker.book1D("emtfTrackPt", "EMTF Track p_{T}", 256, 1, 257);
0177   emtfTrackPt->setAxisTitle("Track p_{T} [GeV]", 1);
0178 
0179   emtfTrackUnconstrainedPt = ibooker.book1D("emtfTrackUnconstrainedPt", "EMTF Track Unconstrained p_{T}", 256, 1, 257);
0180   emtfTrackUnconstrainedPt->setAxisTitle("Track Unconstrained p_{T} [GeV]", 1);
0181 
0182   emtfTrackDxy = ibooker.book1D("emtfTrackDxy", "EMTF Track d_{xy}", 3, 0, 3);
0183   emtfTrackDxy->setAxisTitle("Track d_{xy}", 1);
0184 
0185   emtfTrackEta = ibooker.book1D("emtfTrackEta", "EMTF Track #eta", 100, -2.5, 2.5);
0186   emtfTrackEta->setAxisTitle("Track #eta", 1);
0187 
0188   emtfTrackPhi = ibooker.book1D("emtfTrackPhi", "EMTF Track #phi", 126, -3.15, 3.15);
0189   emtfTrackPhi->setAxisTitle("Track #phi", 1);
0190 
0191   emtfTrackOccupancy = ibooker.book2D("emtfTrackOccupancy", "EMTF Track Occupancy", 100, -2.5, 2.5, 126, -3.15, 3.15);
0192   emtfTrackOccupancy->setAxisTitle("#eta", 1);
0193   emtfTrackOccupancy->setAxisTitle("#phi", 2);
0194 
0195   emtfTrackMode = ibooker.book1D("emtfTrackMode", "EMTF Track Mode", 16, 0, 16);
0196   emtfTrackMode->setAxisTitle("Mode", 1);
0197 
0198   emtfTrackQuality = ibooker.book1D("emtfTrackQuality", "EMTF Track Quality", 16, 0, 16);
0199   emtfTrackQuality->setAxisTitle("Quality", 1);
0200 
0201   emtfTrackQualityVsMode = ibooker.book2D("emtfTrackQualityVsMode", "EMTF Track Quality vs Mode", 16, 0, 16, 16, 0, 16);
0202   emtfTrackQualityVsMode->setAxisTitle("Mode", 1);
0203   emtfTrackQualityVsMode->setAxisTitle("Quality", 2);
0204 
0205   RPCvsEMTFTrackMode = ibooker.book2D("RPCvsEMTFTrackMode", "RPC Mode vs EMTF TrackMode", 16, 0, 16, 16, 0, 16);
0206   RPCvsEMTFTrackMode->setAxisTitle("EMTF Mode", 1);
0207   RPCvsEMTFTrackMode->setAxisTitle("RPC Mode", 2);
0208 
0209   for (int bin = 1; bin <= 16; ++bin) {
0210     emtfTrackMode->setBinLabel(bin, std::to_string(bin - 1), 1);
0211     emtfTrackQuality->setBinLabel(bin, std::to_string(bin - 1), 1);
0212     emtfTrackQualityVsMode->setBinLabel(bin, std::to_string(bin - 1), 1);
0213     emtfTrackQualityVsMode->setBinLabel(bin, std::to_string(bin - 1), 2);
0214     RPCvsEMTFTrackMode->setBinLabel(bin, std::to_string(bin - 1), 1);
0215     RPCvsEMTFTrackMode->setBinLabel(bin, std::to_string(bin - 1), 2);
0216   }
0217 
0218   //Chad Freer May 8 2018 (Selected Tracks)
0219   ibooker.setCurrentFolder(monitorDir + "/SelectedTracks");
0220 
0221   //Chad Freer May 8 2018 (High Quality Track Plots)
0222   emtfTrackPtHighQuality = ibooker.book1D("emtfTrackPtHighQuality", "EMTF High Quality Track p_{T}", 256, 1, 257);
0223   emtfTrackPtHighQuality->setAxisTitle("Track p_{T} [GeV] (Quality #geq 12)", 1);
0224 
0225   emtfTrackUnconstrainedPtHighQuality =
0226       ibooker.book1D("emtfTrackUnconstrainedPtHighQuality", "EMTF High Quality Track Unconstrained p_{T}", 256, 1, 257);
0227   emtfTrackUnconstrainedPtHighQuality->setAxisTitle("Track Unconstrained p_{T} [GeV] (Quality #geq 12)", 1);
0228 
0229   emtfTrackEtaHighQuality = ibooker.book1D("emtfTrackEtaHighQuality", "EMTF High Quality Track #eta", 100, -2.5, 2.5);
0230   emtfTrackEtaHighQuality->setAxisTitle("Track #eta (Quality #geq 12)", 1);
0231 
0232   emtfTrackPhiHighQuality = ibooker.book1D("emtfTrackPhiHighQuality", "EMTF High Quality #phi", 126, -3.15, 3.15);
0233   emtfTrackPhiHighQuality->setAxisTitle("Track #phi (Quality #geq 12)", 1);
0234 
0235   emtfTrackOccupancyHighQuality = ibooker.book2D(
0236       "emtfTrackOccupancyHighQuality", "EMTF High Quality Track Occupancy", 100, -2.5, 2.5, 126, -3.15, 3.15);
0237   emtfTrackOccupancyHighQuality->setAxisTitle("#eta", 1);
0238   emtfTrackOccupancyHighQuality->setAxisTitle("#phi", 2);
0239 
0240   //Chad Freer may 8 2018 (High Quality and High PT [22 GeV] Track Plots)
0241   emtfTrackPtHighQualityHighPT =
0242       ibooker.book1D("emtfTrackPtHighQualityHighPT", "EMTF High Quality High PT Track p_{T}", 256, 1, 257);
0243   emtfTrackPtHighQualityHighPT->setAxisTitle("Track p_{T} [GeV] (Quality #geq 12 and pT>22)", 1);
0244 
0245   emtfTrackEtaHighQualityHighPT =
0246       ibooker.book1D("emtfTrackEtaHighQualityHighPT", "EMTF High Quality High PT Track #eta", 100, -2.5, 2.5);
0247   emtfTrackEtaHighQualityHighPT->setAxisTitle("Track #eta (Quality #geq 12 and pT>22)", 1);
0248 
0249   emtfTrackPhiHighQualityHighPT =
0250       ibooker.book1D("emtfTrackPhiHighQualityHighPT", "EMTF High Quality High PT #phi", 126, -3.15, 3.15);
0251   emtfTrackPhiHighQualityHighPT->setAxisTitle("Track #phi (Quality #geq 12 and pT>22)", 1);
0252 
0253   emtfTrackOccupancyHighQualityHighPT = ibooker.book2D("emtfTrackOccupancyHighQualityHighPT",
0254                                                        "EMTF High Quality High PT Track Occupancy",
0255                                                        100,
0256                                                        -2.5,
0257                                                        2.5,
0258                                                        126,
0259                                                        -3.15,
0260                                                        3.15);
0261   emtfTrackOccupancyHighQualityHighPT->setAxisTitle("#eta", 1);
0262   emtfTrackOccupancyHighQualityHighPT->setAxisTitle("#phi", 2);
0263   //Chad Freer May 8 2018 (END new plots)
0264 
0265   emtfTrackUnconstrainedPtHighQualityHighUPT = ibooker.book1D(
0266       "emtfTrackUnconstrainedPtHighQualityHighUPT", "EMTF High Quality High UPT Track Unconstrained p_{T}", 256, 1, 257);
0267   emtfTrackUnconstrainedPtHighQualityHighUPT->setAxisTitle(
0268       "Track Unconstrained p_{T} [GeV] (Quality #geq 12 and UpT>10)", 1);
0269 
0270   emtfTrackEtaHighQualityHighUPT =
0271       ibooker.book1D("emtfTrackEtaHighQualityHighUPT", "EMTF High Quality High UPT Track #eta", 100, -2.5, 2.5);
0272   emtfTrackEtaHighQualityHighUPT->setAxisTitle("Track #eta (Quality #geq 12 and UpT>10)", 1);
0273 
0274   emtfTrackPhiHighQualityHighUPT =
0275       ibooker.book1D("emtfTrackPhiHighQualityHighUPT", "EMTF High Quality High UPT #phi", 126, -3.15, 3.15);
0276   emtfTrackPhiHighQualityHighUPT->setAxisTitle("Track #phi (Quality #geq 12 and UpT>10)", 1);
0277 
0278   emtfTrackOccupancyHighQualityHighUPT = ibooker.book2D("emtfTrackOccupancyHighQualityHighUPT",
0279                                                         "EMTF High Quality High UPT Track Occupancy",
0280                                                         100,
0281                                                         -2.5,
0282                                                         2.5,
0283                                                         126,
0284                                                         -3.15,
0285                                                         3.15);
0286   emtfTrackOccupancyHighQualityHighUPT->setAxisTitle("#eta", 1);
0287   emtfTrackOccupancyHighQualityHighUPT->setAxisTitle("#phi", 2);
0288 
0289   // CSC Input
0290   ibooker.setCurrentFolder(monitorDir + "/CSCInput");
0291 
0292   for (int hist = 0, i = 0; hist < 20; ++hist, i = hist % 10) {
0293     if (hist < 10) {
0294       name = "MENeg" + suffix_name[i];
0295       label = "ME-" + suffix_label[i];
0296     } else {
0297       name = "MEPos" + suffix_name[9 - i];
0298       label = "ME+" + suffix_label[9 - i];
0299     }
0300 
0301     if (hist < 6) {
0302       nChambs = (i % 2) ? 18 : 36;
0303     } else if (hist > 13) {
0304       nChambs = (i % 2) ? 36 : 18;
0305     } else {
0306       nChambs = 36;
0307     }
0308 
0309     const std::array<int, 10> wiregroups{{64, 96, 64, 96, 64, 112, 32, 64, 48, 48}};
0310     const std::array<int, 10> halfstrips{{160, 160, 160, 160, 160, 160, 128, 160, 128, 96}};
0311 
0312     if (hist < 10) {
0313       nWires = wiregroups[hist];
0314       nStrips = halfstrips[hist];
0315     } else {
0316       nWires = wiregroups[19 - hist];
0317       nStrips = halfstrips[19 - hist];
0318     }
0319 
0320     cscLCTStrip[hist] = ibooker.book1D("cscLCTStrip" + name, "CSC Halfstrip " + label, nStrips, 0, nStrips);
0321     cscLCTStrip[hist]->setAxisTitle("Cathode Halfstrip, " + label, 1);
0322 
0323     cscLCTWire[hist] = ibooker.book1D("cscLCTWire" + name, "CSC Wiregroup " + label, nWires, 0, nWires);
0324     cscLCTWire[hist]->setAxisTitle("Anode Wiregroup, " + label, 1);
0325 
0326     cscChamberStrip[hist] = ibooker.book2D(
0327         "cscChamberStrip" + name, "CSC Halfstrip " + label, nChambs, 1, 1 + nChambs, nStrips, 0, nStrips);
0328     cscChamberStrip[hist]->setAxisTitle("Chamber, " + label, 1);
0329     cscChamberStrip[hist]->setAxisTitle("Cathode Halfstrip", 2);
0330 
0331     cscChamberWire[hist] =
0332         ibooker.book2D("cscChamberWire" + name, "CSC Wiregroup " + label, nChambs, 1, 1 + nChambs, nWires, 0, nWires);
0333     cscChamberWire[hist]->setAxisTitle("Chamber, " + label, 1);
0334     cscChamberWire[hist]->setAxisTitle("Anode Wiregroup", 2);
0335 
0336     for (int bin = 1; bin <= nChambs; ++bin) {
0337       cscChamberStrip[hist]->setBinLabel(bin, std::to_string(bin), 1);
0338       cscChamberWire[hist]->setBinLabel(bin, std::to_string(bin), 1);
0339     }
0340   }
0341 
0342   // RPC Input
0343   ibooker.setCurrentFolder(monitorDir + "/RPCInput");
0344 
0345   for (int hist = 0, i = 0; hist < 12; ++hist, i = hist % 6) {
0346     if (hist < 6) {
0347       name = "RENeg" + rpc_name[i];
0348       label = "RE-" + rpc_label[i];
0349     } else {
0350       name = "REPos" + rpc_name[5 - i];
0351       label = "RE+" + rpc_label[5 - i];
0352     }
0353     rpcHitPhi[hist] = ibooker.book1D("rpcHitPhi" + name, "RPC Hit Phi " + label, 1250, 0, 1250);
0354     rpcHitPhi[hist]->setAxisTitle("#phi", 1);
0355     rpcHitTheta[hist] = ibooker.book1D("rpcHitTheta" + name, "RPC Hit Theta " + label, 32, 0, 32);
0356     rpcHitTheta[hist]->setAxisTitle("#theta", 1);
0357     rpcChamberPhi[hist] = ibooker.book2D("rpcChamberPhi" + name, "RPC Chamber Phi " + label, 36, 1, 37, 1250, 0, 1250);
0358     rpcChamberPhi[hist]->setAxisTitle("Chamber", 1);
0359     rpcChamberPhi[hist]->setAxisTitle("#phi", 2);
0360     rpcChamberTheta[hist] =
0361         ibooker.book2D("rpcChamberTheta" + name, "RPC Chamber Theta " + label, 36, 1, 37, 32, 0, 32);
0362     rpcChamberTheta[hist]->setAxisTitle("Chamber", 1);
0363     rpcChamberTheta[hist]->setAxisTitle("#theta", 2);
0364     for (int xbin = 1; xbin < 37; ++xbin) {
0365       rpcChamberPhi[hist]->setBinLabel(xbin, std::to_string(xbin), 1);
0366       rpcChamberTheta[hist]->setBinLabel(xbin, std::to_string(xbin), 1);
0367     }
0368   }
0369 
0370   // GEM Input Nov 03 2020
0371   ibooker.setCurrentFolder(monitorDir + "/GEMInput");
0372   for (int hist = 0; hist < 2; hist++) {
0373     if (hist == 0) {
0374       name = "GENeg11";
0375       label = "GE-1/1";
0376     }
0377     if (hist == 1) {
0378       name = "GEPos11";
0379       label = "GE+1/1";
0380     }
0381     nChambs = 36;
0382     nStrips = 192;  //use nStrips for number of pads
0383     gemChamberPad[hist] = ibooker.book2D(
0384         "gemChamberPad" + name, "GEM Chamber Pad " + label, nChambs, 1, 1 + nChambs, nStrips, 0, nStrips);  // pads 0-191
0385     gemChamberPad[hist]->setAxisTitle("Chamber, " + label, 1);
0386     gemChamberPad[hist]->setAxisTitle("Pad", 2);
0387     gemChamberPartition[hist] =
0388         ibooker.book2D("gemChamberPartition" + name,
0389                        "GEM Chamber Partition " + label,
0390                        nChambs,
0391                        1,
0392                        1 + nChambs,
0393                        9,
0394                        0,
0395                        9);  // partitions 1-8 or 0-7. There have been changes in different firmware/unpacker verisions.
0396     gemChamberPartition[hist]->setAxisTitle("Chamber, " + label, 1);
0397     gemChamberPartition[hist]->setAxisTitle("Partition", 2);
0398     for (int bin = 1; bin <= nChambs; ++bin) {
0399       gemChamberPad[hist]->setBinLabel(bin, std::to_string(bin), 1);
0400       gemChamberPartition[hist]->setBinLabel(bin, std::to_string(bin), 1);
0401     }
0402     //Added 07-21-22 **
0403     for (int ch = 0; ch < 36; ch++) {
0404       for (int lyr = 0; lyr < 2; lyr++) {
0405         gemVFATBXPerChamber[ch][hist][lyr] = ibooker.book2D(
0406             "gemVFATBXPerChamber_" + std::to_string(ch) + "_" + std::to_string(hist) + "_" + std::to_string(lyr + 1),
0407             "GEM BX vs VFAT in Chamber " + std::to_string(ch + 1) + " " + label + " Layer " + std::to_string(lyr + 1),
0408             7,
0409             -3,
0410             4,
0411             24,
0412             0,
0413             24);
0414         gemVFATBXPerChamber[ch][hist][lyr]->setAxisTitle("BX", 1);
0415         gemVFATBXPerChamber[ch][hist][lyr]->setAxisTitle("VFAT #", 2);
0416 
0417         for (int bin = 1; bin <= 24; ++bin) {
0418           gemVFATBXPerChamber[ch][hist][lyr]->setBinLabel(bin, std::to_string(bin - 1), 2);
0419         }
0420         for (int bx = 1; bx <= 7; ++bx) {
0421           gemVFATBXPerChamber[ch][hist][lyr]->setBinLabel(bx, std::to_string(bx - 4), 1);
0422         }
0423       }
0424     }
0425 
0426     //changed gemChamberVFATBX to be indexed by BX 07-21-2022
0427     string bx_string;
0428     for (int bx = 1; bx <= 7; ++bx) {
0429       //Assign (m)inus or (p)us to plot name
0430       if (bx < 4)
0431         bx_string = "Neg" + std::to_string(-1 * (bx - 4));
0432       else if (bx > 4)
0433         bx_string = "Pos" + std::to_string(bx - 4);
0434       else
0435         bx_string = "0";
0436 
0437       gemChamberVFATBX[hist][bx - 1] =
0438           ibooker.book2D("gemChamberVFATBX" + bx_string + name,
0439                          "GEM Chamber vs VFAT at BX = " + std::to_string(bx - 4) + ", " + label,
0440                          42,
0441                          1,
0442                          43,
0443                          24,
0444                          0,
0445                          24);  // 8* (0-2) phi part + (0-7) eta part
0446       gemChamberVFATBX[hist][bx - 1]->setAxisTitle("Chamber, (Ni = Neighbor of Sector i), " + label, 1);
0447       gemChamberVFATBX[hist][bx - 1]->setAxisTitle("VFAT #", 2);
0448 
0449       for (int bin = 1; bin <= 24; bin++)
0450         gemChamberVFATBX[hist][bx - 1]->setBinLabel(bin, std::to_string(bin - 1), 2);
0451 
0452       int count = 0;
0453       for (int bin = 1; bin <= 42; ++bin) {
0454         gemChamberVFATBX[hist][bx - 1]->setBinLabel(bin, std::to_string(bin - count), 1);
0455         if (bin == 2 || bin == 9 || bin == 16 || bin == 23 || bin == 30 || bin == 37) {
0456           ++bin;
0457           ++count;
0458           gemChamberVFATBX[hist][bx - 1]->setBinLabel(bin, "N" + std::to_string(count), 1);
0459         }
0460       }
0461       gemChamberVFATBX[hist][bx - 1]->getTH2F()->GetXaxis()->SetCanExtend(false);
0462     }
0463   }
0464   // CSC LCT and RPC Hit Timing
0465   ibooker.setCurrentFolder(monitorDir + "/Timing");
0466 
0467   cscTimingTot = ibooker.book2D("cscTimingTotal", "CSC Total BX ", 42, 1, 43, 20, 0, 20);
0468   cscTimingTot->setAxisTitle("10#circ Chamber, (Ni = Neighbor of Sector i)", 1);
0469 
0470   rpcHitTimingTot = ibooker.book2D("rpcHitTimingTot", "RPC Chamber Occupancy ", 42, 1, 43, 12, 0, 12);
0471   rpcHitTimingTot->setAxisTitle("Sector (Ni = Neighbor of Sector i)", 1);
0472 
0473   gemHitTimingTot =
0474       ibooker.book2D("gemHitTimingTot", "GEM Chamber Occupancy ", 42, 1, 43, 2, 0, 2);  // Add GEM Timing Oct 27 2020
0475   gemHitTimingTot->setAxisTitle("10#circ Chamber (Ni = Neighbor of Sector i)", 1);
0476   const std::array<std::string, 5> nameBX{{"BXNeg1", "BXPos1", "BXNeg2", "BXPos2", "BX0"}};
0477   const std::array<std::string, 5> labelBX{{"BX -1", "BX +1", "BX -2", "BX +2", "BX 0"}};
0478 
0479   for (int hist = 0; hist < 5; ++hist) {
0480     count = 0;
0481     cscLCTTiming[hist] =
0482         ibooker.book2D("cscLCTTiming" + nameBX[hist], "CSC Chamber Occupancy " + labelBX[hist], 42, 1, 43, 20, 0, 20);
0483     cscLCTTiming[hist]->setAxisTitle("10#circ Chamber, (Ni = Neighbor of Sector i)", 1);
0484 
0485     for (int xbin = 1; xbin < 43; ++xbin) {
0486       cscLCTTiming[hist]->setBinLabel(xbin, std::to_string(xbin - count), 1);
0487       if (hist == 0)
0488         cscTimingTot->setBinLabel(xbin, std::to_string(xbin - count), 1);  //only fill once in the loop
0489       if (xbin == 2 || xbin == 9 || xbin == 16 || xbin == 23 || xbin == 30 || xbin == 37) {
0490         ++xbin;
0491         ++count;
0492         cscLCTTiming[hist]->setBinLabel(xbin, "N" + std::to_string(count), 1);
0493         if (hist == 0)
0494           cscTimingTot->setBinLabel(xbin, "N" + std::to_string(count), 1);
0495       }
0496     }
0497 
0498     for (int ybin = 1; ybin <= 10; ++ybin) {
0499       cscLCTTiming[hist]->setBinLabel(ybin, "ME-" + suffix_label[ybin - 1], 2);
0500       cscLCTTiming[hist]->setBinLabel(21 - ybin, "ME+" + suffix_label[ybin - 1], 2);
0501       if (hist == 0)
0502         cscTimingTot->setBinLabel(ybin, "ME-" + suffix_label[ybin - 1], 2);
0503       if (hist == 0)
0504         cscTimingTot->setBinLabel(21 - ybin, "ME+" + suffix_label[ybin - 1], 2);
0505     }
0506     if (hist == 0)
0507       cscTimingTot->getTH2F()->GetXaxis()->SetCanExtend(false);      // Needed to stop multi-thread summing
0508     cscLCTTiming[hist]->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0509 
0510     rpcHitTiming[hist] =
0511         ibooker.book2D("rpcHitTiming" + nameBX[hist], "RPC Chamber Occupancy " + labelBX[hist], 42, 1, 43, 12, 0, 12);
0512     rpcHitTiming[hist]->setAxisTitle("Sector, (Ni=Neighbor of Sector i )", 1);
0513     for (int bin = 1; bin < 7; ++bin) {
0514       rpcHitTiming[hist]->setBinLabel(bin * 7 - 6, std::to_string(bin), 1);
0515       rpcHitTiming[hist]->setBinLabel(bin * 7, "N" + std::to_string(bin), 1);
0516       rpcHitTiming[hist]->setBinLabel(bin, "RE-" + rpc_label[bin - 1], 2);
0517       rpcHitTiming[hist]->setBinLabel(13 - bin, "RE+" + rpc_label[bin - 1], 2);
0518     }
0519     rpcHitTiming[hist]->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0520     if (hist == 0) {
0521       for (int bin = 1; bin < 7; ++bin) {
0522         rpcHitTimingTot->setBinLabel(bin * 7 - 6, std::to_string(bin), 1);
0523         rpcHitTimingTot->setBinLabel(bin * 7, "N" + std::to_string(bin), 1);
0524         rpcHitTimingTot->setBinLabel(bin, "RE-" + rpc_label[bin - 1], 2);
0525         rpcHitTimingTot->setBinLabel(13 - bin, "RE+" + rpc_label[bin - 1], 2);
0526       }
0527       rpcHitTimingTot->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0528     }
0529     //if (hist == 4) continue; // Don't book for BX = 0
0530 
0531     // Add GEM Timing Oct 27 2020
0532     gemHitTiming[hist] =
0533         ibooker.book2D("gemHitTiming" + nameBX[hist], "GEM Chamber Occupancy " + labelBX[hist], 42, 1, 43, 2, 0, 2);
0534     gemHitTiming[hist]->setAxisTitle("10#circ Chamber, (Ni = Neighbor of Sector i)", 1);
0535     count = 0;
0536     for (int xbin = 1; xbin < 43; ++xbin) {
0537       gemHitTiming[hist]->setBinLabel(xbin, std::to_string(xbin - count), 1);
0538       if (hist == 0)
0539         gemHitTimingTot->setBinLabel(xbin, std::to_string(xbin - count), 1);  //only fill once in the loop
0540       if (xbin == 2 || xbin == 9 || xbin == 16 || xbin == 23 || xbin == 30 || xbin == 37) {
0541         ++xbin;
0542         ++count;
0543         gemHitTiming[hist]->setBinLabel(xbin, "N" + std::to_string(count), 1);
0544         if (hist == 0)
0545           gemHitTimingTot->setBinLabel(xbin, "N" + std::to_string(count), 1);
0546       }
0547     }
0548 
0549     gemHitTiming[hist]->setBinLabel(1, "GE-1/1", 2);
0550     gemHitTiming[hist]->setBinLabel(2, "GE+1/1", 2);
0551     if (hist == 0) {
0552       gemHitTimingTot->setBinLabel(1, "GE-1/1", 2);
0553       gemHitTimingTot->setBinLabel(2, "GE+1/1", 2);
0554       gemHitTimingTot->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0555     }
0556     gemHitTiming[hist]->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0557 
0558     count = 0;
0559     cscLCTTimingFrac[hist] = ibooker.book2D(
0560         "cscLCTTimingFrac" + nameBX[hist], "CSC Chamber Occupancy " + labelBX[hist], 42, 1, 43, 20, 0, 20);
0561     cscLCTTimingFrac[hist]->setAxisTitle("10#circ Chambers, (Ni = Neighbor of Sector i)", 1);
0562     for (int xbin = 1; xbin < 43; ++xbin) {
0563       cscLCTTimingFrac[hist]->setBinLabel(xbin, std::to_string(xbin - count), 1);
0564       if (xbin == 2 || xbin == 9 || xbin == 16 || xbin == 23 || xbin == 30 || xbin == 37) {
0565         ++xbin;
0566         ++count;
0567         cscLCTTimingFrac[hist]->setBinLabel(xbin, "N" + std::to_string(count), 1);
0568       }
0569     }
0570     for (int ybin = 1; ybin <= 10; ++ybin) {
0571       cscLCTTimingFrac[hist]->setBinLabel(ybin, "ME-" + suffix_label[ybin - 1], 2);
0572       cscLCTTimingFrac[hist]->setBinLabel(21 - ybin, "ME+" + suffix_label[ybin - 1], 2);
0573     }
0574     cscLCTTimingFrac[hist]->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0575 
0576     rpcHitTimingFrac[hist] = ibooker.book2D(
0577         "rpcHitTimingFrac" + nameBX[hist], "RPC Chamber Fraction in " + labelBX[hist], 42, 1, 43, 12, 0, 12);
0578     rpcHitTimingFrac[hist]->setAxisTitle("Sector, (Ni = Neighbor of Sector i)", 1);
0579     for (int bin = 1; bin < 7; ++bin) {
0580       rpcHitTimingFrac[hist]->setBinLabel(bin * 7 - 6, std::to_string(bin), 1);
0581       rpcHitTimingFrac[hist]->setBinLabel(bin * 7, "N" + std::to_string(bin), 1);
0582       rpcHitTimingFrac[hist]->setBinLabel(bin, "RE-" + rpc_label[bin - 1], 2);
0583       rpcHitTimingFrac[hist]->setBinLabel(13 - bin, "RE+" + rpc_label[bin - 1], 2);
0584     }
0585 
0586     // Add GEM Timing Oct 27 2020
0587     gemHitTimingFrac[hist] =
0588         ibooker.book2D("gemHitTimingFrac" + nameBX[hist], "GEM Chamber Occupancy " + labelBX[hist], 42, 1, 43, 2, 0, 2);
0589     gemHitTimingFrac[hist]->setAxisTitle("10#circ Chambers, (Ni = Neighbor of Sector i)", 1);
0590     count = 0;
0591     for (int xbin = 1; xbin < 43; ++xbin) {
0592       gemHitTimingFrac[hist]->setBinLabel(xbin, std::to_string(xbin - count), 1);
0593       if (xbin == 2 || xbin == 9 || xbin == 16 || xbin == 23 || xbin == 30 || xbin == 37) {
0594         ++xbin;
0595         ++count;
0596         gemHitTimingFrac[hist]->setBinLabel(xbin, "N" + std::to_string(count), 1);
0597       }
0598     }
0599     gemHitTimingFrac[hist]->setBinLabel(1, "GE-1/1", 2);
0600     gemHitTimingFrac[hist]->setBinLabel(2, "GE+1/1", 2);
0601     gemHitTimingFrac[hist]->getTH2F()->GetXaxis()->SetCanExtend(false);  // Needed to stop multi-thread summing
0602   }
0603 
0604   rpcHitTimingInTrack =
0605       ibooker.book2D("rpcHitTimingInTrack", "RPC Hit Timing (matched to track in BX 0)", 7, -3, 4, 12, 0, 12);
0606   rpcHitTimingInTrack->setAxisTitle("BX", 1);
0607   for (int xbin = 1, xbin_label = -3; xbin <= 7; ++xbin, ++xbin_label) {
0608     rpcHitTimingInTrack->setBinLabel(xbin, std::to_string(xbin_label), 1);
0609   }
0610   for (int ybin = 1; ybin <= 6; ++ybin) {
0611     rpcHitTimingInTrack->setBinLabel(ybin, "RE-" + rpc_label[ybin - 1], 2);
0612     rpcHitTimingInTrack->setBinLabel(13 - ybin, "RE+" + rpc_label[ybin - 1], 2);
0613   }
0614 
0615   const std::array<std::string, 3> nameNumStation{{"4Station", "3Station", "2Station"}};
0616   const std::array<std::string, 3> labelNumStation{{"4 Station Track", "3 Station Track", "2 Station Track"}};
0617 
0618   for (int hist = 0; hist < 3; ++hist) {
0619     emtfTrackBXVsCSCLCT[hist] = ibooker.book2D("emtfTrackBXVsCSCLCT" + nameNumStation[hist],
0620                                                "EMTF " + labelNumStation[hist] + " BX vs CSC LCT BX",
0621                                                7,
0622                                                -3,
0623                                                4,
0624                                                7,
0625                                                -3,
0626                                                4);
0627     emtfTrackBXVsCSCLCT[hist]->setAxisTitle("LCT BX", 1);
0628     emtfTrackBXVsCSCLCT[hist]->setAxisTitle("Track BX", 2);
0629     for (int bin = 1, bin_label = -3; bin <= 7; ++bin, ++bin_label) {
0630       emtfTrackBXVsCSCLCT[hist]->setBinLabel(bin, std::to_string(bin_label), 1);
0631       emtfTrackBXVsCSCLCT[hist]->setBinLabel(bin, std::to_string(bin_label), 2);
0632     }
0633     emtfTrackBXVsRPCHit[hist] = ibooker.book2D("emtfTrackBXVsRPCHit" + nameNumStation[hist],
0634                                                "EMTF " + labelNumStation[hist] + " BX vs RPC Hit BX",
0635                                                7,
0636                                                -3,
0637                                                4,
0638                                                7,
0639                                                -3,
0640                                                4);
0641     emtfTrackBXVsRPCHit[hist]->setAxisTitle("Hit BX", 1);
0642     emtfTrackBXVsRPCHit[hist]->setAxisTitle("Track BX", 2);
0643     for (int bin = 1, bin_label = -3; bin <= 7; ++bin, ++bin_label) {
0644       emtfTrackBXVsRPCHit[hist]->setBinLabel(bin, std::to_string(bin_label), 1);
0645       emtfTrackBXVsRPCHit[hist]->setBinLabel(bin, std::to_string(bin_label), 2);
0646     }
0647     // Add GEM vs track BX Dec 05 2020
0648     emtfTrackBXVsGEMHit[hist] = ibooker.book2D("emtfTrackBXVsGEMHit" + nameNumStation[hist],
0649                                                "EMTF " + labelNumStation[hist] + " BX vs GEM Hit BX",
0650                                                7,
0651                                                -3,
0652                                                4,
0653                                                7,
0654                                                -3,
0655                                                4);
0656     emtfTrackBXVsGEMHit[hist]->setAxisTitle("Hit BX", 1);
0657     emtfTrackBXVsGEMHit[hist]->setAxisTitle("Track BX", 2);
0658     for (int bin = 1, bin_label = -3; bin <= 7; ++bin, ++bin_label) {
0659       emtfTrackBXVsGEMHit[hist]->setBinLabel(bin, std::to_string(bin_label), 1);
0660       emtfTrackBXVsGEMHit[hist]->setBinLabel(bin, std::to_string(bin_label), 2);
0661     }
0662   }  // End loop: for (int hist = 0; hist < 3; ++hist)
0663 
0664   // Add mode vs BXdiff comparison Dec 07 2020
0665   const std::array<std::string, 2> nameGEMStation{{"GENeg11", "GEPos11"}};
0666   const std::array<std::string, 2> labelGEMStation{{"GE-1/1", "GE+1/1"}};
0667   const std::array<std::string, 8> nameCSCStation{
0668       {"MENeg4", "MENeg3", "MENeg2", "MENeg1", "MEPos1", "MEPos2", "MEPos3", "MEPos4"}};
0669   const std::array<std::string, 8> labelCSCStation{{"ME-4", "ME-3", "ME-2", "ME-1", "ME+1", "ME+2", "ME+3", "ME+4"}};
0670   const std::array<std::string, 8> nameRPCStation{
0671       {"RENeg4", "RENeg3", "RENeg2", "RENeg1", "REPos1", "REPos2", "REPos3", "REPos4"}};
0672   const std::array<std::string, 8> labelRPCStation{{"RE-4", "RE-3", "RE-2", "RE-1", "RE+1", "RE+2", "RE+3", "RE+4"}};
0673 
0674   for (int iGEM = 0; iGEM < 2; iGEM++) {
0675     emtfTrackModeVsGEMBXDiff[iGEM] = ibooker.book2D("emtfTrackModeVsGEMBXDiff" + nameGEMStation[iGEM],
0676                                                     "EMTF Track Mode vs (Track BX - GEM BX) " + labelGEMStation[iGEM],
0677                                                     9,
0678                                                     -4,
0679                                                     5,
0680                                                     16,
0681                                                     0,
0682                                                     16);
0683     emtfTrackModeVsGEMBXDiff[iGEM]->setAxisTitle("Track BX - GEM BX", 1);
0684     emtfTrackModeVsGEMBXDiff[iGEM]->setAxisTitle("Track Mode", 2);
0685   }
0686   for (int iCSC = 0; iCSC < 8; iCSC++) {
0687     emtfTrackModeVsCSCBXDiff[iCSC] = ibooker.book2D("emtfTrackModeVsCSCBXDiff" + nameCSCStation[iCSC],
0688                                                     "EMTF Track Mode vs (Track BX - LCT BX) " + labelCSCStation[iCSC],
0689                                                     9,
0690                                                     -4,
0691                                                     5,
0692                                                     16,
0693                                                     0,
0694                                                     16);
0695     emtfTrackModeVsCSCBXDiff[iCSC]->setAxisTitle("Track BX - LCT BX", 1);
0696     emtfTrackModeVsCSCBXDiff[iCSC]->setAxisTitle("Track Mode", 2);
0697   }
0698   for (int iRPC = 0; iRPC < 8; iRPC++) {
0699     emtfTrackModeVsRPCBXDiff[iRPC] = ibooker.book2D("emtfTrackModeVsRPCBXDiff" + nameRPCStation[iRPC],
0700                                                     "EMTF Track Mode vs (Track BX - RPC BX) " + labelRPCStation[iRPC],
0701                                                     9,
0702                                                     -4,
0703                                                     5,
0704                                                     16,
0705                                                     0,
0706                                                     16);
0707     emtfTrackModeVsRPCBXDiff[iRPC]->setAxisTitle("Track BX - RPC BX", 1);
0708     emtfTrackModeVsRPCBXDiff[iRPC]->setAxisTitle("Track Mode", 2);
0709   }
0710 
0711   // GEM vs CSC Dec 06 2020
0712   ibooker.setCurrentFolder(monitorDir + "/GEMVsCSC");
0713   for (int hist = 0; hist < 2; hist++) {
0714     gemHitPhi[hist] = ibooker.book2D(
0715         "gemHitPhi" + nameGEMStation[hist], "GEM Hit Phi " + labelGEMStation[hist], 4921, 0, 4921, 6, 1, 7);
0716     gemHitTheta[hist] = ibooker.book2D(
0717         "gemHitTheta" + nameGEMStation[hist], "GEM Hit Theta " + labelGEMStation[hist], 128, 0, 128, 6, 1, 7);
0718     gemHitVScscLCTPhi[hist] = ibooker.book2D("gemHitVScscLCTPhi" + nameGEMStation[hist],
0719                                              "GEM Hit Phi - CSC LCT Phi " + labelGEMStation[hist],
0720                                              1200,
0721                                              -600,
0722                                              600,
0723                                              36,
0724                                              1,
0725                                              37);  // one chamber is 10 degrees, 60 integer phi per degree
0726     gemHitVScscLCTTheta[hist] = ibooker.book2D("gemHitVScscLCTTheta" + nameGEMStation[hist],
0727                                                "GEM Hit Theta - CSC LCT Theta " + labelGEMStation[hist],
0728                                                20,
0729                                                -10,
0730                                                10,
0731                                                36,
0732                                                1,
0733                                                37);  // 0.1 eta is at most 9.5 integer theta (between eta 1.5 and 1.6)
0734     gemHitVScscLCTBX[hist] = ibooker.book2D("gemHitVScscLCTBX" + nameGEMStation[hist],
0735                                             "GEM Hit BX - CSC LCT BX " + labelGEMStation[hist],
0736                                             9,
0737                                             -4,
0738                                             5,
0739                                             36,
0740                                             1,
0741                                             37);
0742 
0743     gemHitPhi[hist]->setAxisTitle("Integer #phi", 1);
0744     gemHitTheta[hist]->setAxisTitle("Integer #theta", 1);
0745     gemHitVScscLCTPhi[hist]->setAxisTitle("Integer #phi", 1);
0746     gemHitVScscLCTTheta[hist]->setAxisTitle("Integer #theta", 1);
0747     gemHitVScscLCTBX[hist]->setAxisTitle("GEM BX - CSC BX", 1);
0748 
0749     gemHitPhi[hist]->setAxisTitle("Sector", 2);
0750     gemHitTheta[hist]->setAxisTitle("Sector", 2);
0751     gemHitVScscLCTPhi[hist]->setAxisTitle("Chamber", 2);
0752     gemHitVScscLCTTheta[hist]->setAxisTitle("Chamber", 2);
0753     gemHitVScscLCTBX[hist]->setAxisTitle("Chamber", 2);
0754   }
0755 
0756   // Muon Cand
0757   ibooker.setCurrentFolder(monitorDir + "/MuonCand");
0758 
0759   // Regional Muon Candidate Monitor Elements
0760   emtfMuonBX = ibooker.book1D("emtfMuonBX", "EMTF Muon Cand BX", 7, -3, 4);
0761   emtfMuonBX->setAxisTitle("BX", 1);
0762   for (int xbin = 1, bin_label = -3; xbin <= 7; ++xbin, ++bin_label) {
0763     emtfMuonBX->setBinLabel(xbin, std::to_string(bin_label), 1);
0764   }
0765 
0766   emtfMuonhwPt = ibooker.book1D("emtfMuonhwPt", "EMTF Muon Cand p_{T}", 512, 0, 512);
0767   emtfMuonhwPt->setAxisTitle("Hardware p_{T}", 1);
0768 
0769   emtfMuonhwPtUnconstrained =
0770       ibooker.book1D("emtfMuonhwPtUnconstrained", "EMTF Muon Cand Unconstrained p_{T}", 256, 0, 256);
0771   emtfMuonhwPtUnconstrained->setAxisTitle("Hardware Unconstrained p_{T}", 1);
0772 
0773   emtfMuonhwDxy = ibooker.book1D("emtfMuonhwDxy", "EMTF Muon Cand d_{xy}", 3, 0, 3);
0774   emtfMuonhwDxy->setAxisTitle("Hardware d_{xy}", 1);
0775 
0776   emtfMuonhwEta = ibooker.book1D("emtfMuonhwEta", "EMTF Muon Cand #eta", 460, -230, 230);
0777   emtfMuonhwEta->setAxisTitle("Hardware #eta", 1);
0778 
0779   emtfMuonhwPhi = ibooker.book1D("emtfMuonhwPhi", "EMTF Muon Cand #phi", 145, -40, 105);
0780   emtfMuonhwPhi->setAxisTitle("Hardware #phi", 1);
0781 
0782   emtfMuonhwQual = ibooker.book1D("emtfMuonhwQual", "EMTF Muon Cand Quality", 16, 0, 16);
0783   emtfMuonhwQual->setAxisTitle("Quality", 1);
0784   for (int xbin = 1; xbin <= 16; ++xbin) {
0785     emtfMuonhwQual->setBinLabel(xbin, std::to_string(xbin - 1), 1);
0786   }
0787 }
0788 
0789 // CSCOccupancy chamber mapping for neighbor inclusive plots
0790 int chamber_bin(int station, int ring, int chamber) {
0791   int chamber_bin_index = 0;
0792   if (station > 1 && (ring % 2) == 1) {
0793     chamber_bin_index = (chamber * 2) + ((chamber + 1) / 3);
0794   } else {
0795     chamber_bin_index = chamber + ((chamber + 3) / 6);
0796   }
0797   return chamber_bin_index;
0798 };
0799 
0800 void L1TStage2EMTF::analyze(const edm::Event& e, const edm::EventSetup& c) {
0801   if (verbose)
0802     edm::LogInfo("L1TStage2EMTF") << "L1TStage2EMTF: analyze..." << std::endl;
0803 
0804   // DAQ Output
0805   edm::Handle<l1t::EMTFDaqOutCollection> DaqOutCollection;
0806   e.getByToken(daqToken, DaqOutCollection);
0807 
0808   for (auto DaqOut = DaqOutCollection->begin(); DaqOut != DaqOutCollection->end(); ++DaqOut) {
0809     const l1t::emtf::MECollection* MECollection = DaqOut->PtrMECollection();
0810     for (auto ME = MECollection->begin(); ME != MECollection->end(); ++ME) {
0811       if (ME->SE())
0812         emtfErrors->Fill(1);
0813       if (ME->SM())
0814         emtfErrors->Fill(2);
0815       if (ME->BXE())
0816         emtfErrors->Fill(3);
0817       if (ME->AF())
0818         emtfErrors->Fill(4);
0819     }
0820 
0821     const l1t::emtf::EventHeader* EventHeader = DaqOut->PtrEventHeader();
0822     if (!EventHeader->Rdy())
0823       emtfErrors->Fill(5);
0824 
0825     // Fill MPC input link errors
0826     int offset = (EventHeader->Sector() - 1) * 9;
0827     int endcap = EventHeader->Endcap();
0828     l1t::emtf::Counters CO = DaqOut->GetCounters();
0829     const std::array<std::array<int, 9>, 5> counters{
0830         {{{CO.ME1a_1(),
0831            CO.ME1a_2(),
0832            CO.ME1a_3(),
0833            CO.ME1a_4(),
0834            CO.ME1a_5(),
0835            CO.ME1a_6(),
0836            CO.ME1a_7(),
0837            CO.ME1a_8(),
0838            CO.ME1a_9()}},
0839          {{CO.ME1b_1(),
0840            CO.ME1b_2(),
0841            CO.ME1b_3(),
0842            CO.ME1b_4(),
0843            CO.ME1b_5(),
0844            CO.ME1b_6(),
0845            CO.ME1b_7(),
0846            CO.ME1b_8(),
0847            CO.ME1b_9()}},
0848          {{CO.ME2_1(), CO.ME2_2(), CO.ME2_3(), CO.ME2_4(), CO.ME2_5(), CO.ME2_6(), CO.ME2_7(), CO.ME2_8(), CO.ME2_9()}},
0849          {{CO.ME3_1(), CO.ME3_2(), CO.ME3_3(), CO.ME3_4(), CO.ME3_5(), CO.ME3_6(), CO.ME3_7(), CO.ME3_8(), CO.ME3_9()}},
0850          {{CO.ME4_1(), CO.ME4_2(), CO.ME4_3(), CO.ME4_4(), CO.ME4_5(), CO.ME4_6(), CO.ME4_7(), CO.ME4_8(), CO.ME4_9()}}}};
0851     for (int i = 0; i < 5; i++) {
0852       for (int j = 0; j < 9; j++) {
0853         if (counters.at(i).at(j) != 0)
0854           mpcLinkErrors->Fill(j + 1 + offset, endcap * (i + 0.5), counters.at(i).at(j));
0855         else
0856           mpcLinkGood->Fill(j + 1 + offset, endcap * (i + 0.5));
0857       }
0858     }
0859     if (CO.ME1n_3() == 1)
0860       mpcLinkErrors->Fill(1 + offset, endcap * 5.5);
0861     if (CO.ME1n_6() == 1)
0862       mpcLinkErrors->Fill(2 + offset, endcap * 5.5);
0863     if (CO.ME1n_9() == 1)
0864       mpcLinkErrors->Fill(3 + offset, endcap * 5.5);
0865     if (CO.ME2n_3() == 1)
0866       mpcLinkErrors->Fill(4 + offset, endcap * 5.5);
0867     if (CO.ME2n_9() == 1)
0868       mpcLinkErrors->Fill(5 + offset, endcap * 5.5);
0869     if (CO.ME3n_3() == 1)
0870       mpcLinkErrors->Fill(6 + offset, endcap * 5.5);
0871     if (CO.ME3n_9() == 1)
0872       mpcLinkErrors->Fill(7 + offset, endcap * 5.5);
0873     if (CO.ME4n_3() == 1)
0874       mpcLinkErrors->Fill(8 + offset, endcap * 5.5);
0875     if (CO.ME4n_9() == 1)
0876       mpcLinkErrors->Fill(9 + offset, endcap * 5.5);
0877     if (CO.ME1n_3() == 0)
0878       mpcLinkGood->Fill(1 + offset, endcap * 5.5);
0879     if (CO.ME1n_6() == 0)
0880       mpcLinkGood->Fill(2 + offset, endcap * 5.5);
0881     if (CO.ME1n_9() == 0)
0882       mpcLinkGood->Fill(3 + offset, endcap * 5.5);
0883     if (CO.ME2n_3() == 0)
0884       mpcLinkGood->Fill(4 + offset, endcap * 5.5);
0885     if (CO.ME2n_9() == 0)
0886       mpcLinkGood->Fill(5 + offset, endcap * 5.5);
0887     if (CO.ME3n_3() == 0)
0888       mpcLinkGood->Fill(6 + offset, endcap * 5.5);
0889     if (CO.ME3n_9() == 0)
0890       mpcLinkGood->Fill(7 + offset, endcap * 5.5);
0891     if (CO.ME4n_3() == 0)
0892       mpcLinkGood->Fill(8 + offset, endcap * 5.5);
0893     if (CO.ME4n_9() == 0)
0894       mpcLinkGood->Fill(9 + offset, endcap * 5.5);
0895   }
0896 
0897   // Hits (CSC LCTs and RPC hits)
0898   edm::Handle<l1t::EMTFHitCollection> HitCollection;
0899   e.getByToken(hitToken, HitCollection);
0900 
0901   // Maps CSC station and ring to the monitor element index and uses symmetry of the endcaps
0902   const std::map<std::pair<int, int>, int> histIndexCSC = {{{1, 4}, 9},
0903                                                            {{1, 1}, 8},
0904                                                            {{1, 2}, 7},
0905                                                            {{1, 3}, 6},
0906                                                            {{2, 1}, 5},
0907                                                            {{2, 2}, 4},
0908                                                            {{3, 1}, 3},
0909                                                            {{3, 2}, 2},
0910                                                            {{4, 1}, 1},
0911                                                            {{4, 2}, 0}};
0912 
0913   // Maps CSC BX from -2 to 2 to monitor element cscLCTTIming
0914   const std::map<int, int> histIndexBX = {{0, 4}, {-1, 0}, {1, 1}, {-2, 2}, {2, 3}};
0915 
0916   // Maps RPC station and ring to the monitor element index and uses symmetry of the endcaps
0917   const std::map<std::pair<int, int>, int> histIndexRPC = {
0918       {{4, 3}, 0}, {{4, 2}, 1}, {{3, 3}, 2}, {{3, 2}, 3}, {{2, 2}, 4}, {{1, 2}, 5}};
0919 
0920   // Reverse 'rotate by 2' for RPC subsector
0921   auto get_subsector_rpc_cppf = [](int subsector_rpc) { return ((subsector_rpc + 3) % 6) + 1; };
0922 
0923   for (auto Hit = HitCollection->begin(); Hit != HitCollection->end(); ++Hit) {
0924     int endcap = Hit->Endcap();
0925     int sector = Hit->Sector();
0926     int station = Hit->Station();
0927     int ring = Hit->Ring();
0928     int cscid = Hit->CSC_ID();
0929     int chamber = Hit->Chamber();
0930     int strip = Hit->Strip();
0931     int wire = Hit->Wire();
0932     int cscid_offset = (sector - 1) * 9;
0933 
0934     int hist_index = 0;
0935     if (ring == 4 && strip >= 128)
0936       strip -= 128;
0937 
0938     if (Hit->Is_CSC() == true) {
0939       hist_index = histIndexCSC.at({station, ring});
0940       if (endcap > 0)
0941         hist_index = 19 - hist_index;
0942       cscLCTBX->Fill(Hit->BX(), hist_index);
0943       float evt_wgt = (Hit->Station() > 1 && Hit->Ring() == 1) ? 0.5 : 1.0;
0944       if (Hit->Neighbor() == false) {
0945         //Map for cscDQMOccupancy plot
0946         cscDQMOccupancy->Fill(chamber_bin(station, ring, chamber), hist_index, evt_wgt);
0947         if (station > 1 && (ring % 2) == 1) {
0948           cscDQMOccupancy->Fill(chamber_bin(station, ring, chamber) - 1, hist_index, evt_wgt);
0949         }
0950         cscLCTStrip[hist_index]->Fill(strip);
0951         cscLCTWire[hist_index]->Fill(wire);
0952         cscChamberStrip[hist_index]->Fill(chamber, strip);
0953         cscChamberWire[hist_index]->Fill(chamber, wire);
0954         if (Hit->Subsector() == 1) {
0955           cscLCTOccupancy->Fill(cscid + cscid_offset, endcap * (station - 0.5));
0956         } else {
0957           cscLCTOccupancy->Fill(cscid + cscid_offset, endcap * (station + 0.5));
0958         }
0959       } else {
0960         // Map neighbor chambers to "fake" CSC IDs: 1/3 --> 1, 1/6 --> 2, 1/9 --> 3, 2/3 --> 4, 2/9 --> 5, etc.
0961         int cscid_n = (station == 1 ? (cscid / 3) : (station * 2) + ((cscid - 3) / 6));
0962         cscLCTOccupancy->Fill(cscid_n + cscid_offset, endcap * 5.5);
0963       }
0964       if (Hit->Neighbor() == true) {
0965         cscDQMOccupancy->Fill((sector % 6 + 1) * 7 - 4, hist_index, evt_wgt);
0966       }
0967     }
0968 
0969     if (Hit->Is_RPC() == true) {
0970       hist_index = histIndexRPC.at({station, ring});
0971       if (endcap > 0)
0972         hist_index = 11 - hist_index;
0973 
0974       rpcHitBX->Fill(Hit->BX(), hist_index);
0975 
0976       if (Hit->Neighbor() == false) {
0977         rpcHitPhi[hist_index]->Fill(Hit->Phi_fp() / 4);
0978         rpcHitTheta[hist_index]->Fill(Hit->Theta_fp() / 4);
0979         rpcChamberPhi[hist_index]->Fill(chamber, Hit->Phi_fp() / 4);
0980         rpcChamberTheta[hist_index]->Fill(chamber, Hit->Theta_fp() / 4);
0981         rpcHitOccupancy->Fill((Hit->Sector_RPC() - 1) * 7 + get_subsector_rpc_cppf(Hit->Subsector_RPC()),
0982                               hist_index + 0.5);
0983       } else if (Hit->Neighbor() == true) {
0984         rpcHitOccupancy->Fill((Hit->Sector_RPC() - 1) * 7 + 7, hist_index + 0.5);
0985       }
0986     }
0987 
0988     // Add GEM Oct 27 2020
0989     hitTypeBX->Fill(4, Hit->BX());
0990     if (Hit->Is_CSC() == true)
0991       hitTypeBX->Fill(1, Hit->BX());
0992     else if (Hit->Is_RPC() == true)
0993       hitTypeBX->Fill(2, Hit->BX());
0994     else if (Hit->Is_GEM() == true)
0995       hitTypeBX->Fill(3, Hit->BX());
0996 
0997     if (Hit->Is_GEM() == true) {
0998       gemHitBX->Fill(Hit->BX(), (endcap > 0) ? 1.5 : 0.5);
0999       hist_index = (endcap > 0) ? 1 : 0;
1000       //Added def of layer
1001       int layer = Hit->Layer();
1002       int phi_part = Hit->Pad() / 64;  // 0-2
1003       int vfat = phi_part * 8 + Hit->Partition();
1004       if (Hit->Neighbor() == false) {
1005         gemChamberPad[hist_index]->Fill(chamber, Hit->Pad());
1006         gemChamberPartition[hist_index]->Fill(chamber, Hit->Partition());
1007         gemHitOccupancy->Fill(chamber_bin(1, 1, chamber), (endcap > 0) ? 1.5 : 0.5);  // follow CSC convention
1008         //Added plots 07-21-22 ***
1009         gemVFATBXPerChamber[chamber - 1][hist_index][layer]->Fill(Hit->BX(), vfat);
1010         //indexed plots by BX 07-21-22
1011         gemChamberVFATBX[hist_index][Hit->BX() + 3]->Fill(chamber_bin(1, 1, chamber), vfat);
1012       }
1013       //Added plots 06-07-22
1014 
1015       else {
1016         gemChamberPad[hist_index]->Fill((Hit->Sector() % 6) * 6 + 2, Hit->Pad());
1017         gemChamberPartition[hist_index]->Fill((Hit->Sector() % 6) * 6 + 2, Hit->Partition());
1018         gemHitOccupancy->Fill((Hit->Sector() % 6 + 1) * 7 - 4, (endcap > 0) ? 1.5 : 0.5);  // follow CSC convention
1019       }
1020     }  // End of if (Hit->Is_GEM() == true)
1021   }    // End of for (auto Hit = HitCollection->begin(); Hit != HitCollection->end(); ++Hit)
1022 
1023   // Tracks
1024   edm::Handle<l1t::EMTFTrackCollection> TrackCollection;
1025   e.getByToken(trackToken, TrackCollection);
1026 
1027   int nTracks = TrackCollection->size();
1028 
1029   emtfnTracks->Fill(std::min(nTracks, emtfnTracks->getTH1F()->GetNbinsX() - 1));
1030 
1031   constexpr int singleMuQuality = 12;
1032   constexpr float singleMuPT = 22;
1033   constexpr float singleMuUPT = 10;
1034 
1035   for (auto Track = TrackCollection->begin(); Track != TrackCollection->end(); ++Track) {
1036     int endcap = Track->Endcap();
1037     int sector = Track->Sector();
1038     float eta = Track->Eta();
1039     float phi_glob_rad = Track->Phi_glob() * M_PI / 180.;
1040     int mode = Track->Mode();
1041     int quality = Track->GMT_quality();
1042     int numHits = Track->NumHits();
1043     int modeNeighbor = Track->Mode_neighbor();
1044     int modeRPC = Track->Mode_RPC();
1045 
1046     // Only plot if there are <= 1 neighbor hits in the track to avoid spikes at sector boundaries
1047     if (modeNeighbor >= 2 && modeNeighbor != 4 && modeNeighbor != 8)
1048       continue;
1049 
1050     emtfTracknHits->Fill(numHits);
1051     emtfTrackBX->Fill(endcap * (sector - 0.5), Track->BX());
1052     emtfTrackPt->Fill(Track->Pt());
1053     emtfTrackDxy->Fill(Track->GMT_dxy());
1054     emtfTrackUnconstrainedPt->Fill(Track->Pt_dxy());
1055     emtfTrackEta->Fill(eta);
1056 
1057     emtfTrackOccupancy->Fill(eta, phi_glob_rad);
1058     emtfTrackMode->Fill(mode);
1059     emtfTrackQuality->Fill(quality);
1060     emtfTrackQualityVsMode->Fill(mode, quality);
1061     RPCvsEMTFTrackMode->Fill(mode, modeRPC);
1062     emtfTrackPhi->Fill(phi_glob_rad);
1063 
1064     if (quality >= singleMuQuality) {
1065       emtfTrackPtHighQuality->Fill(Track->Pt());
1066       emtfTrackUnconstrainedPtHighQuality->Fill(Track->Pt_dxy());
1067       emtfTrackEtaHighQuality->Fill(eta);
1068       emtfTrackPhiHighQuality->Fill(phi_glob_rad);
1069       emtfTrackOccupancyHighQuality->Fill(eta, phi_glob_rad);
1070       if (Track->Pt() >= singleMuPT) {
1071         emtfTrackPtHighQualityHighPT->Fill(Track->Pt());
1072         emtfTrackEtaHighQualityHighPT->Fill(eta);
1073         emtfTrackPhiHighQualityHighPT->Fill(phi_glob_rad);
1074         emtfTrackOccupancyHighQualityHighPT->Fill(eta, phi_glob_rad);
1075       }
1076       if (Track->Pt_dxy() >= singleMuUPT) {
1077         emtfTrackUnconstrainedPtHighQualityHighUPT->Fill(Track->Pt_dxy());
1078         emtfTrackEtaHighQualityHighUPT->Fill(eta);
1079         emtfTrackPhiHighQualityHighUPT->Fill(phi_glob_rad);
1080         emtfTrackOccupancyHighQualityHighUPT->Fill(eta, phi_glob_rad);
1081       }
1082     }
1083 
1084     ////////////////////////////////////////////////////
1085     ///  Begin block for CSC LCT and RPC hit timing  ///
1086     ////////////////////////////////////////////////////
1087     {
1088       // LCT and RPC Timing
1089       if (numHits < 2 || numHits > 4)
1090         continue;
1091       int numHitsInTrack_BX0 = 0;
1092       unsigned int hist_index2 = 4 - numHits;
1093 
1094       for (const auto& iTrkHit : Track->Hits()) {
1095         if (iTrkHit.Is_CSC()) {
1096           emtfTrackBXVsCSCLCT[hist_index2]->Fill(iTrkHit.BX(), Track->BX());
1097           int iCSC = (endcap > 0) ? (iTrkHit.Station() + 3) : (4 - iTrkHit.Station());
1098           emtfTrackModeVsCSCBXDiff[iCSC]->Fill(Track->BX() - iTrkHit.BX(),
1099                                                mode);  // Add mode vs BXdiff comparison Dec 07 2020
1100         } else if (iTrkHit.Is_RPC()) {
1101           emtfTrackBXVsRPCHit[hist_index2]->Fill(iTrkHit.BX(), Track->BX());
1102           int iRPC = (endcap > 0) ? (iTrkHit.Station() + 3) : (4 - iTrkHit.Station());
1103           emtfTrackModeVsRPCBXDiff[iRPC]->Fill(Track->BX() - iTrkHit.BX(),
1104                                                mode);  // Add mode vs BXdiff comparison Dec 07 2020
1105         } else if (iTrkHit.Is_GEM()) {
1106           emtfTrackBXVsGEMHit[hist_index2]->Fill(iTrkHit.BX(), Track->BX());
1107           int iGEM = (endcap > 0) ? 1 : 0;
1108           emtfTrackModeVsGEMBXDiff[iGEM]->Fill(Track->BX() - iTrkHit.BX(),
1109                                                mode);  // Add mode vs BXdiff comparison Dec 07 2020
1110         }
1111       }
1112 
1113       // Select well-timed tracks: >= 3 hits, with <= 1 in BX != 0
1114       if (numHits < 3)
1115         continue;
1116       for (const auto& jTrkHit : Track->Hits()) {
1117         if (jTrkHit.BX() == 0)
1118           numHitsInTrack_BX0++;
1119       }
1120       if (numHitsInTrack_BX0 < numHits - 1)
1121         continue;
1122 
1123       for (const auto& TrkHit : Track->Hits()) {
1124         int trackHitBX = TrkHit.BX();
1125         if (std::abs(trackHitBX) > 2)
1126           continue;  // Should never happen, but just to be safe ...
1127         //int cscid        = TrkHit.CSC_ID();
1128         int ring = TrkHit.Ring();
1129         int station = TrkHit.Station();
1130         int sector = TrkHit.Sector();
1131         //int cscid_offset = (sector - 1) * 9;//no longer needed after new time plots (maybe useful for future plots)
1132         int neighbor = TrkHit.Neighbor();
1133         int endcap = TrkHit.Endcap();
1134         int chamber = TrkHit.Chamber();
1135 
1136         int hist_index = 0;
1137         float evt_wgt = (TrkHit.Station() > 1 && TrkHit.Ring() == 1) ? 0.5 : 1.0;
1138 
1139         if (TrkHit.Is_CSC() == true) {
1140           hist_index = histIndexCSC.at({station, ring});
1141           if (endcap > 0)
1142             hist_index = 19 - hist_index;
1143           if (neighbor == false) {
1144             cscLCTTiming[histIndexBX.at(trackHitBX)]->Fill(chamber_bin(station, ring, chamber), hist_index, evt_wgt);
1145             cscTimingTot->Fill(chamber_bin(station, ring, chamber), hist_index, evt_wgt);
1146             if (station > 1 && (ring % 2) == 1) {
1147               cscLCTTiming[histIndexBX.at(trackHitBX)]->Fill(
1148                   chamber_bin(station, ring, chamber) - 1, hist_index, evt_wgt);
1149               cscTimingTot->Fill(chamber_bin(station, ring, chamber) - 1, hist_index, evt_wgt);
1150             }
1151           } else {
1152             // Map neighbor chambers to "fake" CSC IDs: 1/3 --> 1, 1/6 --> 2, 1/9 --> 3, 2/3 --> 4, 2/9 --> 5, etc.
1153             //int cscid_n = (station == 1 ? (cscid / 3) : (station * 2) + ((cscid - 3) / 6) );
1154             cscLCTTiming[histIndexBX.at(trackHitBX)]->Fill((sector % 6 + 1) * 7 - 4, hist_index, evt_wgt);
1155             cscTimingTot->Fill((sector % 6 + 1) * 7 - 4, hist_index, evt_wgt);
1156           }
1157 
1158           // Fill RPC timing with matched CSC LCTs
1159           if (trackHitBX == 0 && ring == 2) {
1160             for (auto Hit = HitCollection->begin(); Hit != HitCollection->end(); ++Hit) {
1161               if (Hit->Is_RPC() == false || neighbor == true)
1162                 continue;
1163               if (std::abs(Track->Eta() - Hit->Eta()) > 0.1)
1164                 continue;
1165               if (Hit->Endcap() != endcap || Hit->Station() != station || Hit->Chamber() != chamber)
1166                 continue;
1167               if (std::abs(Hit->BX()) > 2)
1168                 continue;
1169 
1170               hist_index = histIndexRPC.at({Hit->Station(), Hit->Ring()});
1171               if (Hit->Endcap() > 0)
1172                 hist_index = 11 - hist_index;
1173               rpcHitTimingInTrack->Fill(Hit->BX(), hist_index + 0.5);
1174               rpcHitTiming[histIndexBX.at(Hit->BX())]->Fill(
1175                   (Hit->Sector_RPC() - 1) * 7 + get_subsector_rpc_cppf(Hit->Subsector_RPC()), hist_index + 0.5);
1176               rpcHitTimingTot->Fill((Hit->Sector_RPC() - 1) * 7 + get_subsector_rpc_cppf(Hit->Subsector_RPC()),
1177                                     hist_index + 0.5);
1178             }  // End loop: for (auto Hit = HitCollection->begin(); Hit != HitCollection->end(); ++Hit)
1179           }    // End conditional: if (trackHitBX == 0 && ring == 2)
1180 
1181           // Fill GEM timing with matched CSC LCTs
1182           if (trackHitBX == 0 && station == 1 && ring == 1) {  // GEM only in station 1
1183             for (auto Hit = HitCollection->begin(); Hit != HitCollection->end(); ++Hit) {
1184               if (Hit->Is_GEM() == false)
1185                 continue;
1186               if (std::abs(Track->Eta() - Hit->Eta()) > 0.1)
1187                 continue;
1188               if (Hit->Endcap() != endcap || Hit->Station() != station || Hit->Chamber() != chamber ||
1189                   Hit->Neighbor() != neighbor)  //different neighbor requirement from RPC
1190                 continue;
1191               if (std::abs(Hit->BX()) > 2)
1192                 continue;
1193 
1194               if (neighbor == false) {
1195                 gemHitTiming[histIndexBX.at(Hit->BX())]->Fill(chamber_bin(1, 1, chamber), (endcap > 0) ? 1.5 : 0.5);
1196                 gemHitTimingTot->Fill(chamber_bin(1, 1, chamber), (endcap > 0) ? 1.5 : 0.5);
1197                 int ihist = (endcap > 0) ? 1 : 0;
1198                 gemHitPhi[ihist]->Fill(Hit->Phi_fp(), sector);
1199                 gemHitTheta[ihist]->Fill(Hit->Theta_fp(), sector);
1200                 gemHitVScscLCTPhi[ihist]->Fill(Hit->Phi_fp() - TrkHit.Phi_fp(), chamber);  // GEM vs CSC Dec 06 2020
1201                 gemHitVScscLCTTheta[ihist]->Fill(Hit->Theta_fp() - TrkHit.Theta_fp(), chamber);
1202                 gemHitVScscLCTBX[ihist]->Fill(Hit->BX() - TrkHit.BX(), chamber);
1203               } else {
1204                 gemHitTiming[histIndexBX.at(Hit->BX())]->Fill((sector % 6 + 1) * 7 - 4, (endcap > 0) ? 1.5 : 0.5);
1205                 gemHitTimingTot->Fill((sector % 6 + 1) * 7 - 4, (endcap > 0) ? 1.5 : 0.5);
1206               }
1207 
1208             }  // End loop: for (auto Hit = HitCollection->begin(); Hit != HitCollection->end(); ++Hit)
1209           }    // End conditional: if (trackHitBX == 0 && station == 1 && ring == 1)
1210         }      // End conditional: if (TrkHit.Is_CSC() == true)
1211 
1212         if (TrkHit.Is_RPC() == true && neighbor == false) {
1213           hist_index = histIndexRPC.at({station, ring});
1214           if (endcap > 0)
1215             hist_index = 11 - hist_index;
1216 
1217           rpcHitTimingInTrack->Fill(trackHitBX, hist_index + 0.5);
1218           rpcHitTiming[histIndexBX.at(trackHitBX)]->Fill(
1219               (TrkHit.Sector_RPC() - 1) * 7 + get_subsector_rpc_cppf(TrkHit.Subsector_RPC()), hist_index + 0.5);
1220           rpcHitTimingTot->Fill((TrkHit.Sector_RPC() - 1) * 7 + get_subsector_rpc_cppf(TrkHit.Subsector_RPC()),
1221                                 hist_index + 0.5);
1222         }  // End conditional: if (TrkHit.Is_RPC() == true && neighbor == false)
1223         if (TrkHit.Is_RPC() == true && neighbor == true) {
1224           hist_index = histIndexRPC.at({station, ring});
1225           if (endcap > 0)
1226             hist_index = 11 - hist_index;
1227           rpcHitTiming[histIndexBX.at(trackHitBX)]->Fill((TrkHit.Sector_RPC() - 1) * 7, hist_index + 0.5);
1228         }
1229         // Add GEM Timing Oct 27 2020
1230         if (TrkHit.Is_GEM() == true) {
1231           if (neighbor == false) {
1232             gemHitTiming[histIndexBX.at(trackHitBX)]->Fill(chamber_bin(1, 1, chamber), (endcap > 0) ? 1.5 : 0.5);
1233             gemHitTimingTot->Fill(chamber_bin(1, 1, chamber), (endcap > 0) ? 1.5 : 0.5);
1234           } else {
1235             gemHitTiming[histIndexBX.at(trackHitBX)]->Fill((sector % 6 + 1) * 7 - 4, (endcap > 0) ? 1.5 : 0.5);
1236             gemHitTimingTot->Fill((sector % 6 + 1) * 7 - 4, (endcap > 0) ? 1.5 : 0.5);
1237           }
1238         }  // End condition: if (TrkHit.Is_GEM() == true)
1239       }    // End loop: for (int iHit = 0; iHit < numHits; ++iHit)
1240     }
1241     //////////////////////////////////////////////////
1242     ///  End block for CSC LCT and RPC hit timing  ///
1243     //////////////////////////////////////////////////
1244 
1245   }  // End loop: for (auto Track = TrackCollection->begin(); Track != TrackCollection->end(); ++Track)
1246 
1247   // CSC LCT and RPC Hit Timing Efficieny
1248   for (int hist_index = 0; hist_index < 5; ++hist_index) {
1249     cscLCTTimingFrac[hist_index]->getTH2F()->Divide(cscLCTTiming[hist_index]->getTH2F(), cscTimingTot->getTH2F());
1250     rpcHitTimingFrac[hist_index]->getTH2F()->Divide(rpcHitTiming[hist_index]->getTH2F(), rpcHitTimingTot->getTH2F());
1251     gemHitTimingFrac[hist_index]->getTH2F()->Divide(gemHitTiming[hist_index]->getTH2F(), gemHitTimingTot->getTH2F());
1252   }
1253 
1254   // Regional Muon Candidates
1255   edm::Handle<l1t::RegionalMuonCandBxCollection> MuonBxCollection;
1256   e.getByToken(muonToken, MuonBxCollection);
1257 
1258   for (int itBX = MuonBxCollection->getFirstBX(); itBX <= MuonBxCollection->getLastBX(); ++itBX) {
1259     for (l1t::RegionalMuonCandBxCollection::const_iterator Muon = MuonBxCollection->begin(itBX);
1260          Muon != MuonBxCollection->end(itBX);
1261          ++Muon) {
1262       emtfMuonBX->Fill(itBX);
1263       emtfMuonhwPt->Fill(Muon->hwPt());
1264       emtfMuonhwPtUnconstrained->Fill(Muon->hwPtUnconstrained());
1265       emtfMuonhwDxy->Fill(Muon->hwDXY());
1266       emtfMuonhwEta->Fill(Muon->hwEta());
1267       emtfMuonhwPhi->Fill(Muon->hwPhi());
1268       emtfMuonhwQual->Fill(Muon->hwQual());
1269     }
1270   }
1271 }