Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:01:06

0001 // system include files
0002 #include <memory>
0003 #include <string>
0004 
0005 // user include files
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/stream/EDProducer.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009 
0010 #include "FWCore/Framework/interface/Event.h"
0011 #include "FWCore/Framework/interface/MakerMacros.h"
0012 
0013 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0014 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0015 #include "DataFormats/Common/interface/Association.h"
0016 
0017 #include "CommonTools/ParticleFlow/interface/PFPileUpAlgo.h"
0018 
0019 #include "DataFormats/ParticleFlowCandidate/interface/PileUpPFCandidate.h"
0020 #include "DataFormats/ParticleFlowCandidate/interface/PileUpPFCandidateFwd.h"
0021 #include "DataFormats/VertexReco/interface/Vertex.h"
0022 
0023 #include "FWCore/Framework/interface/ESHandle.h"
0024 
0025 #include "FWCore/Utilities/interface/Exception.h"
0026 #include "FWCore/Framework/interface/EventSetup.h"
0027 
0028 using namespace std;
0029 using namespace edm;
0030 using namespace reco;
0031 
0032 /**\class PFPileUp
0033 \brief Identifies pile-up candidates from a collection of PFCandidates, and
0034 produces the corresponding collection of PileUpCandidates.
0035 
0036 \author Colin Bernet
0037 \date   february 2008
0038 \updated Florian Beaudette 30/03/2012
0039 
0040 */
0041 
0042 class PFPileUp : public edm::stream::EDProducer<> {
0043 public:
0044   typedef std::vector<edm::FwdPtr<reco::PFCandidate>> PFCollection;
0045   typedef edm::View<reco::PFCandidate> PFView;
0046   typedef std::vector<reco::PFCandidate> PFCollectionByValue;
0047   typedef edm::Association<reco::VertexCollection> CandToVertex;
0048 
0049   explicit PFPileUp(const edm::ParameterSet&);
0050 
0051   ~PFPileUp() override;
0052 
0053   void produce(edm::Event&, const edm::EventSetup&) override;
0054 
0055   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0056 
0057 private:
0058   PFPileUpAlgo pileUpAlgo_;
0059 
0060   /// PFCandidates to be analyzed
0061   edm::EDGetTokenT<PFCollection> tokenPFCandidates_;
0062   /// fall-back token
0063   edm::EDGetTokenT<PFView> tokenPFCandidatesView_;
0064 
0065   /// vertices
0066   edm::EDGetTokenT<reco::VertexCollection> tokenVertices_;
0067 
0068   /// enable PFPileUp selection
0069   bool enable_;
0070 
0071   /// verbose ?
0072   bool verbose_;
0073 
0074   /// use the closest z vertex if a track is not in a vertex
0075   bool checkClosestZVertex_;
0076 
0077   edm::EDGetTokenT<CandToVertex> tokenVertexAssociation_;
0078   edm::EDGetTokenT<edm::ValueMap<int>> tokenVertexAssociationQuality_;
0079   bool fUseVertexAssociation;
0080   int vertexAssociationQuality_;
0081   unsigned int fNumOfPUVtxsForCharged_;
0082   double fDzCutForChargedFromPUVtxs_;
0083 };
0084 
0085 PFPileUp::PFPileUp(const edm::ParameterSet& iConfig) {
0086   tokenPFCandidates_ = consumes<PFCollection>(iConfig.getParameter<InputTag>("PFCandidates"));
0087   tokenPFCandidatesView_ = mayConsume<PFView>(iConfig.getParameter<InputTag>("PFCandidates"));
0088 
0089   tokenVertices_ = consumes<VertexCollection>(iConfig.getParameter<InputTag>("Vertices"));
0090 
0091   fUseVertexAssociation = iConfig.getParameter<bool>("useVertexAssociation");
0092   vertexAssociationQuality_ = iConfig.getParameter<int>("vertexAssociationQuality");
0093   if (fUseVertexAssociation) {
0094     tokenVertexAssociation_ = consumes<CandToVertex>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
0095     tokenVertexAssociationQuality_ =
0096         consumes<edm::ValueMap<int>>(iConfig.getParameter<edm::InputTag>("vertexAssociation"));
0097   }
0098   fNumOfPUVtxsForCharged_ = iConfig.getParameter<unsigned int>("NumOfPUVtxsForCharged");
0099   fDzCutForChargedFromPUVtxs_ = iConfig.getParameter<double>("DzCutForChargedFromPUVtxs");
0100 
0101   enable_ = iConfig.getParameter<bool>("enable");
0102 
0103   verbose_ = iConfig.getUntrackedParameter<bool>("verbose", false);
0104 
0105   checkClosestZVertex_ = iConfig.getParameter<bool>("checkClosestZVertex");
0106 
0107   // Configure the algo
0108   pileUpAlgo_.setVerbose(verbose_);
0109   pileUpAlgo_.setCheckClosestZVertex(checkClosestZVertex_);
0110   pileUpAlgo_.setNumOfPUVtxsForCharged(fNumOfPUVtxsForCharged_);
0111   pileUpAlgo_.setDzCutForChargedFromPUVtxs(fDzCutForChargedFromPUVtxs_);
0112 
0113   //produces<reco::PFCandidateCollection>();
0114   produces<PFCollection>();
0115   // produces< PFCollectionByValue > ();
0116 }
0117 
0118 PFPileUp::~PFPileUp() {}
0119 
0120 void PFPileUp::produce(Event& iEvent, const EventSetup& iSetup) {
0121   //   LogDebug("PFPileUp")<<"START event: "<<iEvent.id().event()
0122   //             <<" in run "<<iEvent.id().run()<<endl;
0123 
0124   // get PFCandidates
0125 
0126   unique_ptr<PFCollection> pOutput(new PFCollection);
0127 
0128   unique_ptr<PFCollectionByValue> pOutputByValue(new PFCollectionByValue);
0129 
0130   if (enable_) {
0131     // get vertices
0132     Handle<VertexCollection> vertices;
0133     iEvent.getByToken(tokenVertices_, vertices);
0134 
0135     // get PF Candidates
0136     Handle<PFCollection> pfCandidates;
0137     PFCollection const* pfCandidatesRef = nullptr;
0138     PFCollection usedIfNoFwdPtrs;
0139     bool getFromFwdPtr = iEvent.getByToken(tokenPFCandidates_, pfCandidates);
0140     if (getFromFwdPtr) {
0141       pfCandidatesRef = pfCandidates.product();
0142     }
0143     // Maintain backwards-compatibility.
0144     // If there is no vector of FwdPtr<PFCandidate> found, then
0145     // make a dummy vector<FwdPtr<PFCandidate> > for the PFPileupAlgo,
0146     // set the pointer "pfCandidatesRef" to point to it, and
0147     // then we can pass it to the PFPileupAlgo.
0148     else {
0149       Handle<PFView> pfView;
0150       bool getFromView = iEvent.getByToken(tokenPFCandidatesView_, pfView);
0151       if (!getFromView) {
0152         throw cms::Exception(
0153             "PFPileUp is misconfigured. This needs to be either vector<FwdPtr<PFCandidate> >, or View<PFCandidate>");
0154       }
0155       for (edm::View<reco::PFCandidate>::const_iterator viewBegin = pfView->begin(),
0156                                                         viewEnd = pfView->end(),
0157                                                         iview = viewBegin;
0158            iview != viewEnd;
0159            ++iview) {
0160         usedIfNoFwdPtrs.push_back(
0161             edm::FwdPtr<reco::PFCandidate>(pfView->ptrAt(iview - viewBegin), pfView->ptrAt(iview - viewBegin)));
0162       }
0163       pfCandidatesRef = &usedIfNoFwdPtrs;
0164     }
0165 
0166     if (pfCandidatesRef == nullptr) {
0167       throw cms::Exception(
0168           "Something went dreadfully wrong with PFPileUp. pfCandidatesRef should never be zero, so this is a logic "
0169           "error.");
0170     }
0171 
0172     if (fUseVertexAssociation) {
0173       const edm::Association<reco::VertexCollection>& associatedPV = iEvent.get(tokenVertexAssociation_);
0174       const edm::ValueMap<int>& associationQuality = iEvent.get(tokenVertexAssociationQuality_);
0175       PFCollection pfCandidatesFromPU;
0176       for (auto& p : (*pfCandidatesRef)) {
0177         const reco::VertexRef& PVOrig = associatedPV[p];
0178         int quality = associationQuality[p];
0179         if (PVOrig.isNonnull() && (PVOrig.key() > 0) && (quality >= vertexAssociationQuality_))
0180           pfCandidatesFromPU.push_back(p);
0181       }
0182       pOutput->insert(pOutput->end(), pfCandidatesFromPU.begin(), pfCandidatesFromPU.end());
0183     } else {
0184       pileUpAlgo_.process(*pfCandidatesRef, *vertices);
0185       pOutput->insert(
0186           pOutput->end(), pileUpAlgo_.getPFCandidatesFromPU().begin(), pileUpAlgo_.getPFCandidatesFromPU().end());
0187     }
0188     // for ( PFCollection::const_iterator byValueBegin = pileUpAlgo_.getPFCandidatesFromPU().begin(),
0189     //      byValueEnd = pileUpAlgo_.getPFCandidatesFromPU().end(), ibyValue = byValueBegin;
0190     //    ibyValue != byValueEnd; ++ibyValue ) {
0191     //   pOutputByValue->push_back( **ibyValue );
0192     // }
0193 
0194   }  // end if enabled
0195   // outsize of the loop to fill the collection anyway even when disabled
0196   iEvent.put(std::move(pOutput));
0197   // iEvent.put(std::move(pOutputByValue));
0198 }
0199 
0200 void PFPileUp::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0201   edm::ParameterSetDescription desc;
0202   desc.add<edm::InputTag>("PFCandidates", edm::InputTag("particleFlowTmpPtrs"));
0203   desc.add<edm::InputTag>("Vertices", edm::InputTag("offlinePrimaryVertices"));
0204   desc.add<bool>("enable", true);
0205   desc.addUntracked<bool>("verbose", false);
0206   desc.add<bool>("checkClosestZVertex", true);
0207   desc.add<bool>("useVertexAssociation", false);
0208   desc.add<int>("vertexAssociationQuality", 0);
0209   desc.add<edm::InputTag>("vertexAssociation", edm::InputTag(""));
0210   desc.add<unsigned int>("NumOfPUVtxsForCharged", 0);
0211   desc.add<double>("DzCutForChargedFromPUVtxs", .2);
0212   descriptions.add("pfPileUp", desc);
0213 }
0214 
0215 DEFINE_FWK_MODULE(PFPileUp);