CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 // $Id: FFTJetVertexAdder.cc,v 1.1 2012/06/22 07:23:21 igv Exp $
17 //
18 //
19 
20 #include <iostream>
21 #include "CLHEP/Random/RandGauss.h"
22 
23 // framework include files
32 
35 
38 
40 
41 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
42 
43 //
44 // class declaration
45 //
47 {
48 public:
49  explicit FFTJetVertexAdder(const edm::ParameterSet&);
51 
52 protected:
53  // methods
54  void beginJob();
55  void produce(edm::Event&, const edm::EventSetup&);
56  void endJob();
57 
58 private:
62 
66 
67  const bool useBeamSpot;
68  const bool addExistingVertices;
69 
70  const double fixedX;
71  const double fixedY;
72  const double fixedZ;
73 
74  const double sigmaX;
75  const double sigmaY;
76  const double sigmaZ;
77 
78  const double nDof;
79  const double chi2;
80  const double errX;
81  const double errY;
82  const double errZ;
83 
84  const unsigned nVerticesToMake;
85 
86  CLHEP::RandGauss* rGauss_;
87 };
88 
89 //
90 // constructors and destructor
91 //
93  : init_param(edm::InputTag, beamSpotLabel),
94  init_param(edm::InputTag, existingVerticesLabel),
95  init_param(std::string, outputLabel),
96  init_param(bool, useBeamSpot),
97  init_param(bool, addExistingVertices),
98  init_param(double, fixedX),
99  init_param(double, fixedY),
100  init_param(double, fixedZ),
101  init_param(double, sigmaX),
102  init_param(double, sigmaY),
103  init_param(double, sigmaZ),
104  init_param(double, nDof),
105  init_param(double, chi2),
106  init_param(double, errX),
107  init_param(double, errY),
108  init_param(double, errZ),
109  init_param(unsigned, nVerticesToMake),
110  rGauss_(0)
111 {
112  produces<reco::VertexCollection>(outputLabel);
113 }
114 
115 
117 {
118  delete rGauss_;
119 }
120 
121 
122 // ------------ method called to produce the data ------------
124  edm::Event& iEvent, const edm::EventSetup& iSetup)
125 {
126  // get PFCandidates
127  std::auto_ptr<reco::VertexCollection> pOutput(new reco::VertexCollection);
128 
129  double xmean = fixedX;
130  double ymean = fixedY;
131  double zmean = fixedZ;
132 
133  double xwidth = sigmaX;
134  double ywidth = sigmaY;
135  double zwidth = sigmaZ;
136 
137  if (useBeamSpot)
138  {
139  edm::Handle<reco::BeamSpot> beamSpotHandle;
140  iEvent.getByLabel(beamSpotLabel, beamSpotHandle);
141  if (!beamSpotHandle.isValid())
142  throw cms::Exception("FFTJetBadConfig")
143  << "ERROR in FFTJetVertexAdder:"
144  " could not find beam spot information"
145  << std::endl;
146 
147  xmean = beamSpotHandle->x0();
148  ymean = beamSpotHandle->y0();
149  zmean = beamSpotHandle->z0();
150 
151  xwidth = beamSpotHandle->BeamWidthX();
152  ywidth = beamSpotHandle->BeamWidthY();
153  zwidth = beamSpotHandle->sigmaZ();
154  }
155 
157  for (unsigned i=0; i<3; ++i)
158  for (unsigned j=0; j<3; ++j)
159  err[i][j] = 0.0;
160  err[0][0] = errX*errX;
161  err[1][1] = errY*errY;
162  err[2][2] = errZ*errZ;
163 
164  for (unsigned iv=0; iv<nVerticesToMake; ++iv)
165  {
166  const double x0 = (*rGauss_)(xmean, xwidth);
167  const double y0 = (*rGauss_)(ymean, ywidth);
168  const double z0 = (*rGauss_)(zmean, zwidth);
169  const reco::Vertex::Point position(x0, y0, z0);
170  pOutput->push_back(reco::Vertex(position, err, chi2, nDof, 0));
171  }
172 
174  {
175  typedef reco::VertexCollection::const_iterator IV;
176 
178  iEvent.getByLabel(existingVerticesLabel, vertices);
179  if (!vertices.isValid())
180  throw cms::Exception("FFTJetBadConfig")
181  << "ERROR in FFTJetVertexAdder:"
182  " could not find existing collection of vertices"
183  << std::endl;
184 
185  const IV vertend(vertices->end());
186  for (IV iv=vertices->begin(); iv!=vertend; ++iv)
187  pOutput->push_back(*iv);
188  }
189 
190  iEvent.put(pOutput, outputLabel);
191 }
192 
193 
194 // ------------ method called once each job just before starting event loop
196 {
197  if (!rGauss_)
198  {
200  if ( !rng.isAvailable() )
201  throw cms::Exception("FFTJetBadConfig")
202  << "ERROR in FFTJetVertexAdder:"
203  " failed to initialize the random number generator"
204  << std::endl;
205  rGauss_ = new CLHEP::RandGauss(rng->getEngine());
206  }
207 }
208 
209 
210 // ------------ method called once each job just after ending the event loop
212 {
213 }
214 
215 
216 //define this as a plug-in
int i
Definition: DBlmapReader.cc:9
const std::string outputLabel
const bool addExistingVertices
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
std::vector< Vertex > VertexCollection
collection of Vertex objects
Definition: VertexFwd.h:9
#define init_param(type, varname)
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
const unsigned nVerticesToMake
const edm::InputTag beamSpotLabel
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
bool isAvailable() const
Definition: Service.h:47
int j
Definition: DBlmapReader.cc:9
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
const edm::InputTag existingVerticesLabel
virtual CLHEP::HepRandomEngine & getEngine() const =0
Use this to get the random number engine, this is the only function most users should call...
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
FFTJetVertexAdder & operator=(const FFTJetVertexAdder &)
CLHEP::RandGauss * rGauss_
void produce(edm::Event &, const edm::EventSetup &)