CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATPackedGenParticleProducer.cc
Go to the documentation of this file.
1 #include <string>
2 
3 
22 
26 
33 #include "CLHEP/Vector/ThreeVector.h"
34 #include "CLHEP/Vector/LorentzVector.h"
35 #include "CLHEP/Matrix/Vector.h"
36 #include <string>
37 
38 
39 namespace pat {
41  public:
44 
45  virtual void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const;
46 
47  private:
52  const double maxRapidity_;
53  };
54 }
55 
57  Cands_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("inputCollection"))),
58  GenOrigs_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("inputOriginal"))),
59  Asso_(consumes<edm::Association<reco::GenParticleCollection> >(iConfig.getParameter<edm::InputTag>("map"))),
60  AssoOriginal_(consumes<edm::Association<reco::GenParticleCollection> >(iConfig.getParameter<edm::InputTag>("inputCollection"))),
61  maxRapidity_(iConfig.getParameter<double>("maxRapidity"))
62 {
63  produces< std::vector<pat::PackedGenParticle> > ();
64  produces< edm::Association< std::vector<pat::PackedGenParticle> > >();
65 }
66 
68 
70 
71 
73  iEvent.getByToken( Cands_, cands );
74  std::vector<reco::Candidate>::const_iterator cand;
75 
77  iEvent.getByToken( Asso_, asso );
78 
80  iEvent.getByToken( AssoOriginal_, assoOriginal);
81 
83  iEvent.getByToken( GenOrigs_, genOrigs);
84  std::vector<int> mapping(genOrigs->size(), -1);
85 
86  //invert the value map from Orig2New to New2Orig
87  std::map< edm::Ref<reco::GenParticleCollection> , edm::Ref<reco::GenParticleCollection> > reverseMap;
88  for(unsigned int ic=0, nc = genOrigs->size(); ic < nc; ++ic)
89  {
91  edm::Ref<reco::GenParticleCollection> newRef = (*assoOriginal)[originalRef];
92  reverseMap.insert(std::pair<edm::Ref<reco::GenParticleCollection>,edm::Ref<reco::GenParticleCollection>>(newRef,originalRef));
93  }
94 
95  std::auto_ptr< std::vector<pat::PackedGenParticle> > outPtrP( new std::vector<pat::PackedGenParticle> );
96 
97  unsigned int packed=0;
98  for(unsigned int ic=0, nc = cands->size(); ic < nc; ++ic) {
99  const reco::GenParticle &cand=(*cands)[ic];
100  if(cand.status() ==1 && std::abs(cand.y()) < maxRapidity_)
101  {
102  // Obtain original gen particle collection reference from input reference and map
104  edm::Ref<reco::GenParticleCollection> originalRef=reverseMap[inputRef];
105  mapping[originalRef.key()]=packed;
106  packed++;
107  if(cand.numberOfMothers() > 0) {
108  edm::Ref<reco::GenParticleCollection> newRef=(*asso)[cand.motherRef(0)];
109  outPtrP->push_back( pat::PackedGenParticle(cand,newRef));
110  } else {
112 
113  }
114 
115  }
116  }
117 
118 
120 
121  std::auto_ptr<edm::Association< std::vector<pat::PackedGenParticle> > > gp2pgp(new edm::Association< std::vector<pat::PackedGenParticle> > (oh ));
122  edm::Association< std::vector<pat::PackedGenParticle> >::Filler gp2pgpFiller(*gp2pgp);
123  gp2pgpFiller.insert(genOrigs, mapping.begin(), mapping.end());
124  gp2pgpFiller.fill();
125  iEvent.put(gp2pgp);
126 
127 
128 }
129 
130 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
PATPackedGenParticleProducer(const edm::ParameterSet &)
virtual double y() const final
rapidity
key_type key() const
Accessor for product key.
Definition: Ref.h:264
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const
virtual int status() const final
status word
const edm::EDGetTokenT< reco::GenParticleCollection > Cands_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
virtual size_t numberOfMothers() const
number of mothers
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > Asso_
daughters::value_type motherRef(size_type i=0) const
reference to mother at given position
const edm::EDGetTokenT< reco::GenParticleCollection > GenOrigs_
const edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > AssoOriginal_