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::Event&, const edm::EventSetup&);
46 
47  private:
53  double maxRapidity_;
54  };
55 }
56 
58  Cands_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("inputCollection"))),
59  GenOrigs_(consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("inputOriginal"))),
60  Asso_(consumes<edm::Association<reco::GenParticleCollection> >(iConfig.getParameter<edm::InputTag>("map"))),
61  AssoOriginal_(consumes<edm::Association<reco::GenParticleCollection> >(iConfig.getParameter<edm::InputTag>("inputCollection"))),
62  PVs_(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("inputVertices"))),
63  maxRapidity_(iConfig.getParameter<double>("maxRapidity"))
64 {
65  produces< std::vector<pat::PackedGenParticle> > ();
66  produces< edm::Association< std::vector<pat::PackedGenParticle> > >();
67 
68 }
69 
71 
73 
74 
76  iEvent.getByToken( Cands_, cands );
77  std::vector<reco::Candidate>::const_iterator cand;
78 
80  iEvent.getByToken( Asso_, asso );
81 
83  iEvent.getByToken( AssoOriginal_, assoOriginal);
84 
86  iEvent.getByToken( GenOrigs_, genOrigs);
87  std::vector<int> mapping(genOrigs->size(), -1);
88 
89 
91  iEvent.getByToken( PVs_, PVs );
92  reco::VertexRef PV(PVs.id());
93  math::XYZPoint PVpos;
94  if (!PVs->empty()) {
95  PV = reco::VertexRef(PVs, 0);
96  PVpos = PV->position();
97  }
98 
99  //invert the value map from Orig2New to New2Orig
100  std::map< edm::Ref<reco::GenParticleCollection> , edm::Ref<reco::GenParticleCollection> > reverseMap;
101  for(unsigned int ic=0, nc = genOrigs->size(); ic < nc; ++ic)
102  {
104  edm::Ref<reco::GenParticleCollection> newRef = (*assoOriginal)[originalRef];
105  reverseMap.insert(std::pair<edm::Ref<reco::GenParticleCollection>,edm::Ref<reco::GenParticleCollection>>(newRef,originalRef));
106  }
107 
108  std::auto_ptr< std::vector<pat::PackedGenParticle> > outPtrP( new std::vector<pat::PackedGenParticle> );
109 
110  unsigned int packed=0;
111  for(unsigned int ic=0, nc = cands->size(); ic < nc; ++ic) {
112  const reco::GenParticle &cand=(*cands)[ic];
113  if(cand.status() ==1 && std::abs(cand.y()) < maxRapidity_)
114  {
115  // Obtain original gen particle collection reference from input reference and map
117  edm::Ref<reco::GenParticleCollection> originalRef=reverseMap[inputRef];
118  mapping[originalRef.key()]=packed;
119  packed++;
120  if(cand.numberOfMothers() > 0) {
121  edm::Ref<reco::GenParticleCollection> newRef=(*asso)[cand.motherRef(0)];
122  outPtrP->push_back( pat::PackedGenParticle(cand,newRef));
123  } else {
125 
126  }
127 
128  }
129  }
130 
131 
133 
134  std::auto_ptr<edm::Association< std::vector<pat::PackedGenParticle> > > gp2pgp(new edm::Association< std::vector<pat::PackedGenParticle> > (oh ));
135  edm::Association< std::vector<pat::PackedGenParticle> >::Filler gp2pgpFiller(*gp2pgp);
136  gp2pgpFiller.insert(genOrigs, mapping.begin(), mapping.end());
137  gp2pgpFiller.fill();
138  iEvent.put(gp2pgp);
139 
140 
141 }
142 
143 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
edm::EDGetTokenT< reco::VertexCollection > PVs_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual int status() const
status word
PATPackedGenParticleProducer(const edm::ParameterSet &)
virtual double y() const
rapidity
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
key_type key() const
Accessor for product key.
Definition: Ref.h:266
virtual void produce(edm::Event &, const edm::EventSetup &)
edm::EDGetTokenT< reco::GenParticleCollection > GenOrigs_
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< reco::GenParticleCollection > Cands_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
virtual size_t numberOfMothers() const
number of mothers
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::Ref< VertexCollection > VertexRef
persistent reference to a Vertex
Definition: VertexFwd.h:13
daughters::value_type motherRef(size_type i=0) const
reference to mother at given position
edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > Asso_
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > AssoOriginal_