File indexing completed on 2024-04-06 12:01:02
0001
0002
0003
0004
0005
0006
0007
0008 #include "FWCore/Framework/interface/Frameworkfwd.h"
0009 #include "FWCore/Framework/interface/global/EDProducer.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Utilities/interface/InputTag.h"
0012 #include "DataFormats/Candidate/interface/Candidate.h"
0013 #include <vector>
0014
0015 class ParticleDecayProducer : public edm::global::EDProducer<> {
0016 public:
0017 explicit ParticleDecayProducer(const edm::ParameterSet&);
0018 ~ParticleDecayProducer() override = default;
0019
0020 private:
0021 void produce(edm::StreamID, edm::Event& e, edm::EventSetup const& c) const override;
0022 edm::EDGetTokenT<reco::CandidateCollection> genCandidatesToken_;
0023 const int motherPdgId_;
0024 const std::vector<int> daughtersPdgId_;
0025 const std::string decayChain_;
0026 std::vector<std::string> valias;
0027 };
0028
0029
0030 #include "DataFormats/Candidate/interface/CandMatchMap.h"
0031 #include "DataFormats/Common/interface/AssociationVector.h"
0032 #include "CommonTools/Utils/interface/PtComparator.h"
0033 #include "FWCore/Framework/interface/Event.h"
0034 #include <sstream>
0035 using namespace edm;
0036 using namespace std;
0037 using namespace reco;
0038
0039
0040 ParticleDecayProducer::ParticleDecayProducer(const edm::ParameterSet& iConfig)
0041 : genCandidatesToken_(consumes<CandidateCollection>(iConfig.getParameter<InputTag>("src"))),
0042 motherPdgId_(iConfig.getParameter<int>("motherPdgId")),
0043 daughtersPdgId_(iConfig.getParameter<vector<int> >("daughtersPdgId")),
0044 decayChain_(iConfig.getParameter<std::string>("decayChain")) {
0045 string alias;
0046 produces<CandidateCollection>(alias = decayChain_ + "Mother").setBranchAlias(alias);
0047 for (unsigned int j = 0; j < daughtersPdgId_.size(); ++j) {
0048 ostringstream index, collection;
0049 index << j;
0050 collection << decayChain_ << "Lepton" << index.str();
0051 valias.push_back(collection.str());
0052 produces<CandidateCollection>(valias.at(j)).setBranchAlias(valias.at(j));
0053 }
0054 }
0055
0056 void ParticleDecayProducer::produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const&) const {
0057
0058 const edm::Handle<CandidateCollection> genCandidatesCollection = iEvent.getHandle(genCandidatesToken_);
0059
0060 unique_ptr<CandidateCollection> mothercands(new CandidateCollection);
0061 unique_ptr<CandidateCollection> daughterscands(new CandidateCollection);
0062 size_t daughtersize = daughtersPdgId_.size();
0063 for (CandidateCollection::const_iterator p = genCandidatesCollection->begin(); p != genCandidatesCollection->end();
0064 ++p) {
0065 if (p->pdgId() == motherPdgId_ && p->status() == 3) {
0066 mothercands->push_back(p->clone());
0067 size_t ndau = p->numberOfDaughters();
0068 for (size_t i = 0; i < ndau; ++i) {
0069 for (size_t j = 0; j < daughtersize; ++j) {
0070 if (p->daughter(i)->pdgId() == daughtersPdgId_[j] && p->daughter(i)->status() == 3) {
0071 daughterscands->push_back(p->daughter(i)->clone());
0072 }
0073 }
0074 }
0075 }
0076 }
0077
0078 iEvent.put(std::move(mothercands), decayChain_ + "Mother");
0079 daughterscands->sort(GreaterByPt<reco::Candidate>());
0080
0081 for (unsigned int row = 0; row < daughtersize; ++row) {
0082 unique_ptr<CandidateCollection> leptonscands_(new CandidateCollection);
0083 leptonscands_->push_back((daughterscands->begin() + row)->clone());
0084 iEvent.put(std::move(leptonscands_), valias.at(row));
0085 }
0086 }
0087
0088 #include "FWCore/Framework/interface/MakerMacros.h"
0089 DEFINE_FWK_MODULE(ParticleDecayProducer);