CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
ParticleDecayProducer.cc
Go to the documentation of this file.
1 
7 // user include files
13 #include <vector>
14 
16 public:
18  ~ParticleDecayProducer() override;
19 
20 private:
21  void produce(edm::Event&, const edm::EventSetup&) override;
24  std::vector<int> daughtersPdgId_;
26  std::vector<std::string> valias;
27 };
28 
29 // Candidate handling
34 #include <sstream>
35 using namespace edm;
36 using namespace std;
37 using namespace reco;
38 
39 // constructors
41  : genCandidatesToken_(consumes<CandidateCollection>(iConfig.getParameter<InputTag>("src"))),
42  motherPdgId_(iConfig.getParameter<int>("motherPdgId")),
43  daughtersPdgId_(iConfig.getParameter<vector<int> >("daughtersPdgId")),
44  decayChain_(iConfig.getParameter<std::string>("decayChain")) {
45  string alias;
46  produces<CandidateCollection>(alias = decayChain_ + "Mother").setBranchAlias(alias);
47  for (unsigned int j = 0; j < daughtersPdgId_.size(); ++j) {
48  ostringstream index, collection;
49  index << j;
50  collection << decayChain_ << "Lepton" << index.str();
51  valias.push_back(collection.str());
52  produces<CandidateCollection>(valias.at(j)).setBranchAlias(valias.at(j));
53  }
54 }
55 
56 // destructor
58 
60  // get gen particle candidates
61  edm::Handle<CandidateCollection> genCandidatesCollection;
62  iEvent.getByToken(genCandidatesToken_, genCandidatesCollection);
63 
64  unique_ptr<CandidateCollection> mothercands(new CandidateCollection);
65  unique_ptr<CandidateCollection> daughterscands(new CandidateCollection);
66  size_t daughtersize = daughtersPdgId_.size();
67  for (CandidateCollection::const_iterator p = genCandidatesCollection->begin(); p != genCandidatesCollection->end();
68  ++p) {
69  if (p->pdgId() == motherPdgId_ && p->status() == 3) {
70  mothercands->push_back(p->clone());
71  size_t ndau = p->numberOfDaughters();
72  for (size_t i = 0; i < ndau; ++i) {
73  for (size_t j = 0; j < daughtersize; ++j) {
74  if (p->daughter(i)->pdgId() == daughtersPdgId_[j] && p->daughter(i)->status() == 3) {
75  daughterscands->push_back(p->daughter(i)->clone());
76  }
77  }
78  }
79  }
80  }
81 
82  iEvent.put(std::move(mothercands), decayChain_ + "Mother");
83  daughterscands->sort(GreaterByPt<reco::Candidate>());
84 
85  for (unsigned int row = 0; row < daughtersize; ++row) {
86  unique_ptr<CandidateCollection> leptonscands_(new CandidateCollection);
87  leptonscands_->push_back((daughterscands->begin() + row)->clone());
88  iEvent.put(std::move(leptonscands_), valias.at(row));
89  }
90 }
91 
93 
void produce(edm::Event &, const edm::EventSetup &) override
ParticleDecayProducer(const edm::ParameterSet &)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::CandidateCollection > genCandidatesToken_
std::vector< std::string > valias
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
std::vector< int > daughtersPdgId_
TEveGeoShape * clone(const TEveElement *element, TEveElement *parent)
Definition: eve_macros.cc:135