Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:13:21

0001 #include "L1Trigger/RegionalCaloTrigger/interface/L1RCTProducer.h"
0002 
0003 // RunInfo stuff
0004 #include "CondFormats/RunInfo/interface/RunInfo.h"
0005 #include "FWCore/Framework/interface/LuminosityBlock.h"
0006 #include "FWCore/Framework/interface/Run.h"
0007 
0008 #include <vector>
0009 using std::vector;
0010 #include <iostream>
0011 
0012 using std::cout;
0013 using std::endl;
0014 namespace {
0015   constexpr int crateFED[18][6] = {{613, 614, 603, 702, 718, 1118},
0016                                    {611, 612, 602, 700, 718, 1118},
0017                                    {627, 610, 601, 716, 722, 1122},
0018                                    {625, 626, 609, 714, 722, 1122},
0019                                    {623, 624, 608, 712, 722, 1122},
0020                                    {621, 622, 607, 710, 720, 1120},
0021                                    {619, 620, 606, 708, 720, 1120},
0022                                    {617, 618, 605, 706, 720, 1120},
0023                                    {615, 616, 604, 704, 718, 1118},
0024                                    {631, 632, 648, 703, 719, 1118},
0025                                    {629, 630, 647, 701, 719, 1118},
0026                                    {645, 628, 646, 717, 723, 1122},
0027                                    {643, 644, 654, 715, 723, 1122},
0028                                    {641, 642, 653, 713, 723, 1122},
0029                                    {639, 640, 652, 711, 721, 1120},
0030                                    {637, 638, 651, 709, 721, 1120},
0031                                    {635, 636, 650, 707, 721, 1120},
0032                                    {633, 634, 649, 705, 719, 1118}};
0033 }
0034 L1RCTProducer::L1RCTProducer(const edm::ParameterSet &conf)
0035     : rctLookupTables(new L1RCTLookupTables),
0036       rct(new L1RCT(rctLookupTables.get())),
0037       useEcal(conf.getParameter<bool>("useEcal")),
0038       useHcal(conf.getParameter<bool>("useHcal")),
0039       ecalDigis(conf.getParameter<std::vector<edm::InputTag>>("ecalDigis")),
0040       hcalDigis(conf.getParameter<std::vector<edm::InputTag>>("hcalDigis")),
0041       bunchCrossings(conf.getParameter<std::vector<int>>("BunchCrossings")),
0042       getFedsFromOmds(conf.getParameter<bool>("getFedsFromOmds")),
0043       queryDelayInLS(conf.getParameter<unsigned int>("queryDelayInLS")),
0044       queryIntervalInLS(conf.getParameter<unsigned int>("queryIntervalInLS")),
0045       conditionsLabel(conf.getParameter<std::string>("conditionsLabel")),
0046       fedUpdatedMask(nullptr),
0047 
0048       rctParamsToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", conditionsLabel))),
0049       emScaleToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", conditionsLabel))),
0050       ecalScaleToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", conditionsLabel))),
0051       hcalScaleToken_(esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", conditionsLabel))),
0052 
0053       beginRunRunInfoToken_(esConsumes<edm::Transition::BeginRun>()),
0054 
0055       beginRunChannelMaskToken_(esConsumes<edm::Transition::BeginRun>()),
0056       beginRunHotChannelMaskToken_(esConsumes<edm::Transition::BeginRun>()) {
0057   produces<L1CaloEmCollection>();
0058   produces<L1CaloRegionCollection>();
0059 
0060   if (getFedsFromOmds) {
0061     beginLumiChannelMaskToken_ = esConsumes<edm::Transition::BeginLuminosityBlock>();
0062     beginLumiHotChannelMaskToken_ = esConsumes<edm::Transition::BeginLuminosityBlock>();
0063     beginLumiRunInfoToken_ = esConsumes<edm::Transition::BeginLuminosityBlock>();
0064     omdsRunInfoToken_ = esConsumes<edm::Transition::BeginLuminosityBlock>(edm::ESInputTag("", "OmdsFedVector"));
0065   }
0066 
0067   for (unsigned int ihc = 0; ihc < hcalDigis.size(); ihc++) {
0068     consumes<edm::SortedCollection<HcalTriggerPrimitiveDigi, edm::StrictWeakOrdering<HcalTriggerPrimitiveDigi>>>(
0069         hcalDigis[ihc]);
0070   }
0071 
0072   for (unsigned int iec = 0; iec < ecalDigis.size(); iec++) {
0073     consumes<edm::SortedCollection<EcalTriggerPrimitiveDigi, edm::StrictWeakOrdering<EcalTriggerPrimitiveDigi>>>(
0074         ecalDigis[iec]);
0075   }
0076 }
0077 
0078 void L1RCTProducer::beginRun(edm::Run const &run, const edm::EventSetup &eventSetup) {
0079   //  std::cout << "getFedsFromOmds is " << getFedsFromOmds << std::endl;
0080 
0081   updateConfiguration(eventSetup);
0082 
0083   // list of RCT channels to mask
0084   L1RCTChannelMask const &channelMask = eventSetup.getData(beginRunChannelMaskToken_);
0085 
0086   // list of Noisy RCT channels to mask
0087   L1RCTNoisyChannelMask const &hotChannelMask = eventSetup.getData(beginRunHotChannelMaskToken_);
0088 
0089   updateFedVector(channelMask, hotChannelMask, getFedVectorFromRunInfo(beginRunRunInfoToken_, eventSetup));
0090 }
0091 
0092 void L1RCTProducer::beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &context) {
0093   // check LS number every LS, if the checkOMDS flag is set AND it's the right
0094   // LS, update the FED vector from OMDS can pass the flag as the bool??  but
0095   // only check LS number if flag is true anyhow
0096   if (getFedsFromOmds) {
0097     unsigned int nLumi = lumiSeg.luminosityBlock();  // doesn't even need the (unsigned int) cast
0098                                                      // because LuminosityBlockNumber_t is already
0099                                                      // an unsigned int
0100     // LS count starts at 1, want to be able to delay 0 LS's intuitively
0101     if (((nLumi - 1) == queryDelayInLS) ||
0102         (queryIntervalInLS > 0 &&
0103          nLumi % queryIntervalInLS == 0))  // to guard against problems if online DQM crashes; every 100
0104                                            // LS is ~20-30 minutes, not too big a load, hopefully not too
0105                                            // long between
0106     {
0107       //      std::cout << "Lumi section for this FED vector update is " <<
0108       // nLumi << std::endl;
0109 
0110       // list of RCT channels to mask
0111       L1RCTChannelMask const &channelMask = context.getData(beginLumiChannelMaskToken_);
0112 
0113       // list of Noisy RCT channels to mask
0114       L1RCTNoisyChannelMask const &hotChannelMask = context.getData(beginLumiHotChannelMaskToken_);
0115 
0116       updateFedVector(channelMask, hotChannelMask, getFedVectorFromOmds(context));
0117     } else if (queryIntervalInLS <= 0) {
0118       // don't do interval checking... cout message??
0119     }
0120   }
0121 }
0122 
0123 void L1RCTProducer::updateConfiguration(const edm::EventSetup &eventSetup) {
0124   // Refresh configuration information every event
0125   // Hopefully, this does not take too much time
0126   // There should be a call back function in future to
0127   // handle changes in configuration
0128   // parameters to configure RCT (thresholds, etc)
0129   const L1RCTParameters *r = &eventSetup.getData(rctParamsToken_);
0130 
0131   // SCALES
0132 
0133   // energy scale to convert eGamma output
0134   const L1CaloEtScale *s = &eventSetup.getData(emScaleToken_);
0135 
0136   // get energy scale to convert input from ECAL
0137   const L1CaloEcalScale *e = &eventSetup.getData(ecalScaleToken_);
0138 
0139   // get energy scale to convert input from HCAL
0140   const L1CaloHcalScale *h = &eventSetup.getData(hcalScaleToken_);
0141 
0142   // set scales
0143   rctLookupTables->setEcalScale(e);
0144   rctLookupTables->setHcalScale(h);
0145 
0146   rctLookupTables->setRCTParameters(r);
0147   rctLookupTables->setL1CaloEtScale(s);
0148 }
0149 
0150 void L1RCTProducer::updateFedVector(const L1RCTChannelMask &channelMask,
0151                                     const L1RCTNoisyChannelMask &hotChannelMask,
0152                                     const std::vector<int> &Feds)
0153 // http://cmslxr.fnal.gov/lxr/source/FWCore/Framework/interface/EventSetup.h
0154 {
0155   rctLookupTables->setNoisyChannelMask(&hotChannelMask);
0156 
0157   // Update the channel mask according to the FED VECTOR
0158   // This is the beginning of run. We delete the old
0159   // create the new and set it in the LUTs
0160 
0161   fedUpdatedMask = std::make_unique<L1RCTChannelMask>();
0162   // copy a constant object
0163   for (int i = 0; i < 18; i++) {
0164     for (int j = 0; j < 2; j++) {
0165       for (int k = 0; k < 28; k++) {
0166         fedUpdatedMask->ecalMask[i][j][k] = channelMask.ecalMask[i][j][k];
0167         fedUpdatedMask->hcalMask[i][j][k] = channelMask.hcalMask[i][j][k];
0168       }
0169       for (int k = 0; k < 4; k++) {
0170         fedUpdatedMask->hfMask[i][j][k] = channelMask.hfMask[i][j][k];
0171       }
0172     }
0173   }
0174 
0175   // so can create/initialize/assign const quantity in one line accounting for
0176   // if statement wikipedia says this is exactly what it's for:
0177   // http://en.wikipedia.org/wiki/%3F:#C.2B.2B
0178 
0179   //   std::cout << "Contents of ";
0180   //   std::cout << (getFromOmds ? "OMDS RunInfo" : "standard RunInfo");
0181   //   std::cout << " FED vector" << std::endl;
0182   //   printFedVector(Feds);
0183 
0184   bool useUpgradedHF = false;
0185 
0186   std::vector<int> caloFeds;  // pare down the feds to the interesting ones
0187   // is this unneccesary?
0188   // Mike B : This will decrease the find speed so better do it
0189   for (std::vector<int>::const_iterator cf = Feds.begin(); cf != Feds.end(); ++cf) {
0190     int fedNum = *cf;
0191     if ((fedNum > 600 && fedNum < 724) || fedNum == 1118 || fedNum == 1120 || fedNum == 1122)
0192       caloFeds.push_back(fedNum);
0193 
0194     if (fedNum == 1118 || fedNum == 1120 || fedNum == 1122)
0195       useUpgradedHF = true;
0196   }
0197 
0198   for (int cr = 0; cr < 18; ++cr) {
0199     for (crateSection cs = c_min; cs <= c_max; cs = crateSection(cs + 1)) {
0200       bool fedFound = false;
0201 
0202       // Try to find the FED
0203       std::vector<int>::iterator fv = std::find(caloFeds.begin(), caloFeds.end(), crateFED[cr][cs]);
0204       if (fv != caloFeds.end())
0205         fedFound = true;
0206 
0207       if (!fedFound) {
0208         int eta_min = 0;
0209         int eta_max = 0;
0210         bool phi_even[2] = {false};  //, phi_odd = false;
0211         bool ecal = false;
0212 
0213         switch (cs) {
0214           case ebEvenFed:
0215             eta_min = minBarrel;
0216             eta_max = maxBarrel;
0217             phi_even[0] = true;
0218             ecal = true;
0219             break;
0220 
0221           case ebOddFed:
0222             eta_min = minBarrel;
0223             eta_max = maxBarrel;
0224             phi_even[1] = true;
0225             ecal = true;
0226             break;
0227 
0228           case eeFed:
0229             eta_min = minEndcap;
0230             eta_max = maxEndcap;
0231             phi_even[0] = true;
0232             phi_even[1] = true;
0233             ecal = true;
0234             break;
0235 
0236           case hbheFed:
0237             eta_min = minBarrel;
0238             eta_max = maxEndcap;
0239             phi_even[0] = true;
0240             phi_even[1] = true;
0241             ecal = false;
0242             break;
0243 
0244           case hfFed:
0245             if (useUpgradedHF)
0246               break;
0247 
0248             eta_min = minHF;
0249             eta_max = maxHF;
0250 
0251             phi_even[0] = true;
0252             phi_even[1] = true;
0253             ecal = false;
0254             break;
0255 
0256           case hfFedUp:
0257             if (!useUpgradedHF)
0258               break;
0259 
0260             eta_min = minHF;
0261             eta_max = maxHF;
0262 
0263             phi_even[0] = true;
0264             phi_even[1] = true;
0265             ecal = false;
0266             break;
0267 
0268           default:
0269             break;
0270         }
0271         for (int ieta = eta_min; ieta <= eta_max; ++ieta) {
0272           if (ieta <= 28)  // barrel and endcap
0273             for (int even = 0; even <= 1; even++) {
0274               if (phi_even[even]) {
0275                 if (ecal)
0276                   fedUpdatedMask->ecalMask[cr][even][ieta - 1] = true;
0277                 else
0278                   fedUpdatedMask->hcalMask[cr][even][ieta - 1] = true;
0279               }
0280             }
0281           else
0282             for (int even = 0; even <= 1; even++)
0283               if (phi_even[even])
0284                 fedUpdatedMask->hfMask[cr][even][ieta - 29] = true;
0285         }
0286       }
0287     }
0288   }
0289 
0290   rctLookupTables->setChannelMask(fedUpdatedMask.get());
0291 }
0292 
0293 const std::vector<int> L1RCTProducer::getFedVectorFromRunInfo(const edm::ESGetToken<RunInfo, RunInfoRcd> &token,
0294                                                               const edm::EventSetup &eventSetup) const {
0295   //  std::cout << "Getting FED vector from standard RunInfo object" <<
0296   //  std::endl;
0297   // get FULL FED vector from RUNINFO
0298   return eventSetup.getData(token).m_fed_in;
0299 }
0300 
0301 const std::vector<int> L1RCTProducer::getFedVectorFromOmds(const edm::EventSetup &eventSetup) const {
0302   //  std::cout << "Getting FED vector from my specific ES RunInfo object" <<
0303   //  std::endl;
0304 
0305   // get FULL FED vector from RunInfo object specifically created to have OMDS
0306   // fed vector
0307   edm::ESHandle<RunInfo> sum = eventSetup.getHandle(omdsRunInfoToken_);
0308   if (sum.isValid()) {
0309     return sum->m_fed_in;
0310   } else {
0311     return getFedVectorFromRunInfo(beginLumiRunInfoToken_, eventSetup);
0312   }
0313 }
0314 
0315 void L1RCTProducer::produce(edm::Event &event, const edm::EventSetup &eventSetup) {
0316   std::unique_ptr<L1CaloEmCollection> rctEmCands(new L1CaloEmCollection);
0317   std::unique_ptr<L1CaloRegionCollection> rctRegions(new L1CaloRegionCollection);
0318 
0319   if (!(ecalDigis.size() == hcalDigis.size() && hcalDigis.size() == bunchCrossings.size()))
0320     throw cms::Exception("BadInput") << "From what I see the number of your your ECAL input digi "
0321                                         "collections.\n"
0322                                      << "is different from the size of your HCAL digi input collections\n"
0323                                      << "or the size of your BX factor collection"
0324                                      << "They must be the same to correspond to the same Bxs\n"
0325                                      << "It does not matter if one of them is empty\n";
0326 
0327   // loop through and process each bx
0328   for (unsigned short sample = 0; sample < bunchCrossings.size(); sample++) {
0329     edm::Handle<EcalTrigPrimDigiCollection> ecal;
0330     edm::Handle<HcalTrigPrimDigiCollection> hcal;
0331 
0332     EcalTrigPrimDigiCollection ecalIn;
0333     HcalTrigPrimDigiCollection hcalIn;
0334 
0335     if (useHcal && event.getByLabel(hcalDigis[sample], hcal))
0336       hcalIn = *hcal;
0337 
0338     if (useEcal && event.getByLabel(ecalDigis[sample], ecal))
0339       ecalIn = *ecal;
0340 
0341     rct->digiInput(ecalIn, hcalIn);
0342     rct->processEvent();
0343 
0344     // Stuff to create
0345     for (int j = 0; j < 18; j++) {
0346       L1CaloEmCollection isolatedEGObjects = rct->getIsolatedEGObjects(j);
0347       L1CaloEmCollection nonisolatedEGObjects = rct->getNonisolatedEGObjects(j);
0348       for (int i = 0; i < 4; i++) {
0349         isolatedEGObjects.at(i).setBx(bunchCrossings[sample]);
0350         nonisolatedEGObjects.at(i).setBx(bunchCrossings[sample]);
0351         rctEmCands->push_back(isolatedEGObjects.at(i));
0352         rctEmCands->push_back(nonisolatedEGObjects.at(i));
0353       }
0354     }
0355 
0356     for (int i = 0; i < 18; i++) {
0357       std::vector<L1CaloRegion> regions = rct->getRegions(i);
0358       for (int j = 0; j < 22; j++) {
0359         regions.at(j).setBx(bunchCrossings[sample]);
0360         rctRegions->push_back(regions.at(j));
0361       }
0362     }
0363   }
0364 
0365   // putting stuff back into event
0366   event.put(std::move(rctEmCands));
0367   event.put(std::move(rctRegions));
0368 }
0369 
0370 // print contents of (FULL) FED vector
0371 void L1RCTProducer::printFedVector(const std::vector<int> &fedVector) {
0372   std::cout << "Contents of given fedVector: ";
0373   std::copy(fedVector.begin(), fedVector.end(), std::ostream_iterator<int>(std::cout, ", "));
0374   std::cout << std::endl;
0375 }
0376 
0377 // print contents of RCT channel mask fedUpdatedMask
0378 void L1RCTProducer::printUpdatedFedMask() {
0379   if (fedUpdatedMask != nullptr) {
0380     fedUpdatedMask->print(std::cout);
0381   } else {
0382     std::cout << "Trying to print contents of fedUpdatedMask, but it doesn't exist!" << std::endl;
0383   }
0384 }
0385 
0386 // print contents of RCT channel mask fedUpdatedMask
0387 void L1RCTProducer::printUpdatedFedMaskVerbose() {
0388   if (fedUpdatedMask != nullptr) {
0389     // print contents of fedvector
0390     std::cout << "Contents of fedUpdatedMask: ";
0391     //       std::copy(fedUpdatedMask.begin(), fedUpdatedMask.end(),
0392     //       std::ostream_iterator<int>(std::cout, ", "));
0393     std::cout << "--> ECAL mask: " << std::endl;
0394     for (int i = 0; i < 18; i++) {
0395       for (int j = 0; j < 2; j++) {
0396         for (int k = 0; k < 28; k++) {
0397           std::cout << fedUpdatedMask->ecalMask[i][j][k] << ", ";
0398         }
0399       }
0400     }
0401     std::cout << "--> HCAL mask: " << std::endl;
0402     for (int i = 0; i < 18; i++) {
0403       for (int j = 0; j < 2; j++) {
0404         for (int k = 0; k < 28; k++) {
0405           std::cout << fedUpdatedMask->hcalMask[i][j][k] << ", ";
0406         }
0407       }
0408     }
0409     std::cout << "--> HF mask: " << std::endl;
0410     for (int i = 0; i < 18; i++) {
0411       for (int j = 0; j < 2; j++) {
0412         for (int k = 0; k < 4; k++) {
0413           std::cout << fedUpdatedMask->hfMask[i][j][k] << ", ";
0414         }
0415       }
0416     }
0417 
0418     std::cout << std::endl;
0419   } else {
0420     // print error message
0421     std::cout << "Trying to print contents of fedUpdatedMask, but it doesn't exist!" << std::endl;
0422   }
0423 }