Skip to content

Commit

Permalink
Don't fail if there are no LHE weights
Browse files Browse the repository at this point in the history
  • Loading branch information
kdlong committed Aug 1, 2020
1 parent 35b32c4 commit 3180b59
Show file tree
Hide file tree
Showing 3 changed files with 32 additions and 20 deletions.
45 changes: 27 additions & 18 deletions PhysicsTools/NanoAOD/plugins/LHEWeightsTableProducer.cc
Original file line number Diff line number Diff line change
Expand Up @@ -48,25 +48,32 @@ class LHEWeightsTableProducer : public edm::global::EDProducer<edm::LuminosityBl

void produce(edm::StreamID id, edm::Event& iEvent, const edm::EventSetup& iSetup) const override {
edm::Handle<GenWeightProduct> lheWeightHandle;
for (auto& token : lheWeightTokens_) {
iEvent.getByToken(token, lheWeightHandle);
if (lheWeightHandle.isValid()) {
break;
bool foundLheWeights = false;
if (foundLheWeights) {
for (auto& token : lheWeightTokens_) {
iEvent.getByToken(token, lheWeightHandle);
if (lheWeightHandle.isValid()) {
break;
}
}
}

const GenWeightProduct* lheWeightProduct = lheWeightHandle.product();
WeightsContainer lheWeights = lheWeightProduct->weights();
WeightsContainer lheWeights;
if (foundLheWeights) {
const GenWeightProduct* lheWeightProduct = lheWeightHandle.product();
lheWeights = lheWeightProduct->weights();
}

edm::Handle<GenWeightProduct> genWeightHandle;
iEvent.getByToken(genWeightToken_, genWeightHandle);
const GenWeightProduct* genWeightProduct = genWeightHandle.product();
WeightsContainer genWeights = genWeightProduct->weights();

auto lheWeightTables = std::make_unique<std::vector<nanoaod::FlatTable>>();
auto const& weightInfos = *luminosityBlockCache(iEvent.getLuminosityBlock().index());

addWeightGroupToTable(lheWeightTables, "LHE", weightInfos.at(inLHE), lheWeights);
auto lheWeightTables = std::make_unique<std::vector<nanoaod::FlatTable>>();
if (foundLheWeights)
addWeightGroupToTable(lheWeightTables, "LHE", weightInfos.at(inLHE), lheWeights);
addWeightGroupToTable(lheWeightTables, "Gen", weightInfos.at(inGen), genWeights);

iEvent.put(std::move(lheWeightTables));
Expand All @@ -76,14 +83,13 @@ class LHEWeightsTableProducer : public edm::global::EDProducer<edm::LuminosityBl
const char* typeName,
const WeightGroupDataContainer& weightInfos,
WeightsContainer& allWeights) const {
size_t typeCount = 0;
gen::WeightType previousType = gen::WeightType::kUnknownWeights;
std::unordered_map<gen::WeightType, int> typeCount = {};
for (auto& type : gen::allWeightTypes)
typeCount[type] = 0;

for (const auto& groupInfo : weightInfos) {
std::string entryName = typeName;
gen::WeightType weightType = groupInfo.group->weightType();
if (previousType != weightType)
typeCount = 0;

std::string name = weightTypeNames_.at(weightType);
std::string label = groupInfo.group->name();
Expand All @@ -102,28 +108,29 @@ class LHEWeightsTableProducer : public edm::global::EDProducer<edm::LuminosityBl

entryName.append(name);
entryName.append("Weight");
if (typeCount > 0) {
if (typeCount[weightType] > 0) {
entryName.append("AltSet");
entryName.append(std::to_string(typeCount));
entryName.append(std::to_string(typeCount[weightType]));
}

lheWeightTables->emplace_back(weights.size(), entryName, false);
lheWeightTables->back().addColumn<float>(
"", weights, label, nanoaod::FlatTable::FloatColumn, lheWeightPrecision_);

typeCount++;
previousType = weightType;
typeCount[weightType]++;
}
}

std::shared_ptr<WeightGroupsToStore> globalBeginLuminosityBlock(edm::LuminosityBlock const& iLumi,
edm::EventSetup const&) const override {
// Set equal to the max number of groups
// subtrack 1 for each weight group you find
bool foundLheWeights = false;
edm::Handle<GenWeightInfoProduct> lheWeightInfoHandle;
for (auto& token : lheWeightInfoTokens_) {
iLumi.getByToken(token, lheWeightInfoHandle);
if (lheWeightInfoHandle.isValid()) {
foundLheWeights = true;
break;
}
}
Expand All @@ -137,8 +144,10 @@ class LHEWeightsTableProducer : public edm::global::EDProducer<edm::LuminosityBl

WeightGroupsToStore weightsToStore;
for (auto weightType : gen::allWeightTypes) {
auto lheWeights = weightDataPerType(lheWeightInfoHandle, weightType, storePerType[weightType]);
weightsToStore.at(inLHE).insert(weightsToStore.at(inLHE).end(), lheWeights.begin(), lheWeights.end());
if (foundLheWeights) {
auto lheWeights = weightDataPerType(lheWeightInfoHandle, weightType, storePerType[weightType]);
weightsToStore.at(inLHE).insert(weightsToStore.at(inLHE).end(), lheWeights.begin(), lheWeights.end());
}

auto genWeights = weightDataPerType(genWeightInfoHandle, weightType, storePerType[weightType]);
weightsToStore.at(inGen).insert(weightsToStore.at(inGen).end(), genWeights.begin(), genWeights.end());
Expand Down
2 changes: 1 addition & 1 deletion PhysicsTools/NanoAOD/python/nanogen_cff.py
Original file line number Diff line number Diff line change
Expand Up @@ -78,7 +78,7 @@

nanogenMiniSequence = cms.Sequence(
genWeights+
lheWeights+
#lheWeights+
nanoMetadata+
mergedGenParticles+
genParticles2HepMC+
Expand Down
5 changes: 4 additions & 1 deletion PhysicsTools/NanoAOD/test/testNanoWeights.py
Original file line number Diff line number Diff line change
Expand Up @@ -19,7 +19,10 @@ def variableAndNumber(varName, tree):
rtfile = ROOT.TFile(args.inputFile)
tree = rtfile.Get("Events")
tree.GetEntry(0)
variables = ["LHEScaleWeight", "LHEPdfWeight", "LHEMEParamWeight", "GenPartonShowerWeight", "LHEUnknownWeight", ]
types = ["ScaleWeight", "PdfWeight", "MEParamWeight", "UnknownWeight", ]
variables = ["LHE"+t for t in types]
variables.append("GenPartonShowerWeight")
variables.extend(["Gen"+t for t in types])

for varName in variables:
variableAndNumber(varName, tree)
Expand Down

0 comments on commit 3180b59

Please sign in to comment.