Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-08-30 02:10:34

0001 #include "DQM/L1TMonitor/interface/L1TCaloLayer1Summary.h"
0002 
0003 L1TCaloLayer1Summary::L1TCaloLayer1Summary(const edm::ParameterSet& iConfig)
0004     : caloLayer1CICADAScoreToken_(
0005           consumes<l1t::CICADABxCollection>(iConfig.getParameter<edm::InputTag>("caloLayer1CICADAScore"))),
0006       gtCICADAScoreToken_(consumes<l1t::CICADABxCollection>(iConfig.getParameter<edm::InputTag>("gtCICADAScore"))),
0007       simCICADAScoreToken_(consumes<l1t::CICADABxCollection>(iConfig.getParameter<edm::InputTag>("simCICADAScore"))),
0008       caloLayer1RegionsToken_(
0009           consumes<L1CaloRegionCollection>(iConfig.getParameter<edm::InputTag>("caloLayer1Regions"))),
0010       simRegionsToken_(consumes<L1CaloRegionCollection>(iConfig.getParameter<edm::InputTag>("simRegions"))),
0011       fedRawData_(consumes<FEDRawDataCollection>(iConfig.getParameter<edm::InputTag>("fedRawDataLabel"))),
0012       histFolder_(iConfig.getParameter<std::string>("histFolder")) {}
0013 
0014 // ------------ method called for each event  ------------
0015 void L1TCaloLayer1Summary::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0016   edm::Handle<FEDRawDataCollection> fedRawDataCollection;
0017   iEvent.getByToken(fedRawData_, fedRawDataCollection);
0018   if (fedRawDataCollection.isValid()) {
0019     for (int iFed = FEDNumbering::MINRCTFEDID + 4; iFed < FEDNumbering::MAXRCTFEDID; iFed += 2) {
0020       const FEDRawData& fedRawData = fedRawDataCollection->FEDData(iFed);
0021       if (fedRawData.size() == 0) {
0022         continue;
0023       }
0024       const uint64_t* fedRawDataArray = (const uint64_t*)fedRawData.data();
0025       UCTDAQRawData daqData(fedRawDataArray);
0026 
0027       if (daqData.nAMCs() == 7) {
0028         UCTAMCRawData amcSlot7(daqData.amcPayload(3));
0029         if (amcSlot7.amcNo() == 7) {
0030           histoSlot7MinusDaqBxid->Fill(amcSlot7.BXID() - daqData.BXID());
0031         }
0032       }
0033     }
0034   }
0035 
0036   L1CaloRegionCollection caloLayer1Regions = iEvent.get(caloLayer1RegionsToken_);
0037   L1CaloRegionCollection simRegions = iEvent.get(simRegionsToken_);
0038   int nRegions = caloLayer1Regions.size();
0039 
0040   unsigned int maxEtaIdx = 0;
0041   for (int iRegion = 0; iRegion < nRegions; iRegion++) {
0042     if (maxEtaIdx < caloLayer1Regions[iRegion].gctEta()) {
0043       maxEtaIdx = caloLayer1Regions[iRegion].gctEta();
0044     }
0045   }
0046   int matrixSize = maxEtaIdx + 1;
0047 
0048   bool foundMatrix[2][matrixSize][matrixSize];
0049   int etMatrix[2][matrixSize][matrixSize];
0050   for (int i = 0; i < 2; i++) {
0051     for (int j = 0; j < matrixSize; j++) {
0052       for (int k = 0; k < matrixSize; k++) {
0053         foundMatrix[i][j][k] = false;
0054         etMatrix[i][j][k] = 0;
0055       }
0056     }
0057   }
0058 
0059   for (int iRegion = 0; iRegion < nRegions; iRegion++) {
0060     L1CaloRegion cRegion = caloLayer1Regions[iRegion];
0061     L1CaloRegion sRegion = simRegions[iRegion];
0062 
0063     foundMatrix[0][cRegion.gctEta()][cRegion.gctPhi()] = true;
0064     etMatrix[0][cRegion.gctEta()][cRegion.gctPhi()] = cRegion.et();
0065     foundMatrix[1][sRegion.gctEta()][sRegion.gctPhi()] = true;
0066     etMatrix[1][sRegion.gctEta()][sRegion.gctPhi()] = sRegion.et();
0067   }
0068   int iRegion = 0;
0069   for (int iEta = 0; iEta < matrixSize; iEta++) {
0070     for (int iPhi = 0; iPhi < matrixSize; iPhi++) {
0071       if (foundMatrix[0][iEta][iPhi] && foundMatrix[1][iEta][iPhi]) {
0072         histoCaloRegions->Fill(iRegion, etMatrix[0][iEta][iPhi]);
0073         histoSimRegions->Fill(iRegion, etMatrix[1][iEta][iPhi]);
0074         histoCaloMinusSimRegions->Fill(iRegion, etMatrix[0][iEta][iPhi] - etMatrix[1][iEta][iPhi]);
0075         iRegion++;
0076       }
0077     }
0078   }
0079 
0080   auto caloCICADAScores = iEvent.get(caloLayer1CICADAScoreToken_);
0081   auto gtCICADAScores = iEvent.get(gtCICADAScoreToken_);
0082   auto simCICADAScores = iEvent.get(simCICADAScoreToken_);
0083 
0084   if (caloCICADAScores.size() > 0) {
0085     histoCaloLayer1CICADAScore->Fill(caloCICADAScores[0]);
0086     if (gtCICADAScores.size() > 0) {
0087       histoGtCICADAScore->Fill(gtCICADAScores.at(0, 0));
0088       histoCaloMinusGt->Fill(caloCICADAScores[0] - gtCICADAScores.at(0, 0));
0089     }
0090     if (simCICADAScores.size() > 0) {
0091       histoSimCICADAScore->Fill(simCICADAScores[0]);
0092       histoCaloMinusSim->Fill(caloCICADAScores[0] - simCICADAScores[0]);
0093     }
0094   }
0095 }
0096 
0097 void L1TCaloLayer1Summary::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const&) {
0098   ibooker.setCurrentFolder(histFolder_);
0099   histoSlot7MinusDaqBxid = ibooker.book1D("slot7BXID", "Slot 7- DAQ BXID", 50, -20, 20);
0100 
0101   ibooker.setCurrentFolder(histFolder_ + "/CICADAScore");
0102   histoCaloLayer1CICADAScore = ibooker.book1D("caloLayer1CICADAScore", "CaloLayer1 CICADAScore", 50, 0, 200);
0103   histoGtCICADAScore = ibooker.book1D("gtCICADAScore", "GT CICADAScore at BX0", 50, 0, 200);
0104   histoCaloMinusGt = ibooker.book1D("caloMinusGtCICADAScore", "CaloLayer1 - GT CICADAScore at BX0", 50, -50, 50);
0105   histoSimCICADAScore =
0106       ibooker.book1D("simCaloLayer1CICADAScore", "simCaloLayer1 CICADAScore (input: DAQ regions)", 50, 0, 200);
0107   histoCaloMinusSim = ibooker.book1D(
0108       "caloMinusSimCICADAScore", "CaloLayer1 - simCaloLayer1 (input: DAQ regions) CICADAScore", 50, -50, 50);
0109 
0110   ibooker.setCurrentFolder(histFolder_ + "/Regions");
0111   histoCaloMinusSimRegions =
0112       ibooker.book2D("caloMinusSumRegions",
0113                      "CaloLayer1 - simCaloLayer1 (input: DAQ trigger primatives) Regions;Region;ET Difference",
0114                      252,
0115                      -0.5,
0116                      252.5,
0117                      100,
0118                      -400,
0119                      400);
0120   histoCaloRegions = ibooker.book2D("caloLayer1Regions", "CaloLayer1 Regions;Region;ET", 252, -0.5, 252.5, 100, 0, 800);
0121   histoSimRegions = ibooker.book2D("simCaloLayer1Regions",
0122                                    "simCaloLayer1 Regions (input: DAQ trigger primatives);Region;ET",
0123                                    252,
0124                                    -0.5,
0125                                    252.5,
0126                                    100,
0127                                    0,
0128                                    800);
0129 }
0130 
0131 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0132 void L1TCaloLayer1Summary::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0133   // l1tCaloLayer1Summary
0134   edm::ParameterSetDescription desc;
0135   desc.add<edm::InputTag>("caloLayer1CICADAScore", edm::InputTag("caloLayer1Digis", "CICADAScore"));
0136   desc.add<edm::InputTag>("gtCICADAScore", edm::InputTag("gtTestcrateStage2Digis", "CICADAScore"));
0137   desc.add<edm::InputTag>("simCICADAScore", edm::InputTag("simCaloStage2Layer1Summary", "CICADAScore"));
0138   desc.add<edm::InputTag>("caloLayer1Regions", edm::InputTag("caloLayer1Digis"));
0139   desc.add<edm::InputTag>("simRegions", edm::InputTag("simCaloStage2Layer1Digis"));
0140   desc.add<edm::InputTag>("fedRawDataLabel", edm::InputTag("rawDataCollector"));
0141   desc.add<std::string>("histFolder", "L1T/L1TCaloLayer1Summary");
0142   descriptions.add("l1tCaloLayer1Summary", desc);
0143 }