CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_6/src/PhysicsTools/PatAlgos/plugins/PATPFParticleProducer.cc

Go to the documentation of this file.
00001 //
00002 // $Id: PATPFParticleProducer.cc,v 1.5 2009/06/25 23:49:35 gpetrucc Exp $
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   // general configurables
00028   pfCandidateSrc_ = iConfig.getParameter<edm::InputTag>( "pfCandidateSource" );
00029  
00030   // MC matching configurables
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   // Efficiency configurables
00042   addEfficiencies_ = iConfig.getParameter<bool>("addEfficiencies");
00043   if (addEfficiencies_) {
00044      efficiencyLoader_ = pat::helper::EfficiencyLoader(iConfig.getParameter<edm::ParameterSet>("efficiencies"));
00045   }
00046 
00047   // Resolution configurables
00048   addResolutions_ = iConfig.getParameter<bool>("addResolutions");
00049   if (addResolutions_) {
00050      resolutionLoader_ = pat::helper::KinResolutionsLoader(iConfig.getParameter<edm::ParameterSet>("resolutions"));
00051   }
00052 
00053 
00054 
00055   // produces vector of muons
00056   produces<std::vector<PFParticle> >();
00057 
00058 }
00059 
00060 
00061 PATPFParticleProducer::~PATPFParticleProducer() {
00062 }
00063 
00064 
00065 void PATPFParticleProducer::produce(edm::Event & iEvent, 
00066                                     const edm::EventSetup & iSetup) {
00067 
00068   // Get the collection of PFCandidates from the event
00069   edm::Handle<edm::View<reco::PFCandidate> > pfCandidates;
00070 
00071   fetchCandidateCollection(pfCandidates, 
00072                            pfCandidateSrc_, 
00073                            iEvent );
00074 
00075   // prepare the MC matching
00076   std::vector<edm::Handle<edm::Association<reco::GenParticleCollection> > > genMatches(genMatchSrc_.size());
00077   if (addGenMatch_) {
00078         for (size_t j = 0, nd = genMatchSrc_.size(); j < nd; ++j) {
00079             iEvent.getByLabel(genMatchSrc_[j], genMatches[j]);
00080         }
00081   }
00082 
00083   if (efficiencyLoader_.enabled()) efficiencyLoader_.newEvent(iEvent);
00084   if (resolutionLoader_.enabled()) resolutionLoader_.newEvent(iEvent, iSetup);
00085 
00086   // loop over PFCandidates
00087   std::vector<PFParticle> * patPFParticles = new std::vector<PFParticle>();
00088   for (edm::View<reco::PFCandidate>::const_iterator 
00089          itPFParticle = pfCandidates->begin(); 
00090        itPFParticle != pfCandidates->end(); 
00091        ++itPFParticle) {
00092 
00093     // construct the PFParticle from the ref -> save ref to original object
00094     unsigned int idx = itPFParticle - pfCandidates->begin();
00095     edm::RefToBase<reco::PFCandidate> pfCandidatesRef = pfCandidates->refAt(idx);
00096 
00097     PFParticle aPFParticle(pfCandidatesRef);
00098 
00099     if (addGenMatch_) {
00100       for(size_t i = 0, n = genMatches.size(); i < n; ++i) {
00101           reco::GenParticleRef genPFParticle = (*genMatches[i])[pfCandidatesRef];
00102           aPFParticle.addGenParticleRef(genPFParticle);
00103       }
00104       if (embedGenMatch_) aPFParticle.embedGenParticle();
00105     }
00106       
00107     if (efficiencyLoader_.enabled()) {
00108         efficiencyLoader_.setEfficiencies( aPFParticle, pfCandidatesRef );
00109     }
00110 
00111     if (resolutionLoader_.enabled()) {
00112         resolutionLoader_.setResolutions(aPFParticle);
00113     }
00114 
00115     // add sel to selected
00116     patPFParticles->push_back(aPFParticle);
00117   }
00118 
00119   // sort pfCandidates in pt
00120   std::sort(patPFParticles->begin(), patPFParticles->end(), pTComparator_);
00121 
00122   // put genEvt object in Event
00123   std::auto_ptr<std::vector<PFParticle> > ptr(patPFParticles);
00124   iEvent.put(ptr);
00125 
00126 }
00127 
00128 void 
00129 PATPFParticleProducer::fetchCandidateCollection( edm::Handle<edm::View<reco::PFCandidate> >& c, 
00130                                                  const edm::InputTag& tag, 
00131                                                  const edm::Event& iEvent) const {
00132   
00133   bool found = iEvent.getByLabel(tag, c);
00134   
00135   if(!found ) {
00136     std::ostringstream  err;
00137     err<<" cannot get PFCandidates: "
00138        <<tag<<std::endl;
00139     edm::LogError("PFCandidates")<<err.str();
00140     throw cms::Exception( "MissingProduct", err.str());
00141   }
00142   
00143 }
00144 
00145 
00146 
00147 #include "FWCore/Framework/interface/MakerMacros.h"
00148 
00149 DEFINE_FWK_MODULE(PATPFParticleProducer);