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
00028 pfCandidateSrc_ = iConfig.getParameter<edm::InputTag>( "pfCandidateSource" );
00029
00030
00031 addGenMatch_ = iConfig.getParameter<bool> ( "addGenMatch" );
00032 if (addGenMatch_) {
00033 embedGenMatch_ = iConfig.getParameter<bool> ( "embedGenMatch" );
00034 if (iConfig.existsAs<edm::InputTag>("genParticleMatch")) {
00035 genMatchSrc_.push_back(iConfig.getParameter<edm::InputTag>( "genParticleMatch" ));
00036 } else {
00037 genMatchSrc_ = iConfig.getParameter<std::vector<edm::InputTag> >( "genParticleMatch" );
00038 }
00039 }
00040
00041
00042
00043 produces<std::vector<PFParticle> >();
00044
00045 }
00046
00047
00048 PATPFParticleProducer::~PATPFParticleProducer() {
00049 }
00050
00051
00052 void PATPFParticleProducer::produce(edm::Event & iEvent,
00053 const edm::EventSetup & iSetup) {
00054
00055
00056 edm::Handle<edm::View<PFParticleType> > pfCandidates;
00057
00058 fetchCandidateCollection(pfCandidates,
00059 pfCandidateSrc_,
00060 iEvent );
00061
00062
00063 std::vector<edm::Handle<edm::Association<reco::GenParticleCollection> > > genMatches(genMatchSrc_.size());
00064 if (addGenMatch_) {
00065 for (size_t j = 0, nd = genMatchSrc_.size(); j < nd; ++j) {
00066 iEvent.getByLabel(genMatchSrc_[j], genMatches[j]);
00067 }
00068 }
00069
00070
00071 std::vector<PFParticle> * patPFParticles = new std::vector<PFParticle>();
00072 for (edm::View<PFParticleType>::const_iterator
00073 itPFParticle = pfCandidates->begin();
00074 itPFParticle != pfCandidates->end();
00075 ++itPFParticle) {
00076
00077
00078 unsigned int idx = itPFParticle - pfCandidates->begin();
00079 edm::RefToBase<PFParticleType> pfCandidatesRef = pfCandidates->refAt(idx);
00080
00081 PFParticle aPFParticle(pfCandidatesRef);
00082
00083 if (addGenMatch_) {
00084 for(size_t i = 0, n = genMatches.size(); i < n; ++i) {
00085 reco::GenParticleRef genPFParticle = (*genMatches[i])[pfCandidatesRef];
00086 aPFParticle.addGenParticleRef(genPFParticle);
00087 }
00088 if (embedGenMatch_) aPFParticle.embedGenParticle();
00089 }
00090
00091
00092 patPFParticles->push_back(aPFParticle);
00093 }
00094
00095
00096 std::sort(patPFParticles->begin(), patPFParticles->end(), pTComparator_);
00097
00098
00099 std::auto_ptr<std::vector<PFParticle> > ptr(patPFParticles);
00100 iEvent.put(ptr);
00101
00102 }
00103
00104 void
00105 PATPFParticleProducer::fetchCandidateCollection( edm::Handle< edm::View<PFParticleType> >& c,
00106 const edm::InputTag& tag,
00107 const edm::Event& iEvent) const {
00108
00109 bool found = iEvent.getByLabel(tag, c);
00110
00111 if(!found ) {
00112 std::ostringstream err;
00113 err<<" cannot get PFCandidates: "
00114 <<tag<<std::endl;
00115 edm::LogError("PFCandidates")<<err.str();
00116 throw cms::Exception( "MissingProduct", err.str());
00117 }
00118
00119 }
00120
00121
00122
00123 #include "FWCore/Framework/interface/MakerMacros.h"
00124
00125 DEFINE_FWK_MODULE(PATPFParticleProducer);