Go to the documentation of this file.00001
00002
00003
00004
00005 #include "PhysicsTools/PatAlgos/plugins/PATPFParticleProducer.h"
00006
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 #include "FWCore/ParameterSet/interface/FileInPath.h"
00009
00010 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00011 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00012 #include "DataFormats/HepMCCandidate/interface/GenParticleFwd.h"
00013 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00014
00015 #include "DataFormats/Common/interface/Association.h"
00016
00017 #include "TMath.h"
00018
00019 #include <vector>
00020 #include <memory>
00021
00022
00023 using namespace pat;
00024
00025
00026 PATPFParticleProducer::PATPFParticleProducer(const edm::ParameterSet & iConfig) :
00027 userDataHelper_ ( iConfig.getParameter<edm::ParameterSet>("userData") )
00028 {
00029
00030 pfCandidateSrc_ = iConfig.getParameter<edm::InputTag>( "pfCandidateSource" );
00031
00032
00033 addGenMatch_ = iConfig.getParameter<bool> ( "addGenMatch" );
00034 if (addGenMatch_) {
00035 embedGenMatch_ = iConfig.getParameter<bool> ( "embedGenMatch" );
00036 if (iConfig.existsAs<edm::InputTag>("genParticleMatch")) {
00037 genMatchSrc_.push_back(iConfig.getParameter<edm::InputTag>( "genParticleMatch" ));
00038 } else {
00039 genMatchSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >( "genParticleMatch" );
00040 }
00041 }
00042
00043
00044 addEfficiencies_ = iConfig.getParameter<bool>("addEfficiencies");
00045 if (addEfficiencies_) {
00046 efficiencyLoader_ = pat::helper::EfficiencyLoader(iConfig.getParameter<edm::ParameterSet>("efficiencies"));
00047 }
00048
00049
00050 addResolutions_ = iConfig.getParameter<bool>("addResolutions");
00051 if (addResolutions_) {
00052 resolutionLoader_ = pat::helper::KinResolutionsLoader(iConfig.getParameter<edm::ParameterSet>("resolutions"));
00053 }
00054
00055
00056 useUserData_ = false;
00057 if ( iConfig.exists("userData") ) {
00058 useUserData_ = true;
00059 }
00060
00061
00062 produces<std::vector<PFParticle> >();
00063
00064 }
00065
00066
00067 PATPFParticleProducer::~PATPFParticleProducer() {
00068 }
00069
00070
00071 void PATPFParticleProducer::produce(edm::Event & iEvent,
00072 const edm::EventSetup & iSetup) {
00073
00074
00075 edm::Handle<edm::View<reco::PFCandidate> > pfCandidates;
00076
00077 fetchCandidateCollection(pfCandidates,
00078 pfCandidateSrc_,
00079 iEvent );
00080
00081
00082 std::vector<edm::Handle<edm::Association<reco::GenParticleCollection> > > genMatches(genMatchSrc_.size());
00083 if (addGenMatch_) {
00084 for (size_t j = 0, nd = genMatchSrc_.size(); j < nd; ++j) {
00085 iEvent.getByLabel(genMatchSrc_[j], genMatches[j]);
00086 }
00087 }
00088
00089 if (efficiencyLoader_.enabled()) efficiencyLoader_.newEvent(iEvent);
00090 if (resolutionLoader_.enabled()) resolutionLoader_.newEvent(iEvent, iSetup);
00091
00092
00093 std::vector<PFParticle> * patPFParticles = new std::vector<PFParticle>();
00094 for (edm::View<reco::PFCandidate>::const_iterator
00095 itPFParticle = pfCandidates->begin();
00096 itPFParticle != pfCandidates->end();
00097 ++itPFParticle) {
00098
00099
00100 unsigned int idx = itPFParticle - pfCandidates->begin();
00101 edm::RefToBase<reco::PFCandidate> pfCandidatesRef = pfCandidates->refAt(idx);
00102
00103 PFParticle aPFParticle(pfCandidatesRef);
00104
00105 if (addGenMatch_) {
00106 for(size_t i = 0, n = genMatches.size(); i < n; ++i) {
00107 reco::GenParticleRef genPFParticle = (*genMatches[i])[pfCandidatesRef];
00108 aPFParticle.addGenParticleRef(genPFParticle);
00109 }
00110 if (embedGenMatch_) aPFParticle.embedGenParticle();
00111 }
00112
00113 if (efficiencyLoader_.enabled()) {
00114 efficiencyLoader_.setEfficiencies( aPFParticle, pfCandidatesRef );
00115 }
00116
00117 if (resolutionLoader_.enabled()) {
00118 resolutionLoader_.setResolutions(aPFParticle);
00119 }
00120
00121 if ( useUserData_ ) {
00122 userDataHelper_.add( aPFParticle, iEvent, iSetup );
00123 }
00124
00125
00126
00127 patPFParticles->push_back(aPFParticle);
00128 }
00129
00130
00131 std::sort(patPFParticles->begin(), patPFParticles->end(), pTComparator_);
00132
00133
00134 std::auto_ptr<std::vector<PFParticle> > ptr(patPFParticles);
00135 iEvent.put(ptr);
00136
00137 }
00138
00139 void
00140 PATPFParticleProducer::fetchCandidateCollection( edm::Handle<edm::View<reco::PFCandidate> >& c,
00141 const edm::InputTag& tag,
00142 const edm::Event& iEvent) const {
00143
00144 bool found = iEvent.getByLabel(tag, c);
00145
00146 if(!found ) {
00147 std::ostringstream err;
00148 err<<" cannot get PFCandidates: "
00149 <<tag<<std::endl;
00150 edm::LogError("PFCandidates")<<err.str();
00151 throw cms::Exception( "MissingProduct", err.str());
00152 }
00153
00154 }
00155
00156
00157
00158 #include "FWCore/Framework/interface/MakerMacros.h"
00159
00160 DEFINE_FWK_MODULE(PATPFParticleProducer);