Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // -*- C++ -*-
0002 //
0003 // Package:    Validation/RecoVertex
0004 // Class:      BSvsPVAnalyzer
0005 //
0006 /**\class BSvsPVAnalyzer BSvsPVAnalyzer.cc Validation/RecoVertex/plugins/BSvsPVAnalyzer.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 
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/BSvsPVHistogramMaker.h"
0041 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0042 #include "DataFormats/VertexReco/interface/Vertex.h"
0043 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0044 
0045 //
0046 // class decleration
0047 //
0048 
0049 class BSvsPVAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0050 public:
0051   explicit BSvsPVAnalyzer(const edm::ParameterSet&);
0052   ~BSvsPVAnalyzer() override;
0053 
0054 private:
0055   void beginJob() override;
0056   void analyze(const edm::Event&, const edm::EventSetup&) override;
0057   void beginRun(const edm::Run&, const edm::EventSetup&) override;
0058   void endRun(const edm::Run&, const edm::EventSetup&) override;
0059   void endJob() override;
0060 
0061   // ----------member data ---------------------------
0062 
0063   BSvsPVHistogramMaker _bspvhm;
0064   edm::EDGetTokenT<reco::VertexCollection> _recoVertexCollectionToken;
0065   edm::EDGetTokenT<reco::BeamSpot> _recoBeamSpotToken;
0066   bool _firstOnly;
0067 };
0068 
0069 //
0070 // constants, enums and typedefs
0071 //
0072 
0073 //
0074 // static data member definitions
0075 //
0076 
0077 //
0078 // constructors and destructor
0079 //
0080 BSvsPVAnalyzer::BSvsPVAnalyzer(const edm::ParameterSet& iConfig)
0081     : _bspvhm(iConfig.getParameter<edm::ParameterSet>("bspvHistogramMakerPSet"), consumesCollector()),
0082       _recoVertexCollectionToken(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("pvCollection"))),
0083       _recoBeamSpotToken(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("bsCollection"))),
0084       _firstOnly(iConfig.getUntrackedParameter<bool>("firstOnly", false)) {
0085   //now do what ever initialization is needed
0086   usesResource(TFileService::kSharedResource);
0087   //
0088 
0089   _bspvhm.book();
0090 }
0091 
0092 BSvsPVAnalyzer::~BSvsPVAnalyzer() {
0093   // do anything here that needs to be done at desctruction time
0094   // (e.g. close files, deallocate resources etc.)
0095 }
0096 
0097 //
0098 // member functions
0099 //
0100 
0101 // ------------ method called to for each event  ------------
0102 void BSvsPVAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0103   // get BS
0104 
0105   edm::Handle<reco::BeamSpot> bs;
0106   iEvent.getByToken(_recoBeamSpotToken, bs);
0107 
0108   // get PV
0109 
0110   edm::Handle<reco::VertexCollection> pvcoll;
0111   iEvent.getByToken(_recoVertexCollectionToken, pvcoll);
0112 
0113   if (_firstOnly) {
0114     reco::VertexCollection firstpv;
0115     if (!pvcoll->empty())
0116       firstpv.push_back((*pvcoll)[0]);
0117     _bspvhm.fill(iEvent, firstpv, *bs);
0118   } else {
0119     _bspvhm.fill(iEvent, *pvcoll, *bs);
0120   }
0121 }
0122 
0123 // ------------ method called once each job just before starting event loop  ------------
0124 void BSvsPVAnalyzer::beginJob() {}
0125 
0126 void BSvsPVAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) { _bspvhm.beginRun(iRun.run()); }
0127 
0128 void BSvsPVAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {}
0129 // ------------ method called once each job just after ending the event loop  ------------
0130 void BSvsPVAnalyzer::endJob() {}
0131 
0132 //define this as a plug-in
0133 DEFINE_FWK_MODULE(BSvsPVAnalyzer);