CMS 3D CMS Logo

PATPFParticleProducer.cc
Go to the documentation of this file.
1 //
2 //
3 
5 
8 
13 
15 
16 #include "TMath.h"
17 
18 #include <vector>
19 #include <memory>
20 
21 
22 using namespace pat;
23 
24 
26  userDataHelper_ ( iConfig.getParameter<edm::ParameterSet>("userData"), consumesCollector() )
27 {
28  // general configurables
29  pfCandidateToken_ = consumes<edm::View<reco::PFCandidate> >(iConfig.getParameter<edm::InputTag>( "pfCandidateSource" ));
30 
31  // MC matching configurables
32  addGenMatch_ = iConfig.getParameter<bool> ( "addGenMatch" );
33  if (addGenMatch_) {
34  embedGenMatch_ = iConfig.getParameter<bool>( "embedGenMatch" );
35  if (iConfig.existsAs<edm::InputTag>("genParticleMatch")) {
37  }
38  else {
39  genMatchTokens_ = edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag> >( "genParticleMatch" ), [this](edm::InputTag const & tag){return consumes<edm::Association<reco::GenParticleCollection> >(tag);});
40  }
41  }
42 
43  // Efficiency configurables
44  addEfficiencies_ = iConfig.getParameter<bool>("addEfficiencies");
45  if (addEfficiencies_) {
47  }
48 
49  // Resolution configurables
50  addResolutions_ = iConfig.getParameter<bool>("addResolutions");
51  if (addResolutions_) {
53  }
54 
55  // Check to see if the user wants to add user data
56  useUserData_ = false;
57  if ( iConfig.exists("userData") ) {
58  useUserData_ = true;
59  }
60 
61  // produces vector of muons
62  produces<std::vector<PFParticle> >();
63 
64 }
65 
66 
68 }
69 
70 
72  const edm::EventSetup & iSetup) {
73 
74  // Get the collection of PFCandidates from the event
76  iEvent.getByToken(pfCandidateToken_, pfCandidates);
77 
78  // prepare the MC matching
79  std::vector<edm::Handle<edm::Association<reco::GenParticleCollection> > > genMatches(genMatchTokens_.size());
80  if (addGenMatch_) {
81  for (size_t j = 0, nd = genMatchTokens_.size(); j < nd; ++j) {
82  iEvent.getByToken(genMatchTokens_[j], genMatches[j]);
83  }
84  }
85 
87  if (resolutionLoader_.enabled()) resolutionLoader_.newEvent(iEvent, iSetup);
88 
89  // loop over PFCandidates
90  std::vector<PFParticle> * patPFParticles = new std::vector<PFParticle>();
92  itPFParticle = pfCandidates->begin();
93  itPFParticle != pfCandidates->end();
94  ++itPFParticle) {
95 
96  // construct the PFParticle from the ref -> save ref to original object
97  unsigned int idx = itPFParticle - pfCandidates->begin();
98  edm::RefToBase<reco::PFCandidate> pfCandidatesRef = pfCandidates->refAt(idx);
99 
100  PFParticle aPFParticle(pfCandidatesRef);
101 
102  if (addGenMatch_) {
103  for(size_t i = 0, n = genMatches.size(); i < n; ++i) {
104  reco::GenParticleRef genPFParticle = (*genMatches[i])[pfCandidatesRef];
105  aPFParticle.addGenParticleRef(genPFParticle);
106  }
107  if (embedGenMatch_) aPFParticle.embedGenParticle();
108  }
109 
110  if (efficiencyLoader_.enabled()) {
111  efficiencyLoader_.setEfficiencies( aPFParticle, pfCandidatesRef );
112  }
113 
114  if (resolutionLoader_.enabled()) {
115  resolutionLoader_.setResolutions(aPFParticle);
116  }
117 
118  if ( useUserData_ ) {
119  userDataHelper_.add( aPFParticle, iEvent, iSetup );
120  }
121 
122 
123  // add sel to selected
124  patPFParticles->push_back(aPFParticle);
125  }
126 
127  // sort pfCandidates in pt
128  std::sort(patPFParticles->begin(), patPFParticles->end(), pTComparator_);
129 
130  // put genEvt object in Event
131  std::unique_ptr<std::vector<PFParticle> > ptr(patPFParticles);
132  iEvent.put(std::move(ptr));
133 
134 }
135 
136 
137 
139 
bool enabled() const
&#39;true&#39; if this there is at least one efficiency configured
T getParameter(std::string const &) const
void newEvent(const edm::Event &event)
To be called for each new event, reads in the ValueMaps for efficiencies.
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:136
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
pat::helper::KinResolutionsLoader resolutionLoader_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< edm::EDGetTokenT< edm::Association< reco::GenParticleCollection > > > genMatchTokens_
bool exists(std::string const &parameterName) const
checks if a parameter exists
PATPFParticleProducer(const edm::ParameterSet &iConfig)
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
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
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
Definition: HeavyIon.h:7
int iEvent
Definition: GenABIO.cc:230
void embedGenParticle()
Definition: PATObject.h:694
void newEvent(const edm::Event &event, const edm::EventSetup &setup)
To be called for each new event, reads in the EventSetup object.
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
pat::helper::EfficiencyLoader efficiencyLoader_
edm::EDGetTokenT< edm::View< reco::PFCandidate > > pfCandidateToken_
Produces pat::PFParticle&#39;s.
void addGenParticleRef(const reco::GenParticleRef &ref)
Definition: PATObject.h:678
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:35
HLT enums.
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:86
pat::PATUserDataHelper< pat::PFParticle > userDataHelper_
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) override
def move(src, dest)
Definition: eostools.py:510
GreaterByPt< PFParticle > pTComparator_