CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_13_patch3/src/PhysicsTools/HepMCCandAlgos/plugins/GenParticleDecaySelector.cc

Go to the documentation of this file.
00001 /* \class GenParticleDecaySelector
00002  *
00003  * \author Luca Lista, INFN
00004  *
00005  * \version $Id: GenParticleDecaySelector.cc,v 1.4 2010/02/11 00:12:56 wmtan Exp $
00006  *
00007  */
00008 
00009 #include "FWCore/Framework/interface/EDProducer.h"
00010 #include "FWCore/Utilities/interface/InputTag.h"
00011 #include "SimGeneral/HepPDTRecord/interface/PdtEntry.h"
00012 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00013 
00014 class GenParticleDecaySelector : public edm::EDProducer {
00015 public:
00017   GenParticleDecaySelector(const edm::ParameterSet&);
00018 private:
00020   void produce(edm::Event& e, const edm::EventSetup&);
00021   bool firstEvent_;
00023   edm::InputTag src_;  
00025   PdtEntry particle_;
00027   int status_;
00029   std::pair<reco::GenParticleRef, reco::GenParticle*>
00030   add(reco::GenParticleCollection&, const reco::GenParticle &, 
00031       reco::GenParticleRefProd);
00032 };
00033 
00034 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00035 #include "DataFormats/Common/interface/Handle.h"
00036 #include "FWCore/Framework/interface/Event.h"
00037 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00038 
00039 using namespace edm;
00040 using namespace reco;
00041 using namespace std;
00042 
00043 GenParticleDecaySelector::GenParticleDecaySelector(const edm::ParameterSet& cfg) :
00044   firstEvent_(true),
00045   src_(cfg.getParameter<InputTag>("src")),
00046   particle_(cfg.getParameter<PdtEntry>("particle")),
00047   status_(cfg.getParameter<int>("status")) {
00048   produces<GenParticleCollection>();
00049 }
00050 
00051 void GenParticleDecaySelector::produce(edm::Event& evt, const edm::EventSetup& es) {
00052   if (firstEvent_) {particle_.setup(es); firstEvent_ = false;}
00053 
00054   Handle<GenParticleCollection> genParticles;
00055   evt.getByLabel(src_, genParticles);
00056   auto_ptr<GenParticleCollection> decay(new GenParticleCollection);
00057   const GenParticleRefProd ref = evt.getRefBeforePut<GenParticleCollection>();
00058   for(GenParticleCollection::const_iterator g = genParticles->begin();
00059       g != genParticles->end(); ++g)
00060     if(g->pdgId() == particle_.pdgId() && g->status() == status_)
00061       add(*decay, *g, ref);
00062   evt.put(decay);
00063 }
00064 
00065 pair<GenParticleRef, GenParticle*> GenParticleDecaySelector::add(GenParticleCollection & decay, const GenParticle & p,
00066                                                                  GenParticleRefProd ref) {
00067   size_t idx = decay.size();
00068   GenParticleRef r(ref, idx);
00069   decay.resize(idx+1);
00070   const LeafCandidate & part = p;
00071   GenParticle g(part);
00072   size_t n = p.numberOfDaughters();
00073   for(size_t i = 0; i < n; ++i) {
00074     pair<GenParticleRef, GenParticle*> d = add(decay, *p.daughterRef(i), ref);
00075     d.second->addMother(r);
00076     g.addDaughter(d.first);
00077   }
00078   GenParticle & gp = decay[idx] = g;
00079   return make_pair(r, &gp);
00080 }
00081 
00082 
00083 #include "FWCore/Framework/interface/MakerMacros.h"
00084 
00085 DEFINE_FWK_MODULE(GenParticleDecaySelector);