CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AnotherPrimaryVertexAnalyzer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: Validation/RecoVertex
4 // Class: AnotherPrimaryVertexAnalyzer
5 //
13 //
14 // Original Author: Andrea Venturi
15 // Created: Mon Oct 27 17:37:53 CET 2008
16 // $Id: AnotherPrimaryVertexAnalyzer.cc,v 1.2 2011/12/11 10:51:49 venturia Exp $
17 //
18 //
19 
20 
21 // system include files
22 #include <memory>
23 
24 // user include files
25 
26 #include <vector>
27 #include <map>
28 #include <limits>
29 #include <string>
30 
33 
37 
39 
41 
45 
47 
48 
49 //
50 // class decleration
51 //
52 
54  public:
57 
58 
59 private:
60  virtual void beginJob() ;
61  virtual void analyze(const edm::Event&, const edm::EventSetup&);
62  virtual void beginRun(const edm::Run&, const edm::EventSetup&);
63  virtual void endRun(const edm::Run&, const edm::EventSetup&);
64  virtual void endJob() ;
65 
66  // ----------member data ---------------------------
67 
70  bool _firstOnly;
71 
73 };
74 
75 //
76 // constants, enums and typedefs
77 //
78 
79 //
80 // static data member definitions
81 //
82 
83 //
84 // constructors and destructor
85 //
87  _vhm(iConfig.getParameter<edm::ParameterSet>("vHistogramMakerPSet")),
88  _pvcollection(iConfig.getParameter<edm::InputTag>("pvCollection")),
89  _firstOnly(iConfig.getUntrackedParameter<bool>("firstOnly",false)),
90  _weightprov(iConfig.getParameter<bool>("usePrescaleWeight") ?
91  new PrescaleWeightProvider(iConfig.getParameter<edm::ParameterSet>("prescaleWeightProviderPSet")) : 0)
92 {
93  //now do what ever initialization is needed
94 
95  //
96 
97  _vhm.book();
98 
99 }
100 
101 
103 {
104 
105  // do anything here that needs to be done at desctruction time
106  // (e.g. close files, deallocate resources etc.)
107 
108  delete _weightprov;
109 
110 }
111 
112 
113 //
114 // member functions
115 //
116 
117 // ------------ method called to for each event ------------
118 void
120 {
121  using namespace edm;
122 
123  // compute event weigth
124 
125  double weight = 1.;
126 
127  if(_weightprov) weight = _weightprov->prescaleWeight(iEvent,iSetup);
128 
129  // get PV
130 
132  iEvent.getByLabel(_pvcollection,pvcoll);
133 
134  if(_firstOnly) {
135  reco::VertexCollection firstpv;
136  if(pvcoll->size()) firstpv.push_back((*pvcoll)[0]);
137  _vhm.fill(iEvent,firstpv,weight);
138  }
139  else {
140  _vhm.fill(iEvent,*pvcoll,weight);
141  }
142 }
143 
144 
145 // ------------ method called once each job just before starting event loop ------------
146 void
148 { }
149 
150 void
152 
153  _vhm.beginRun(iRun);
154 
155  if(_weightprov) _weightprov->initRun(iRun,iSetup);
156 
157 }
158 
159 void
161 
162 }
163 // ------------ method called once each job just after ending the event loop ------------
164 void
166 }
167 
168 
169 //define this as a plug-in
virtual void analyze(const edm::Event &, const edm::EventSetup &)
int prescaleWeight(const edm::Event &event, const edm::EventSetup &setup)
AnotherPrimaryVertexAnalyzer(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
void fill(const edm::Event &iEvent, const reco::VertexCollection &vertices, const double weight=1.)
int iEvent
Definition: GenABIO.cc:243
void initRun(const edm::Run &run, const edm::EventSetup &setup)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
void book(const std::string dirname="")
void beginRun(const edm::Run &iRun)
This class takes a vector of HLT paths and returns a weight based on their HLT and L1 prescales...
virtual void endRun(const edm::Run &, const edm::EventSetup &)
int weight
Definition: histoStyle.py:50
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
Definition: Run.h:36