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 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
24 
25 #include <vector>
26 #include <map>
27 #include <limits>
28 #include <string>
29 
36 
38 
40 
44 
46 
47 
48 //
49 // class decleration
50 //
51 
53  public:
56 
57 
58 private:
59  virtual void beginJob() ;
60  virtual void analyze(const edm::Event&, const edm::EventSetup&);
61  virtual void beginRun(const edm::Run&, const edm::EventSetup&);
62  virtual void endRun(const edm::Run&, const edm::EventSetup&);
63  virtual void endJob() ;
64 
65  // ----------member data ---------------------------
66 
69  bool _firstOnly;
70 
72 };
73 
74 //
75 // constants, enums and typedefs
76 //
77 
78 //
79 // static data member definitions
80 //
81 
82 //
83 // constructors and destructor
84 //
86  : _vhm(iConfig.getParameter<edm::ParameterSet>("vHistogramMakerPSet"), consumesCollector())
87  , _recoVertexCollectionToken(consumes<reco::VertexCollection>(iConfig.getParameter<edm::InputTag>("pvCollection")))
88  , _firstOnly(iConfig.getUntrackedParameter<bool>("firstOnly",false))
89  , _weightprov(iConfig.getParameter<bool>("usePrescaleWeight")
90  ? new PrescaleWeightProvider(iConfig.getParameter<edm::ParameterSet>("prescaleWeightProviderPSet"), consumesCollector())
91  : 0
92  )
93 {
94  //now do what ever initialization is needed
95 
96  //
97 
98  _vhm.book();
99 
100 }
101 
102 
104 {
105 
106  // do anything here that needs to be done at desctruction time
107  // (e.g. close files, deallocate resources etc.)
108 
109  delete _weightprov;
110 
111 }
112 
113 
114 //
115 // member functions
116 //
117 
118 // ------------ method called to for each event ------------
119 void
121 {
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.getByToken(_recoVertexCollectionToken,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 &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#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.)
edm::EDGetTokenT< reco::VertexCollection > _recoVertexCollectionToken
int iEvent
Definition: GenABIO.cc:230
void initRun(const edm::Run &run, const edm::EventSetup &setup)
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 &)
volatile std::atomic< bool > shutdown_flag false
int weight
Definition: histoStyle.py:50
virtual void beginRun(const edm::Run &, const edm::EventSetup &)
Definition: Run.h:41