Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-08-12 02:29:04

0001 #include <memory>
0002 
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/stream/EDProducer.h"
0005 #include "FWCore/Framework/interface/MakerMacros.h"
0006 
0007 #include "DataFormats/L1Scouting/interface/OrbitCollection.h"
0008 #include "DataFormats/NanoAOD/interface/OrbitFlatTable.h"
0009 #include "DataFormats/L1Scouting/interface/L1ScoutingCalo.h"
0010 #include "DataFormats/L1Trigger/interface/EtSum.h"
0011 #include "L1TriggerScouting/Utilities/interface/conversion.h"
0012 
0013 class L1ScoutingEtSumOrbitFlatTableProducer : public edm::stream::EDProducer<> {
0014 public:
0015   explicit L1ScoutingEtSumOrbitFlatTableProducer(const edm::ParameterSet&);
0016   ~L1ScoutingEtSumOrbitFlatTableProducer() override = default;
0017 
0018   static void fillDescriptions(edm::ConfigurationDescriptions&);
0019 
0020   void produce(edm::Event&, edm::EventSetup const&) override;
0021 
0022 private:
0023   std::unique_ptr<l1ScoutingRun3::OrbitFlatTable> produceSingle(l1ScoutingRun3::BxSumsOrbitCollection const&) const;
0024   std::unique_ptr<l1ScoutingRun3::OrbitFlatTable> produceMultiple(l1ScoutingRun3::BxSumsOrbitCollection const&) const;
0025 
0026   edm::EDGetTokenT<l1ScoutingRun3::BxSumsOrbitCollection> src_;
0027 
0028   std::string name_;
0029   std::string doc_;
0030   bool singleton_;
0031   bool writePhysicalValues_;
0032   bool writeHardwareValues_;
0033   bool writeHF_;
0034   bool writeAsym_;
0035   bool writeMinBias_;
0036   bool writeTowerCount_;
0037   bool writeCentrality_;
0038   int ptPrecision_;
0039   int phiPrecision_;
0040 };
0041 
0042 L1ScoutingEtSumOrbitFlatTableProducer::L1ScoutingEtSumOrbitFlatTableProducer(const edm::ParameterSet& params)
0043     : src_(consumes<OrbitCollection<l1ScoutingRun3::BxSums>>(params.getParameter<edm::InputTag>("src"))),
0044       name_(params.getParameter<std::string>("name")),
0045       doc_(params.getParameter<std::string>("doc")),
0046       singleton_(params.getParameter<bool>("singleton")),
0047       writePhysicalValues_(params.getParameter<bool>("writePhysicalValues")),
0048       writeHardwareValues_(params.getParameter<bool>("writeHardwareValues")),
0049       writeHF_(params.getParameter<bool>("writeHF")),
0050       writeAsym_(params.getParameter<bool>("writeAsym")),
0051       writeMinBias_(params.getParameter<bool>("writeMinBias")),
0052       writeTowerCount_(params.getParameter<bool>("writeTowerCount")),
0053       writeCentrality_(params.getParameter<bool>("writeCentrality")),
0054       ptPrecision_(params.getParameter<int>("ptPrecision")),
0055       phiPrecision_(params.getParameter<int>("phiPrecision")) {
0056   if (!writePhysicalValues_ && !writeHardwareValues_) {
0057     throw cms::Exception("L1ScoutingEtSumOrbitFlatTableProducer")
0058         << "writePhysicalValues and writeHardwareValues cannot be false at the same time!";
0059   }
0060   produces<l1ScoutingRun3::OrbitFlatTable>();
0061 }
0062 
0063 void L1ScoutingEtSumOrbitFlatTableProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0064   edm::ParameterSetDescription desc;
0065 
0066   desc.add<edm::InputTag>("src");
0067   desc.add<std::string>("name");
0068   desc.add<std::string>("doc");
0069   desc.add<bool>("singleton", true)
0070       ->setComment("whether to output as singleton (one EtSum per bx) or not (multiple EtSums per bx)");
0071   desc.add<bool>("writePhysicalValues", true);
0072   desc.add<bool>("writeHardwareValues", false);
0073   desc.add<bool>("writeHF", true);
0074   desc.add<bool>("writeAsym", true);
0075   desc.add<bool>("writeMinBias", true);
0076   desc.add<bool>("writeTowerCount", true);
0077   desc.add<bool>("writeCentrality", true);
0078   desc.add<int>("ptPrecision", -1);
0079   desc.add<int>("phiPrecision", -1);
0080 
0081   descriptions.addDefault(desc);
0082 }
0083 
0084 void L1ScoutingEtSumOrbitFlatTableProducer::produce(edm::Event& iEvent, edm::EventSetup const&) {
0085   edm::Handle<l1ScoutingRun3::BxSumsOrbitCollection> src;
0086   iEvent.getByToken(src_, src);
0087 
0088   auto out = singleton_ ? produceSingle(*src) : produceMultiple(*src);
0089   iEvent.put(std::move(out));
0090 }
0091 
0092 std::unique_ptr<l1ScoutingRun3::OrbitFlatTable> L1ScoutingEtSumOrbitFlatTableProducer::produceSingle(
0093     l1ScoutingRun3::BxSumsOrbitCollection const& src) const {
0094   using namespace l1ScoutingRun3;
0095   auto out = std::make_unique<l1ScoutingRun3::OrbitFlatTable>(src.bxOffsets(), name_, /*singleton=*/true);
0096   out->setDoc(doc_);
0097 
0098   unsigned int nobjs = out->size();
0099 
0100   // physical values (float)
0101   std::vector<float> totalEt(nobjs);
0102   std::vector<float> totalEtEm(nobjs);
0103   std::vector<float> missEt(nobjs);
0104   std::vector<float> missEtPhi(nobjs);
0105   std::vector<float> missEtHF(nobjs);
0106   std::vector<float> missEtHFPhi(nobjs);
0107   std::vector<float> totalHt(nobjs);
0108   std::vector<float> missHt(nobjs);
0109   std::vector<float> missHtPhi(nobjs);
0110   std::vector<float> missHtHF(nobjs);
0111   std::vector<float> missHtHFPhi(nobjs);
0112   std::vector<float> asymEt(nobjs);
0113   std::vector<float> asymHt(nobjs);
0114   std::vector<float> asymEtHF(nobjs);
0115   std::vector<float> asymHtHF(nobjs);
0116 
0117   // hardware values (int)
0118   std::vector<int> hwTotalEt(nobjs);
0119   std::vector<int> hwTotalEtEm(nobjs);
0120   std::vector<int> hwMissEt(nobjs);
0121   std::vector<int> hwMissEtPhi(nobjs);
0122   std::vector<int> hwMissEtHF(nobjs);
0123   std::vector<int> hwMissEtHFPhi(nobjs);
0124   std::vector<int> hwTotalHt(nobjs);
0125   std::vector<int> hwMissHt(nobjs);
0126   std::vector<int> hwMissHtPhi(nobjs);
0127   std::vector<int> hwMissHtHF(nobjs);
0128   std::vector<int> hwMissHtHFPhi(nobjs);
0129   std::vector<int> hwAsymEt(nobjs);
0130   std::vector<int> hwAsymHt(nobjs);
0131   std::vector<int> hwAsymEtHF(nobjs);
0132   std::vector<int> hwAsymHtHF(nobjs);
0133 
0134   std::vector<int> minBiasHFP0(nobjs);
0135   std::vector<int> minBiasHFM0(nobjs);
0136   std::vector<int> minBiasHFP1(nobjs);
0137   std::vector<int> minBiasHFM1(nobjs);
0138   std::vector<int> towerCount(nobjs);
0139   std::vector<int> centrality(nobjs);
0140 
0141   for (unsigned int i = 0; i < nobjs; i++) {
0142     const auto& sums = src[i];
0143 
0144     // physical values
0145     totalEt[i] = demux::fEt(sums.hwTotalEt());
0146     totalEtEm[i] = demux::fEt(sums.hwTotalEtEm());
0147     missEt[i] = demux::fEt(sums.hwMissEt());
0148     missEtPhi[i] = demux::fPhi(sums.hwMissEtPhi());
0149     missEtHF[i] = demux::fEt(sums.hwMissEtHF());
0150     missEtHFPhi[i] = demux::fPhi(sums.hwMissEtHFPhi());
0151     totalHt[i] = demux::fEt(sums.hwTotalHt());
0152     missHt[i] = demux::fEt(sums.hwMissHt());
0153     missHtPhi[i] = demux::fPhi(sums.hwMissHtPhi());
0154     missHtHF[i] = demux::fEt(sums.hwMissHtHF());
0155     missHtHFPhi[i] = demux::fPhi(sums.hwMissHtHFPhi());
0156     asymEt[i] = demux::fEt(sums.hwAsymEt());
0157     asymHt[i] = demux::fEt(sums.hwAsymHt());
0158     asymEtHF[i] = demux::fEt(sums.hwAsymEtHF());
0159     asymHtHF[i] = demux::fEt(sums.hwAsymHtHF());
0160 
0161     // hardware values
0162     hwTotalEt[i] = sums.hwTotalEt();
0163     hwTotalEtEm[i] = sums.hwTotalEtEm();
0164     hwMissEt[i] = sums.hwMissEt();
0165     hwMissEtPhi[i] = sums.hwMissEtPhi();
0166     hwMissEtHF[i] = sums.hwMissEtHF();
0167     hwMissEtHFPhi[i] = sums.hwMissEtHFPhi();
0168     hwTotalHt[i] = sums.hwTotalHt();
0169     hwMissHt[i] = sums.hwMissHt();
0170     hwMissHtPhi[i] = sums.hwMissHtPhi();
0171     hwMissHtHF[i] = sums.hwMissHtHF();
0172     hwMissHtHFPhi[i] = sums.hwMissHtHFPhi();
0173     hwAsymEt[i] = sums.hwAsymEt();
0174     hwAsymHt[i] = sums.hwAsymHt();
0175     hwAsymEtHF[i] = sums.hwAsymEtHF();
0176     hwAsymHtHF[i] = sums.hwAsymHtHF();
0177 
0178     minBiasHFP0[i] = sums.minBiasHFP0();
0179     minBiasHFM0[i] = sums.minBiasHFM0();
0180     minBiasHFP1[i] = sums.minBiasHFP1();
0181     minBiasHFM1[i] = sums.minBiasHFM1();
0182     towerCount[i] = sums.towerCount();
0183     centrality[i] = sums.centrality();
0184   }
0185 
0186   // fill table
0187 
0188   if (writePhysicalValues_) {
0189     out->template addColumn<float>("totalEt", totalEt, "totalEt", ptPrecision_);
0190     out->template addColumn<float>("totalEtEm", totalEtEm, "totalEtEm", ptPrecision_);
0191     out->template addColumn<float>("missEt", missEt, "missEt pt", ptPrecision_);
0192     out->template addColumn<float>("missEtPhi", missEtPhi, "missEt phi", phiPrecision_);
0193     out->template addColumn<float>("totalHt", totalHt, "totalHt", ptPrecision_);
0194     out->template addColumn<float>("missHt", missHt, "missHt pt", ptPrecision_);
0195     out->template addColumn<float>("missHtPhi", missHtPhi, "missHt phi", phiPrecision_);
0196   }
0197   if (writeHardwareValues_) {
0198     out->template addColumn<int>("hwTotalEt", hwTotalEt, "hardware totalEt");
0199     out->template addColumn<int>("hwTotalEtEm", hwTotalEtEm, "hardware totalEtEm");
0200     out->template addColumn<int>("hwMissEt", hwMissEt, "hardware missEt pt");
0201     out->template addColumn<int>("hwMissEtPhi", hwMissEtPhi, "hardware missEt phi");
0202     out->template addColumn<int>("hwTotalHt", hwTotalHt, "hardware totalHt");
0203     out->template addColumn<int>("hwMissHt", hwMissHt, "hardware missHt pt");
0204     out->template addColumn<int>("hwMissHtPhi", hwMissHtPhi, "hardware missHt phi");
0205   }
0206 
0207   if (writeHF_) {
0208     if (writePhysicalValues_) {
0209       out->template addColumn<float>("missEtHF", missEtHF, "missEtHF", ptPrecision_);
0210       out->template addColumn<float>("missEtHFPhi", missEtHFPhi, "missEtHF phi", phiPrecision_);
0211       out->template addColumn<float>("missHtHF", missHtHF, "missHtHF pt", ptPrecision_);
0212       out->template addColumn<float>("missHtHFPhi", missHtHFPhi, "missHtHF phi", phiPrecision_);
0213     }
0214     if (writeHardwareValues_) {
0215       out->template addColumn<int>("hwMissEtHF", hwMissEtHF, "hardware missEtHF");
0216       out->template addColumn<int>("hwMissEtHFPhi", hwMissEtHFPhi, "hardware missEtHF phi");
0217       out->template addColumn<int>("hwMissHtHF", hwMissHtHF, "hardware missHtHF");
0218       out->template addColumn<int>("hwMissHtHFPhi", hwMissHtHFPhi, "hardware missHtHF phi");
0219     }
0220   }
0221   if (writeAsym_) {
0222     if (writePhysicalValues_) {
0223       out->template addColumn<float>("asymEt", asymEt, "asymEt", ptPrecision_);
0224       out->template addColumn<float>("asymHt", asymHt, "asymHt", ptPrecision_);
0225     }
0226     if (writeHardwareValues_) {
0227       out->template addColumn<int>("hwAsymEt", hwAsymEt, "hardware asymEt");
0228       out->template addColumn<int>("hwAsymHt", hwAsymHt, "hardware asymHt");
0229     }
0230   }
0231 
0232   if (writeAsym_ && writeHF_) {
0233     if (writePhysicalValues_) {
0234       out->template addColumn<float>("asymEtHF", asymEtHF, "asymEtHF", ptPrecision_);
0235       out->template addColumn<float>("asymHtHF", asymHtHF, "asymHtHF", ptPrecision_);
0236     }
0237     if (writeHardwareValues_) {
0238       out->template addColumn<int>("hwAsymEtHF", hwAsymEtHF, "asymEtHF");
0239       out->template addColumn<int>("hwAsymHtHF", hwAsymHtHF, "asymHtHF");
0240     }
0241   }
0242 
0243   if (writeMinBias_) {
0244     out->template addColumn<int>("minBiasHFP0", minBiasHFP0, "minBiasHFP0");
0245     out->template addColumn<int>("minBiasHFM0", minBiasHFM0, "minBiasHFM0");
0246     out->template addColumn<int>("minBiasHFP1", minBiasHFP1, "minBiasHFP1");
0247     out->template addColumn<int>("minBiasHFM1", minBiasHFM1, "minBiasHFM1");
0248   }
0249 
0250   if (writeTowerCount_) {
0251     out->template addColumn<int>("towerCount", towerCount, "towerCount");
0252   }
0253 
0254   if (writeCentrality_) {
0255     out->template addColumn<int>("centrality", centrality, "centrality");
0256   }
0257 
0258   return out;
0259 }
0260 
0261 std::unique_ptr<l1ScoutingRun3::OrbitFlatTable> L1ScoutingEtSumOrbitFlatTableProducer::produceMultiple(
0262     l1ScoutingRun3::BxSumsOrbitCollection const& src) const {
0263   using namespace l1ScoutingRun3;
0264   // compute number of objects per bx to adjust bxOffsets
0265   unsigned int nitems = 5;  // totalEt, totalEtEm, missEt, totalHt, missHt
0266   if (writeHF_)
0267     nitems += 2;  // missEtHF, missHtHF
0268   if (writeAsym_)
0269     nitems += (writeHF_ ? 4 : 2);  // asymEt, asymHt, asymEtHF, asymHtHF
0270   if (writeMinBias_)
0271     nitems += 4;  // minBiasHFP0, minBiasHFM0, minBiasHFP1, minBiasHFM1
0272   if (writeTowerCount_)
0273     nitems += 1;  // towerCount
0274   if (writeCentrality_)
0275     nitems += 1;  // centrality
0276 
0277   // adjust bxOffsets since each bx now contains multiple objects instead of single object
0278   std::vector<unsigned> offsets(src.bxOffsets());
0279   for (auto& v : offsets)
0280     v *= nitems;
0281 
0282   auto out = std::make_unique<l1ScoutingRun3::OrbitFlatTable>(offsets, name_, /*singleton=*/false);
0283   out->setDoc(doc_);
0284 
0285   unsigned int nobjs = out->size();
0286 
0287   // physical values
0288   std::vector<float> pt(nobjs);
0289   std::vector<float> phi(nobjs, 0.);
0290 
0291   // hardware values
0292   std::vector<int> hwEt(nobjs);
0293   std::vector<int> hwPhi(nobjs, 0);
0294 
0295   std::vector<int> sumType(nobjs);
0296 
0297   unsigned int i = 0;
0298   for (const l1ScoutingRun3::BxSums& sums : src) {
0299     assert(i + nitems <= nobjs && i % nitems == 0);
0300 
0301     // totalEt
0302     pt[i] = demux::fEt(sums.hwTotalEt());
0303     hwEt[i] = sums.hwTotalEt();
0304     sumType[i++] = l1t::EtSum::kTotalEt;
0305     // totalEtEm
0306     pt[i] = demux::fEt(sums.hwTotalEtEm());
0307     hwEt[i] = sums.hwTotalEtEm();
0308     sumType[i++] = l1t::EtSum::kTotalEtEm;
0309     // missEt
0310     pt[i] = demux::fEt(sums.hwMissEt());
0311     phi[i] = demux::fPhi(sums.hwMissEtPhi());
0312     hwEt[i] = sums.hwMissEt();
0313     hwPhi[i] = sums.hwMissEtPhi();
0314     sumType[i++] = l1t::EtSum::kMissingEt;
0315     // totalHt
0316     pt[i] = demux::fEt(sums.hwTotalHt());
0317     hwEt[i] = sums.hwTotalHt();
0318     sumType[i++] = l1t::EtSum::kTotalHt;
0319     // missHt
0320     pt[i] = demux::fEt(sums.hwMissHt());
0321     phi[i] = demux::fPhi(sums.hwMissHtPhi());
0322     hwEt[i] = sums.hwMissHt();
0323     hwPhi[i] = sums.hwMissHtPhi();
0324     sumType[i++] = l1t::EtSum::kMissingHt;
0325 
0326     if (writeHF_) {
0327       // missEtHF
0328       pt[i] = demux::fEt(sums.hwMissEtHF());
0329       phi[i] = demux::fPhi(sums.hwMissEtHFPhi());
0330       hwEt[i] = sums.hwMissEtHF();
0331       hwPhi[i] = sums.hwMissEtHFPhi();
0332       sumType[i++] = l1t::EtSum::kMissingEtHF;
0333       // missHtHF
0334       pt[i] = demux::fEt(sums.hwMissHtHF());
0335       phi[i] = demux::fPhi(sums.hwMissHtHFPhi());
0336       hwEt[i] = sums.hwMissHtHF();
0337       hwPhi[i] = sums.hwMissHtHFPhi();
0338       sumType[i++] = l1t::EtSum::kMissingHtHF;
0339     }
0340 
0341     if (writeAsym_) {
0342       // asymEt
0343       pt[i] = demux::fEt(sums.hwAsymEt());
0344       hwEt[i] = sums.hwAsymEt();
0345       sumType[i++] = l1t::EtSum::kAsymEt;
0346       // asymHt
0347       pt[i] = demux::fEt(sums.hwAsymHt());
0348       hwEt[i] = sums.hwAsymHt();
0349       sumType[i++] = l1t::EtSum::kAsymHt;
0350 
0351       if (writeHF_) {
0352         // asymEtHF
0353         pt[i] = demux::fEt(sums.hwAsymEtHF());
0354         hwEt[i] = sums.hwAsymEtHF();
0355         sumType[i++] = l1t::EtSum::kAsymEtHF;
0356         // asymHtHF
0357         pt[i] = demux::fEt(sums.hwAsymHtHF());
0358         hwEt[i] = sums.hwAsymHtHF();
0359         sumType[i++] = l1t::EtSum::kAsymHtHF;
0360       }
0361     }
0362 
0363     if (writeMinBias_) {
0364       // minBiasHFP0
0365       pt[i] = sums.minBiasHFP0();
0366       hwEt[i] = sums.minBiasHFP0();
0367       sumType[i++] = l1t::EtSum::kMinBiasHFP0;
0368       // minBiasHFM0
0369       pt[i] = sums.minBiasHFM0();
0370       hwEt[i] = sums.minBiasHFM0();
0371       sumType[i++] = l1t::EtSum::kMinBiasHFM0;
0372       // minBiasHFP1
0373       pt[i] = sums.minBiasHFP1();
0374       hwEt[i] = sums.minBiasHFP1();
0375       sumType[i++] = l1t::EtSum::kMinBiasHFP1;
0376       // minBiasHFM1
0377       pt[i] = sums.minBiasHFM1();
0378       hwEt[i] = sums.minBiasHFM1();
0379       sumType[i++] = l1t::EtSum::kMinBiasHFM1;
0380     }
0381 
0382     if (writeTowerCount_) {
0383       // towerCount
0384       pt[i] = sums.towerCount();
0385       hwEt[i] = sums.towerCount();
0386       sumType[i++] = l1t::EtSum::kTowerCount;
0387     }
0388 
0389     if (writeCentrality_) {
0390       // centrality
0391       pt[i] = sums.centrality();
0392       hwEt[i] = sums.centrality();
0393       sumType[i++] = l1t::EtSum::kCentrality;
0394     }
0395   }
0396 
0397   // fill table
0398 
0399   if (writePhysicalValues_) {
0400     out->template addColumn<float>("pt", pt, "pt", ptPrecision_);
0401     out->template addColumn<float>("phi", phi, "phi", phiPrecision_);
0402   }
0403   if (writeHardwareValues_) {
0404     out->template addColumn<int>("hwEt", pt, "hardware Et");
0405     out->template addColumn<int>("hwPhi", phi, "hardware phi");
0406   }
0407 
0408   out->template addColumn<int>(
0409       "etSumType",
0410       sumType,
0411       "the type of the EtSum "
0412       "(https://github.com/cms-sw/cmssw/blob/master/DataFormats/L1Trigger/interface/EtSum.h#L27-L56)");
0413   return out;
0414 }
0415 
0416 DEFINE_FWK_MODULE(L1ScoutingEtSumOrbitFlatTableProducer);