File indexing completed on 2024-04-06 12:21:15
0001 #include "L1Trigger/L1TNtuples/interface/L1AnalysisEvent.h"
0002
0003 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005
0006 #include <string>
0007 #include <iostream>
0008 #include <sys/stat.h>
0009
0010 L1Analysis::L1AnalysisEvent::L1AnalysisEvent(std::string puMCFile,
0011 std::string puMCHist,
0012 std::string puDataFile,
0013 std::string puDataHist,
0014 bool useAvgVtx,
0015 double maxWeight,
0016 edm::ConsumesCollector&& iConsumes)
0017 : fillHLT_(true), doPUWeights_(false), useAvgVtx_(useAvgVtx), maxAllowedWeight_(maxWeight), lumiWeights_() {
0018 pileupSummaryInfoToken_ = iConsumes.consumes<std::vector<PileupSummaryInfo>>(edm::InputTag("addPileupInfo"));
0019
0020 struct stat buf;
0021 if ((stat(puMCFile.c_str(), &buf) != -1) && (stat(puDataFile.c_str(), &buf) != -1)) {
0022 lumiWeights_ = edm::LumiReWeighting(puMCFile, puDataFile, puMCHist, puDataHist);
0023 doPUWeights_ = true;
0024 } else {
0025 edm::LogWarning("L1Prompt") << "No PU reweighting inputs - not going to calculate weights" << std::endl;
0026 }
0027 }
0028
0029 L1Analysis::L1AnalysisEvent::~L1AnalysisEvent() {}
0030
0031 void L1Analysis::L1AnalysisEvent::Set(const edm::Event& e, const edm::EDGetTokenT<edm::TriggerResults>& hlt_) {
0032 event_.run = e.id().run();
0033 event_.event = e.id().event();
0034 event_.time = e.time().value();
0035 event_.bx = e.bunchCrossing();
0036 event_.lumi = e.luminosityBlock();
0037 event_.orbit = e.orbitNumber();
0038
0039 if (!hlt_.isUninitialized()) {
0040 edm::Handle<edm::TriggerResults> hltresults;
0041 e.getByToken(hlt_, hltresults);
0042 const edm::TriggerNames& TrigNames_ = e.triggerNames(*hltresults);
0043 const int ntrigs = hltresults->size();
0044
0045 for (int itr = 0; itr < ntrigs; itr++) {
0046 TString trigName = TrigNames_.triggerName(itr);
0047 if (!hltresults->accept(itr))
0048 continue;
0049 event_.hlt.push_back(trigName);
0050 }
0051 }
0052
0053
0054 double weight = 1.;
0055
0056 if (doPUWeights_ && (!e.eventAuxiliary().isRealData())) {
0057 edm::Handle<std::vector<PileupSummaryInfo>> puInfo;
0058 e.getByLabel(edm::InputTag("addPileupInfo"), puInfo);
0059
0060 if (puInfo.isValid()) {
0061 std::vector<PileupSummaryInfo>::const_iterator pvi;
0062
0063 float npv = -1;
0064 for (pvi = puInfo->begin(); pvi != puInfo->end(); ++pvi) {
0065 int bx = pvi->getBunchCrossing();
0066
0067 if (bx == 0) {
0068 npv = useAvgVtx_ ? pvi->getTrueNumInteractions() : pvi->getPU_NumInteractions();
0069 continue;
0070 }
0071 }
0072
0073 weight = lumiWeights_.weight(npv);
0074 if (maxAllowedWeight_ > 0. && weight > maxAllowedWeight_)
0075 weight = maxAllowedWeight_;
0076 }
0077 }
0078
0079 if (!e.eventAuxiliary().isRealData()) {
0080 edm::Handle<std::vector<PileupSummaryInfo>> puInfo;
0081 e.getByToken(pileupSummaryInfoToken_, puInfo);
0082 if (puInfo.isValid()) {
0083 for (std::vector<PileupSummaryInfo>::const_iterator pvi = puInfo->begin(); pvi != puInfo->end(); pvi++) {
0084 int bx = pvi->getBunchCrossing();
0085 if (bx == 0) {
0086 event_.nPV = pvi->getPU_NumInteractions();
0087 event_.nPV_True = pvi->getTrueNumInteractions();
0088 }
0089 }
0090 }
0091 }
0092
0093 event_.puWeight = weight;
0094 }