Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:34

0001 // -*- C++ -*-
0002 //
0003 // Package:    Validation/RecoVertex
0004 // Class:      AnotherPrimaryVertexAnalyzer
0005 //
0006 /**\class AnotherPrimaryVertexAnalyzer AnotherPrimaryVertexAnalyzer.cc Validation/RecoVertex/plugins/AnotherPrimaryVertexAnalyzer.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Andrea Venturi
0015 //         Created:  Mon Oct 27 17:37:53 CET 2008
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 
0024 #include <vector>
0025 #include <map>
0026 #include <limits>
0027 #include <string>
0028 
0029 #include "FWCore/Framework/interface/Frameworkfwd.h"
0030 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0031 #include "FWCore/Framework/interface/ConsumesCollector.h"
0032 #include "FWCore/Framework/interface/Event.h"
0033 #include "FWCore/Framework/interface/Run.h"
0034 #include "FWCore/Framework/interface/MakerMacros.h"
0035 
0036 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0037 
0038 #include "FWCore/Utilities/interface/InputTag.h"
0039 
0040 #include "Validation/RecoVertex/interface/VertexHistogramMaker.h"
0041 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0042 #include "DataFormats/VertexReco/interface/Vertex.h"
0043 
0044 #include "CommonTools/TriggerUtils/interface/PrescaleWeightProvider.h"
0045 
0046 //
0047 // class decleration
0048 //
0049 
0050 class AnotherPrimaryVertexAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0051 public:
0052   explicit AnotherPrimaryVertexAnalyzer(const edm::ParameterSet&);
0053   ~AnotherPrimaryVertexAnalyzer() override;
0054 
0055 private:
0056   void beginJob() override;
0057   void analyze(const edm::Event&, const edm::EventSetup&) override;
0058   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0059   void endRun(const edm::Run&, const edm::EventSetup&) override;
0060   void endJob() override;
0061 
0062   // ----------member data ---------------------------
0063 
0064   VertexHistogramMaker _vhm;
0065   edm::EDGetTokenT<reco::VertexCollection> _recoVertexCollectionToken;
0066   bool _firstOnly;
0067 
0068   std::unique_ptr<PrescaleWeightProvider> _weightprov;
0069 };
0070 
0071 //
0072 // constants, enums and typedefs
0073 //
0074 
0075 //
0076 // static data member definitions
0077 //
0078 
0079 //
0080 // constructors and destructor
0081 //
0082 AnotherPrimaryVertexAnalyzer::AnotherPrimaryVertexAnalyzer(const edm::ParameterSet& iConfig)
0083     : _vhm(iConfig.getParameter<edm::ParameterSet>("vHistogramMakerPSet"), consumesCollector()),
0084       _recoVertexCollectionToken(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("pvCollection"))),
0085       _firstOnly(iConfig.getUntrackedParameter<bool>("firstOnly", false)),
0086       _weightprov(
0087           iConfig.getParameter<bool>("usePrescaleWeight")
0088               ? new PrescaleWeightProvider(
0089                     iConfig.getParameter<edm::ParameterSet>("prescaleWeightProviderPSet"), consumesCollector(), *this)
0090               : nullptr) {
0091   //now do what ever initialization is needed
0092   usesResource(TFileService::kSharedResource);
0093   //
0094   _vhm.book();
0095 }
0096 
0097 AnotherPrimaryVertexAnalyzer::~AnotherPrimaryVertexAnalyzer() {}
0098 
0099 //
0100 // member functions
0101 //
0102 
0103 // ------------ method called to for each event  ------------
0104 void AnotherPrimaryVertexAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0105   // compute event weight
0106   auto const weight = _weightprov ? _weightprov->prescaleWeight<double>(iEvent, iSetup) : 1.;
0107 
0108   // get PV
0109   edm::Handle<reco::VertexCollection> pvcoll;
0110   iEvent.getByToken(_recoVertexCollectionToken, pvcoll);
0111 
0112   if (_firstOnly) {
0113     reco::VertexCollection firstpv;
0114     if (!pvcoll->empty())
0115       firstpv.push_back((*pvcoll)[0]);
0116     _vhm.fill(iEvent, firstpv, weight);
0117   } else {
0118     _vhm.fill(iEvent, *pvcoll, weight);
0119   }
0120 }
0121 
0122 // ------------ method called once each job just before starting event loop  ------------
0123 void AnotherPrimaryVertexAnalyzer::beginJob() {}
0124 
0125 void AnotherPrimaryVertexAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0126   _vhm.beginRun(iRun);
0127 
0128   if (_weightprov)
0129     _weightprov->initRun(iRun, iSetup);
0130 }
0131 
0132 void AnotherPrimaryVertexAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {}
0133 // ------------ method called once each job just after ending the event loop  ------------
0134 void AnotherPrimaryVertexAnalyzer::endJob() {}
0135 
0136 //define this as a plug-in
0137 DEFINE_FWK_MODULE(AnotherPrimaryVertexAnalyzer);