CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PATPFParticleProducer.cc
Go to the documentation of this file.
1 //
2 // $Id: PATPFParticleProducer.cc,v 1.5 2009/06/25 23:49:35 gpetrucc Exp $
3 //
4 
6 
9 
14 
16 
17 #include "TMath.h"
18 
19 #include <vector>
20 #include <memory>
21 
22 
23 using namespace pat;
24 
25 
27  // general configurables
28  pfCandidateSrc_ = iConfig.getParameter<edm::InputTag>( "pfCandidateSource" );
29 
30  // MC matching configurables
31  addGenMatch_ = iConfig.getParameter<bool> ( "addGenMatch" );
32  if (addGenMatch_) {
33  embedGenMatch_ = iConfig.getParameter<bool> ( "embedGenMatch" );
34  if (iConfig.existsAs<edm::InputTag>("genParticleMatch")) {
35  genMatchSrc_.push_back(iConfig.getParameter<edm::InputTag>( "genParticleMatch" ));
36  } else {
37  genMatchSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >( "genParticleMatch" );
38  }
39  }
40 
41  // Efficiency configurables
42  addEfficiencies_ = iConfig.getParameter<bool>("addEfficiencies");
43  if (addEfficiencies_) {
45  }
46 
47  // Resolution configurables
48  addResolutions_ = iConfig.getParameter<bool>("addResolutions");
49  if (addResolutions_) {
51  }
52 
53 
54 
55  // produces vector of muons
56  produces<std::vector<PFParticle> >();
57 
58 }
59 
60 
62 }
63 
64 
66  const edm::EventSetup & iSetup) {
67 
68  // Get the collection of PFCandidates from the event
70 
71  fetchCandidateCollection(pfCandidates,
73  iEvent );
74 
75  // prepare the MC matching
76  std::vector<edm::Handle<edm::Association<reco::GenParticleCollection> > > genMatches(genMatchSrc_.size());
77  if (addGenMatch_) {
78  for (size_t j = 0, nd = genMatchSrc_.size(); j < nd; ++j) {
79  iEvent.getByLabel(genMatchSrc_[j], genMatches[j]);
80  }
81  }
82 
84  if (resolutionLoader_.enabled()) resolutionLoader_.newEvent(iEvent, iSetup);
85 
86  // loop over PFCandidates
87  std::vector<PFParticle> * patPFParticles = new std::vector<PFParticle>();
89  itPFParticle = pfCandidates->begin();
90  itPFParticle != pfCandidates->end();
91  ++itPFParticle) {
92 
93  // construct the PFParticle from the ref -> save ref to original object
94  unsigned int idx = itPFParticle - pfCandidates->begin();
95  edm::RefToBase<reco::PFCandidate> pfCandidatesRef = pfCandidates->refAt(idx);
96 
97  PFParticle aPFParticle(pfCandidatesRef);
98 
99  if (addGenMatch_) {
100  for(size_t i = 0, n = genMatches.size(); i < n; ++i) {
101  reco::GenParticleRef genPFParticle = (*genMatches[i])[pfCandidatesRef];
102  aPFParticle.addGenParticleRef(genPFParticle);
103  }
104  if (embedGenMatch_) aPFParticle.embedGenParticle();
105  }
106 
107  if (efficiencyLoader_.enabled()) {
108  efficiencyLoader_.setEfficiencies( aPFParticle, pfCandidatesRef );
109  }
110 
111  if (resolutionLoader_.enabled()) {
112  resolutionLoader_.setResolutions(aPFParticle);
113  }
114 
115  // add sel to selected
116  patPFParticles->push_back(aPFParticle);
117  }
118 
119  // sort pfCandidates in pt
120  std::sort(patPFParticles->begin(), patPFParticles->end(), pTComparator_);
121 
122  // put genEvt object in Event
123  std::auto_ptr<std::vector<PFParticle> > ptr(patPFParticles);
124  iEvent.put(ptr);
125 
126 }
127 
128 void
130  const edm::InputTag& tag,
131  const edm::Event& iEvent) const {
132 
133  bool found = iEvent.getByLabel(tag, c);
134 
135  if(!found ) {
136  std::ostringstream err;
137  err<<" cannot get PFCandidates: "
138  <<tag<<std::endl;
139  edm::LogError("PFCandidates")<<err.str();
140  throw cms::Exception( "MissingProduct", err.str());
141  }
142 
143 }
144 
145 
146 
148 
bool enabled() const
&#39;true&#39; if this there is at least one efficiency configured
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:187
pat::helper::KinResolutionsLoader resolutionLoader_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
PATPFParticleProducer(const edm::ParameterSet &iConfig)
void setResolutions(pat::PATObject< T > &obj) const
Sets the efficiencies for this object, using the reference to the original objects.
bool enabled() const
&#39;true&#39; if this there is at least one efficiency configured
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
std::vector< edm::InputTag > genMatchSrc_
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
void embedGenParticle()
Definition: PATObject.h:675
int j
Definition: DBlmapReader.cc:9
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
pat::helper::EfficiencyLoader efficiencyLoader_
Produces pat::PFParticle&#39;s.
void addGenParticleRef(const reco::GenParticleRef &ref)
Definition: PATObject.h:659
void setEfficiencies(pat::PATObject< T > &obj, const R &originalRef) const
Sets the efficiencies for this object, using the reference to the original objects.
Analysis-level class for reconstructed particles.
Definition: PFParticle.h:37
void fetchCandidateCollection(edm::Handle< edm::View< reco::PFCandidate > > &c, const edm::InputTag &tag, const edm::Event &iSetup) const
void newEvent(const edm::Event &event, const edm::EventSetup &setup) const
To be called for each new event, reads in the EventSetup object.
void newEvent(const edm::Event &event) const
To be called for each new event, reads in the ValueMaps for efficiencies.
GreaterByPt< PFParticle > pTComparator_