CMS 3D CMS Logo

BSvsPVAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Validation/RecoVertex
4 // Class: BSvsPVAnalyzer
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Mon Oct 27 17:37:53 CET 2008
16 //
17 //
18 
19 // system include files
20 #include <memory>
21 
22 // user include files
23 
24 #include <vector>
25 #include <map>
26 #include <limits>
27 #include <string>
28 
31 
35 
37 
39 
44 
45 //
46 // class decleration
47 //
48 
49 class BSvsPVAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
50 public:
51  explicit BSvsPVAnalyzer(const edm::ParameterSet&);
52  ~BSvsPVAnalyzer() override;
53 
54 private:
55  void beginJob() override;
56  void analyze(const edm::Event&, const edm::EventSetup&) override;
57  void beginRun(const edm::Run&, const edm::EventSetup&) override;
58  void endRun(const edm::Run&, const edm::EventSetup&) override;
59  void endJob() override;
60 
61  // ----------member data ---------------------------
62 
66  bool _firstOnly;
67 };
68 
69 //
70 // constants, enums and typedefs
71 //
72 
73 //
74 // static data member definitions
75 //
76 
77 //
78 // constructors and destructor
79 //
81  : _bspvhm(iConfig.getParameter<edm::ParameterSet>("bspvHistogramMakerPSet"), consumesCollector()),
82  _recoVertexCollectionToken(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("pvCollection"))),
83  _recoBeamSpotToken(consumes<reco::BeamSpot>(iConfig.getParameter<edm::InputTag>("bsCollection"))),
84  _firstOnly(iConfig.getUntrackedParameter<bool>("firstOnly", false)) {
85  //now do what ever initialization is needed
86  usesResource(TFileService::kSharedResource);
87  //
88 
89  _bspvhm.book();
90 }
91 
93  // do anything here that needs to be done at desctruction time
94  // (e.g. close files, deallocate resources etc.)
95 }
96 
97 //
98 // member functions
99 //
100 
101 // ------------ method called to for each event ------------
103  // get BS
104 
106  iEvent.getByToken(_recoBeamSpotToken, bs);
107 
108  // get PV
109 
111  iEvent.getByToken(_recoVertexCollectionToken, pvcoll);
112 
113  if (_firstOnly) {
114  reco::VertexCollection firstpv;
115  if (!pvcoll->empty())
116  firstpv.push_back((*pvcoll)[0]);
117  _bspvhm.fill(iEvent, firstpv, *bs);
118  } else {
119  _bspvhm.fill(iEvent, *pvcoll, *bs);
120  }
121 }
122 
123 // ------------ method called once each job just before starting event loop ------------
125 
126 void BSvsPVAnalyzer::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) { _bspvhm.beginRun(iRun.run()); }
127 
128 void BSvsPVAnalyzer::endRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {}
129 // ------------ method called once each job just after ending the event loop ------------
131 
132 //define this as a plug-in
static const std::string kSharedResource
Definition: TFileService.h:76
void beginJob() override
void endRun(const edm::Run &, const edm::EventSetup &) override
~BSvsPVAnalyzer() override
void endJob() override
edm::EDGetTokenT< reco::VertexCollection > _recoVertexCollectionToken
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
std::vector< Vertex > VertexCollection
Definition: Vertex.h:31
int iEvent
Definition: GenABIO.cc:224
RunNumber_t run() const
Definition: RunBase.h:40
void beginRun(const edm::Run &, const edm::EventSetup &) override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::BeamSpot > _recoBeamSpotToken
void book(const std::string dirname="")
BSvsPVHistogramMaker _bspvhm
fixed size matrix
HLT enums.
BSvsPVAnalyzer(const edm::ParameterSet &)
void beginRun(const unsigned int nrun)
void analyze(const edm::Event &, const edm::EventSetup &) override
void fill(const unsigned int orbit, const int bx, const reco::VertexCollection &vertices, const reco::BeamSpot &bs)
Definition: Run.h:45