-
Notifications
You must be signed in to change notification settings - Fork 4.2k
/
FFTJetVertexAdder.cc
212 lines (172 loc) · 5.83 KB
/
FFTJetVertexAdder.cc
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
// -*- C++ -*-
//
// Package: FFTJetProducers
// Class: FFTJetVertexAdder
//
/**\class FFTJetVertexAdder FFTJetVertexAdder.cc RecoJets/FFTJetProducers/plugins/FFTJetVertexAdder.cc
Description: adds a collection of fake vertices to the event record
Implementation:
[Notes on implementation]
*/
//
// Original Author: Igor Volobouev
// Created: Thu Jun 21 19:19:40 CDT 2012
//
//
#include <iostream>
#include "CLHEP/Random/RandGauss.h"
// framework include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ServiceRegistry/interface/Service.h"
#include "FWCore/Utilities/interface/RandomNumberGenerator.h"
#include "FWCore/Utilities/interface/Exception.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/Common/interface/View.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/VertexReco/interface/Vertex.h"
#include "DataFormats/VertexReco/interface/VertexFwd.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
//
// class declaration
//
class FFTJetVertexAdder : public edm::EDProducer
{
public:
explicit FFTJetVertexAdder(const edm::ParameterSet&);
~FFTJetVertexAdder();
protected:
// methods
void beginJob() override;
void produce(edm::Event&, const edm::EventSetup&) override;
void endJob() override;
private:
FFTJetVertexAdder();
FFTJetVertexAdder(const FFTJetVertexAdder&);
FFTJetVertexAdder& operator=(const FFTJetVertexAdder&);
const edm::InputTag beamSpotLabel;
const edm::InputTag existingVerticesLabel;
const std::string outputLabel;
const bool useBeamSpot;
const bool addExistingVertices;
const double fixedX;
const double fixedY;
const double fixedZ;
const double sigmaX;
const double sigmaY;
const double sigmaZ;
const double nDof;
const double chi2;
const double errX;
const double errY;
const double errZ;
const unsigned nVerticesToMake;
};
//
// constructors and destructor
//
FFTJetVertexAdder::FFTJetVertexAdder(const edm::ParameterSet& ps)
: init_param(edm::InputTag, beamSpotLabel),
init_param(edm::InputTag, existingVerticesLabel),
init_param(std::string, outputLabel),
init_param(bool, useBeamSpot),
init_param(bool, addExistingVertices),
init_param(double, fixedX),
init_param(double, fixedY),
init_param(double, fixedZ),
init_param(double, sigmaX),
init_param(double, sigmaY),
init_param(double, sigmaZ),
init_param(double, nDof),
init_param(double, chi2),
init_param(double, errX),
init_param(double, errY),
init_param(double, errZ),
init_param(unsigned, nVerticesToMake)
{
produces<reco::VertexCollection>(outputLabel);
}
FFTJetVertexAdder::~FFTJetVertexAdder()
{
}
// ------------ method called to produce the data ------------
void FFTJetVertexAdder::produce(
edm::Event& iEvent, const edm::EventSetup& iSetup)
{
edm::Service<edm::RandomNumberGenerator> rng;
CLHEP::RandGauss rGauss(rng->getEngine(iEvent.streamID()));
// get PFCandidates
std::auto_ptr<reco::VertexCollection> pOutput(new reco::VertexCollection);
double xmean = fixedX;
double ymean = fixedY;
double zmean = fixedZ;
double xwidth = sigmaX;
double ywidth = sigmaY;
double zwidth = sigmaZ;
if (useBeamSpot)
{
edm::Handle<reco::BeamSpot> beamSpotHandle;
iEvent.getByLabel(beamSpotLabel, beamSpotHandle);
if (!beamSpotHandle.isValid())
throw cms::Exception("FFTJetBadConfig")
<< "ERROR in FFTJetVertexAdder:"
" could not find beam spot information"
<< std::endl;
xmean = beamSpotHandle->x0();
ymean = beamSpotHandle->y0();
zmean = beamSpotHandle->z0();
xwidth = beamSpotHandle->BeamWidthX();
ywidth = beamSpotHandle->BeamWidthY();
zwidth = beamSpotHandle->sigmaZ();
}
reco::Vertex::Error err;
for (unsigned i=0; i<3; ++i)
for (unsigned j=0; j<3; ++j)
err[i][j] = 0.0;
err[0][0] = errX*errX;
err[1][1] = errY*errY;
err[2][2] = errZ*errZ;
for (unsigned iv=0; iv<nVerticesToMake; ++iv)
{
const double x0 = rGauss(xmean, xwidth);
const double y0 = rGauss(ymean, ywidth);
const double z0 = rGauss(zmean, zwidth);
const reco::Vertex::Point position(x0, y0, z0);
pOutput->push_back(reco::Vertex(position, err, chi2, nDof, 0));
}
if (addExistingVertices)
{
typedef reco::VertexCollection::const_iterator IV;
edm::Handle<reco::VertexCollection> vertices;
iEvent.getByLabel(existingVerticesLabel, vertices);
if (!vertices.isValid())
throw cms::Exception("FFTJetBadConfig")
<< "ERROR in FFTJetVertexAdder:"
" could not find existing collection of vertices"
<< std::endl;
const IV vertend(vertices->end());
for (IV iv=vertices->begin(); iv!=vertend; ++iv)
pOutput->push_back(*iv);
}
iEvent.put(pOutput, outputLabel);
}
// ------------ method called once each job just before starting event loop
void FFTJetVertexAdder::beginJob()
{
edm::Service<edm::RandomNumberGenerator> rng;
if ( !rng.isAvailable() ) {
throw cms::Exception("FFTJetBadConfig")
<< "ERROR in FFTJetVertexAdder:"
" failed to initialize the random number generator"
<< std::endl;
}
}
// ------------ method called once each job just after ending the event loop
void FFTJetVertexAdder::endJob()
{
}
//define this as a plug-in
DEFINE_FWK_MODULE(FFTJetVertexAdder);