CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GenParticleDecaySelector.cc
Go to the documentation of this file.
1 /* \class GenParticleDecaySelector
2  *
3  * \author Luca Lista, INFN
4  *
5  *
6  */
7 
15 
17 public:
20 
21  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
22 
23 private:
25  void produce(edm::Event& e, const edm::EventSetup&) override;
33  int status_;
35  std::pair<reco::GenParticleRef, reco::GenParticle*> add(reco::GenParticleCollection&,
36  const reco::GenParticle&,
38 };
39 
46 
47 using namespace edm;
48 using namespace reco;
49 using namespace std;
50 
52  : firstEvent_(true),
53  srcToken_(consumes<GenParticleCollection>(cfg.getParameter<InputTag>("src"))),
54  tableToken_(esConsumes()),
55  particle_(cfg.getParameter<PdtEntry>("particle")),
56  status_(cfg.getParameter<int>("status")) {
57  produces<GenParticleCollection>();
58 }
59 
62  desc.add<edm::InputTag>("src");
63  desc.addNode(edm::ParameterDescription<int>("particle", true) xor
65  ->setComment("the PdtEntry can be specified as either an 'int' or via its name using a 'string'");
66  desc.add<int>("status");
67 
68  descriptions.addDefault(desc);
69 }
70 
72  if (firstEvent_) {
73  auto const& pdt = es.getData(tableToken_);
74  particle_.setup(pdt);
75  firstEvent_ = false;
76  }
77 
79  evt.getByToken(srcToken_, genParticles);
80  auto decay = std::make_unique<GenParticleCollection>();
82  for (GenParticleCollection::const_iterator g = genParticles->begin(); g != genParticles->end(); ++g)
83  if (g->pdgId() == particle_.pdgId() && g->status() == status_)
84  add(*decay, *g, ref);
85  evt.put(std::move(decay));
86 }
87 
88 pair<GenParticleRef, GenParticle*> GenParticleDecaySelector::add(GenParticleCollection& decay,
89  const GenParticle& p,
90  GenParticleRefProd ref) {
91  size_t idx = decay.size();
92  GenParticleRef r(ref, idx);
93  decay.resize(idx + 1);
94  const LeafCandidate& part = p;
95  GenParticle g(part);
96  size_t n = p.numberOfDaughters();
97  for (size_t i = 0; i < n; ++i) {
98  pair<GenParticleRef, GenParticle*> d = add(decay, *p.daughterRef(i), ref);
99  d.second->addMother(r);
100  g.addDaughter(d.first);
101  }
102  GenParticle& gp = decay[idx] = g;
103  return make_pair(r, &gp);
104 }
105 
107 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
std::pair< reco::GenParticleRef, reco::GenParticle * > add(reco::GenParticleCollection &, const reco::GenParticle &, reco::GenParticleRefProd)
recursively add a new particle to the output collection
void setComment(std::string const &value)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
tuple cfg
Definition: looper.py:296
daughters::value_type daughterRef(size_type i) const
reference to daughter at given position
void setup(const HepPDT::ParticleDataTable &)
fill data from Event Setup
Definition: PdtEntry.cc:26
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
size_t numberOfDaughters() const override
number of daughters
ParameterDescriptionNode * addNode(ParameterDescriptionNode const &node)
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
void addDaughter(const typename daughters::value_type &)
add a daughter via a reference
bool getData(T &iHolder) const
Definition: EventSetup.h:122
tuple d
Definition: ztail.py:151
GenParticleDecaySelector(const edm::ParameterSet &)
constructor
void addDefault(ParameterSetDescription const &psetDescription)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
def move
Definition: eostools.py:511
ParameterDescriptionBase * add(U const &iLabel, T const &value)
RefProd< PROD > getRefBeforePut()
Definition: Event.h:158
edm::EDGetTokenT< reco::GenParticleCollection > srcToken_
source collection name
part
Definition: HCALResponse.h:20
edm::ESGetToken< HepPDT::ParticleDataTable, edm::DefaultRecord > tableToken_
particle type
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
int pdgId() const
PDG id.
Definition: PdtEntry.cc:5
void produce(edm::Event &e, const edm::EventSetup &) override
process one event