File indexing completed on 2025-03-26 01:51:26
0001 #include "FWCore/Framework/interface/global/EDProducer.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0004 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0006 #include "DataFormats/NanoAOD/interface/FlatTable.h"
0007 #include "DataFormats/Common/interface/View.h"
0008 #include "DataFormats/Candidate/interface/Candidate.h"
0009 #include "DataFormats/PatCandidates/interface/PackedGenParticle.h"
0010 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0011 #include <DataFormats/Math/interface/deltaR.h>
0012 #include "DataFormats/JetReco/interface/GenJetCollection.h"
0013 #include "PhysicsTools/JetMCUtils/interface/CandMCTag.h"
0014
0015 #include <vector>
0016 #include <iostream>
0017
0018 class GenCandMotherTableProducer : public edm::global::EDProducer<> {
0019 public:
0020 GenCandMotherTableProducer(edm::ParameterSet const ¶ms)
0021 : objName_(params.getParameter<std::string>("objName")),
0022 branchName_(params.getParameter<std::string>("branchName")),
0023 src_(consumes<edm::View<pat::PackedGenParticle>>(params.getParameter<edm::InputTag>("src"))),
0024 candMap_(consumes<edm::Association<reco::GenParticleCollection>>(params.getParameter<edm::InputTag>("mcMap"))) {
0025 produces<nanoaod::FlatTable>();
0026 }
0027
0028 ~GenCandMotherTableProducer() override {}
0029
0030 void produce(edm::StreamID id, edm::Event &iEvent, const edm::EventSetup &iSetup) const override {
0031 edm::Handle<edm::View<pat::PackedGenParticle>> cands;
0032 iEvent.getByToken(src_, cands);
0033 unsigned int ncand = cands->size();
0034
0035 auto tab = std::make_unique<nanoaod::FlatTable>(ncand, objName_, false, true);
0036
0037 edm::Handle<edm::Association<reco::GenParticleCollection>> map;
0038 iEvent.getByToken(candMap_, map);
0039
0040 std::vector<int> key(ncand, -1), fromB(ncand, 0), fromC(ncand, 0);
0041 for (unsigned int i = 0; i < ncand; ++i) {
0042 reco::GenParticleRef motherRef = cands->at(i).motherRef();
0043 reco::GenParticleRef match = (*map)[motherRef];
0044
0045 if (match.isNull())
0046 continue;
0047
0048 key[i] = match.key();
0049 fromB[i] = isFromB(cands->at(i));
0050 fromC[i] = isFromC(cands->at(i));
0051 }
0052
0053 tab->addColumn<int>(branchName_ + "MotherIdx", key, "Mother index into GenPart list");
0054 tab->addColumn<uint8_t>("isFromB", fromB, "Is from B hadron: no: 0, any: 1, final: 2");
0055 tab->addColumn<uint8_t>("isFromC", fromC, "Is from C hadron: no: 0, any: 1, final: 2");
0056 iEvent.put(std::move(tab));
0057 }
0058
0059 bool isFinalB(const reco::Candidate &particle) const {
0060 if (!CandMCTagUtils::hasBottom(particle))
0061 return false;
0062
0063
0064 unsigned int npart = particle.numberOfDaughters();
0065
0066 for (size_t i = 0; i < npart; ++i) {
0067 if (CandMCTagUtils::hasBottom(*particle.daughter(i)))
0068 return false;
0069 }
0070
0071 return true;
0072 }
0073
0074 int isFromB(const reco::Candidate &particle) const {
0075 int fromB = 0;
0076
0077 unsigned int npart = particle.numberOfMothers();
0078 for (size_t i = 0; i < npart; ++i) {
0079 const reco::Candidate &mom = *particle.mother(i);
0080 if (CandMCTagUtils::hasBottom(mom)) {
0081 fromB = isFinalB(mom) ? 2 : 1;
0082 break;
0083 } else
0084 fromB = isFromB(mom);
0085 }
0086 return fromB;
0087 }
0088
0089 bool isFinalC(const reco::Candidate &particle) const {
0090 if (!CandMCTagUtils::hasCharm(particle))
0091 return false;
0092
0093
0094 unsigned int npart = particle.numberOfDaughters();
0095
0096 for (size_t i = 0; i < npart; ++i) {
0097 if (CandMCTagUtils::hasCharm(*particle.daughter(i)))
0098 return false;
0099 }
0100
0101 return true;
0102 }
0103
0104 int isFromC(const reco::Candidate &particle) const {
0105 int fromC = 0;
0106
0107 unsigned int npart = particle.numberOfMothers();
0108 for (size_t i = 0; i < npart; ++i) {
0109 const reco::Candidate &mom = *particle.mother(i);
0110 if (CandMCTagUtils::hasCharm(mom)) {
0111 fromC = isFinalC(mom) ? 2 : 1;
0112 break;
0113 } else
0114 fromC = isFromC(mom);
0115 }
0116 return fromC;
0117 }
0118
0119 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0120 edm::ParameterSetDescription desc;
0121 desc.add<std::string>("objName", "GenCands")
0122 ->setComment("name of the nanoaod::FlatTable to extend with this table");
0123 desc.add<std::string>("branchName", "GenPart")
0124 ->setComment(
0125 "name of the column to write (the final branch in the nanoaod will be <objName>_<branchName>Idx and "
0126 "<objName>_<branchName>Flav");
0127 desc.add<edm::InputTag>("src", edm::InputTag("packedGenParticles"))
0128 ->setComment("collection of the packedGenParticles, with association to prunedGenParticles");
0129 desc.add<edm::InputTag>("mcMap", edm::InputTag("finalGenParticles"))
0130 ->setComment(
0131 "tag to an edm::Association<GenParticleCollection> mapping prunedGenParticles to finalGenParticles");
0132 desc.addOptional<edm::InputTag>("genparticles", edm::InputTag("finalGenparticles"))
0133 ->setComment("Collection of genParticles to be mapped.");
0134 descriptions.add("genCandMotherTable", desc);
0135 }
0136
0137 protected:
0138 const std::string objName_, branchName_, doc_;
0139 const edm::EDGetTokenT<edm::View<pat::PackedGenParticle>> src_;
0140 const edm::EDGetTokenT<edm::Association<reco::GenParticleCollection>> candMap_;
0141 edm::EDGetTokenT<reco::GenParticleCollection> genPartsToken_;
0142 };
0143
0144 #include "FWCore/Framework/interface/MakerMacros.h"
0145 DEFINE_FWK_MODULE(GenCandMotherTableProducer);