CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FFTJetVertexAdder.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: FFTJetProducers
4 // Class: FFTJetVertexAdder
5 //
13 //
14 // Original Author: Igor Volobouev
15 // Created: Thu Jun 21 19:19:40 CDT 2012
16 //
17 //
18 
19 #include <iostream>
20 #include "CLHEP/Random/RandGauss.h"
21 
22 // framework include files
31 
34 
37 
39 
40 #define init_param(type, varname) varname(ps.getParameter<type>(#varname))
41 
42 //
43 // class declaration
44 //
46 public:
47  explicit FFTJetVertexAdder(const edm::ParameterSet&);
48  FFTJetVertexAdder() = delete;
49  FFTJetVertexAdder(const FFTJetVertexAdder&) = delete;
51  ~FFTJetVertexAdder() override;
52 
53 protected:
54  // methods
55  void beginJob() override;
56  void produce(edm::Event&, const edm::EventSetup&) override;
57  void endJob() override;
58 
59 private:
62 
65 
67 
68  const bool useBeamSpot;
69  const bool addExistingVertices;
70 
71  const double fixedX;
72  const double fixedY;
73  const double fixedZ;
74 
75  const double sigmaX;
76  const double sigmaY;
77  const double sigmaZ;
78 
79  const double nDof;
80  const double chi2;
81  const double errX;
82  const double errY;
83  const double errZ;
84 
85  const unsigned nVerticesToMake;
86 };
87 
88 //
89 // constructors and destructor
90 //
93  init_param(edm::InputTag, existingVerticesLabel),
94  init_param(std::string, outputLabel),
95  init_param(bool, useBeamSpot),
96  init_param(bool, addExistingVertices),
97  init_param(double, fixedX),
98  init_param(double, fixedY),
99  init_param(double, fixedZ),
100  init_param(double, sigmaX),
101  init_param(double, sigmaY),
102  init_param(double, sigmaZ),
103  init_param(double, nDof),
104  init_param(double, chi2),
105  init_param(double, errX),
106  init_param(double, errY),
107  init_param(double, errZ),
108  init_param(unsigned, nVerticesToMake) {
109  if (useBeamSpot)
110  beamSpotToken = consumes<reco::BeamSpot>(beamSpotLabel);
112  existingVerticesToken = consumes<reco::VertexCollection>(existingVerticesLabel);
113  produces<reco::VertexCollection>(outputLabel);
114 }
115 
117 
118 // ------------ method called to produce the data ------------
121  CLHEP::RandGauss rGauss(rng->getEngine(iEvent.streamID()));
122 
123  // get PFCandidates
124  auto pOutput = std::make_unique<reco::VertexCollection>();
125 
126  double xmean = fixedX;
127  double ymean = fixedY;
128  double zmean = fixedZ;
129 
130  double xwidth = sigmaX;
131  double ywidth = sigmaY;
132  double zwidth = sigmaZ;
133 
134  if (useBeamSpot) {
135  edm::Handle<reco::BeamSpot> beamSpotHandle;
136  iEvent.getByToken(beamSpotToken, beamSpotHandle);
137  if (!beamSpotHandle.isValid())
138  throw cms::Exception("FFTJetBadConfig") << "ERROR in FFTJetVertexAdder:"
139  " could not find beam spot information"
140  << std::endl;
141 
142  xmean = beamSpotHandle->x0();
143  ymean = beamSpotHandle->y0();
144  zmean = beamSpotHandle->z0();
145 
146  xwidth = beamSpotHandle->BeamWidthX();
147  ywidth = beamSpotHandle->BeamWidthY();
148  zwidth = beamSpotHandle->sigmaZ();
149  }
150 
152  for (unsigned i = 0; i < 3; ++i)
153  for (unsigned j = 0; j < 3; ++j)
154  err[i][j] = 0.0;
155  err[0][0] = errX * errX;
156  err[1][1] = errY * errY;
157  err[2][2] = errZ * errZ;
158 
159  for (unsigned iv = 0; iv < nVerticesToMake; ++iv) {
160  const double x0 = rGauss(xmean, xwidth);
161  const double y0 = rGauss(ymean, ywidth);
162  const double z0 = rGauss(zmean, zwidth);
163  const reco::Vertex::Point position(x0, y0, z0);
164  pOutput->push_back(reco::Vertex(position, err, chi2, nDof, 0));
165  }
166 
167  if (addExistingVertices) {
168  typedef reco::VertexCollection::const_iterator IV;
169 
171  iEvent.getByToken(existingVerticesToken, vertices);
172  if (!vertices.isValid())
173  throw cms::Exception("FFTJetBadConfig") << "ERROR in FFTJetVertexAdder:"
174  " could not find existing collection of vertices"
175  << std::endl;
176 
177  const IV vertend(vertices->end());
178  for (IV iv = vertices->begin(); iv != vertend; ++iv)
179  pOutput->push_back(*iv);
180  }
181 
182  iEvent.put(std::move(pOutput), outputLabel);
183 }
184 
185 // ------------ method called once each job just before starting event loop
188  if (!rng.isAvailable()) {
189  throw cms::Exception("FFTJetBadConfig") << "ERROR in FFTJetVertexAdder:"
190  " failed to initialize the random number generator"
191  << std::endl;
192  }
193 }
194 
195 // ------------ method called once each job just after ending the event loop
197 
198 //define this as a plug-in
const std::string outputLabel
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
int32_t *__restrict__ iv
const bool addExistingVertices
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
#define init_param(type, varname)
virtual CLHEP::HepRandomEngine & getEngine(StreamID const &)=0
Use this engine in event methods.
edm::EDGetTokenT< reco::VertexCollection > existingVerticesToken
const unsigned nVerticesToMake
const edm::InputTag beamSpotLabel
int iEvent
Definition: GenABIO.cc:224
def move
Definition: eostools.py:511
FFTJetVertexAdder & operator=(const FFTJetVertexAdder &)=delete
bool isAvailable() const
Definition: Service.h:40
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
const edm::InputTag existingVerticesLabel
bool isValid() const
Definition: HandleBase.h:70
~FFTJetVertexAdder() override
void produce(edm::Event &, const edm::EventSetup &) override
void beginJob() override
static int position[264][3]
Definition: ReadPGInfo.cc:289
StreamID streamID() const
Definition: Event.h:98
FFTJetVertexAdder()=delete
void endJob() override
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken