CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2_patch1/src/RecoJets/FFTJetProducers/plugins/FFTJetVertexAdder.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    FFTJetProducers
00004 // Class:      FFTJetVertexAdder
00005 // 
00013 //
00014 // Original Author:  Igor Volobouev
00015 //         Created:  Thu Jun 21 19:19:40 CDT 2012
00016 // $Id: FFTJetVertexAdder.cc,v 1.1 2012/06/22 07:23:21 igv Exp $
00017 //
00018 //
00019 
00020 #include <iostream>
00021 #include "CLHEP/Random/RandGauss.h"
00022 
00023 // framework include files
00024 #include "FWCore/Framework/interface/Frameworkfwd.h"
00025 #include "FWCore/Framework/interface/EDProducer.h"
00026 #include "FWCore/Framework/interface/Event.h"
00027 #include "FWCore/Framework/interface/MakerMacros.h"
00028 #include "FWCore/ServiceRegistry/interface/Service.h"
00029 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00030 #include "FWCore/Utilities/interface/Exception.h"
00031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00032 
00033 #include "DataFormats/Common/interface/View.h"
00034 #include "DataFormats/Common/interface/Handle.h"
00035 
00036 #include "DataFormats/VertexReco/interface/Vertex.h"
00037 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00038 
00039 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00040 
00041 #define init_param(type, varname) varname (ps.getParameter< type >( #varname ))
00042 
00043 //
00044 // class declaration
00045 //
00046 class FFTJetVertexAdder : public edm::EDProducer
00047 {
00048 public:
00049     explicit FFTJetVertexAdder(const edm::ParameterSet&);
00050     ~FFTJetVertexAdder();
00051 
00052 protected:
00053     // methods
00054     void beginJob();
00055     void produce(edm::Event&, const edm::EventSetup&);
00056     void endJob();
00057 
00058 private:
00059     FFTJetVertexAdder();
00060     FFTJetVertexAdder(const FFTJetVertexAdder&);
00061     FFTJetVertexAdder& operator=(const FFTJetVertexAdder&);
00062 
00063     const edm::InputTag beamSpotLabel;
00064     const edm::InputTag existingVerticesLabel;
00065     const std::string outputLabel;
00066 
00067     const bool useBeamSpot;
00068     const bool addExistingVertices;
00069 
00070     const double fixedX;
00071     const double fixedY;
00072     const double fixedZ;
00073 
00074     const double sigmaX;
00075     const double sigmaY;
00076     const double sigmaZ;
00077 
00078     const double nDof;
00079     const double chi2;
00080     const double errX;
00081     const double errY;
00082     const double errZ;
00083 
00084     const unsigned nVerticesToMake;
00085 
00086     CLHEP::RandGauss* rGauss_;
00087 };
00088 
00089 //
00090 // constructors and destructor
00091 //
00092 FFTJetVertexAdder::FFTJetVertexAdder(const edm::ParameterSet& ps)
00093     : init_param(edm::InputTag, beamSpotLabel),
00094       init_param(edm::InputTag, existingVerticesLabel),
00095       init_param(std::string, outputLabel),
00096       init_param(bool, useBeamSpot),
00097       init_param(bool, addExistingVertices),
00098       init_param(double, fixedX),
00099       init_param(double, fixedY),
00100       init_param(double, fixedZ),
00101       init_param(double, sigmaX),
00102       init_param(double, sigmaY),
00103       init_param(double, sigmaZ),
00104       init_param(double, nDof),
00105       init_param(double, chi2),
00106       init_param(double, errX),
00107       init_param(double, errY),
00108       init_param(double, errZ),
00109       init_param(unsigned, nVerticesToMake),
00110       rGauss_(0)
00111 {
00112     produces<reco::VertexCollection>(outputLabel);
00113 }
00114 
00115 
00116 FFTJetVertexAdder::~FFTJetVertexAdder()
00117 {
00118     delete rGauss_;
00119 }
00120 
00121 
00122 // ------------ method called to produce the data  ------------
00123 void FFTJetVertexAdder::produce(
00124     edm::Event& iEvent, const edm::EventSetup& iSetup)
00125 {
00126     // get PFCandidates
00127     std::auto_ptr<reco::VertexCollection> pOutput(new reco::VertexCollection);
00128 
00129     double xmean = fixedX;
00130     double ymean = fixedY;
00131     double zmean = fixedZ;
00132 
00133     double xwidth = sigmaX;
00134     double ywidth = sigmaY;
00135     double zwidth = sigmaZ;
00136 
00137     if (useBeamSpot)
00138     {
00139         edm::Handle<reco::BeamSpot> beamSpotHandle;
00140         iEvent.getByLabel(beamSpotLabel, beamSpotHandle);
00141         if (!beamSpotHandle.isValid())
00142             throw cms::Exception("FFTJetBadConfig")
00143                 << "ERROR in FFTJetVertexAdder:"
00144                 " could not find beam spot information"
00145                 << std::endl;
00146 
00147         xmean = beamSpotHandle->x0();
00148         ymean = beamSpotHandle->y0();
00149         zmean = beamSpotHandle->z0();
00150 
00151         xwidth = beamSpotHandle->BeamWidthX();
00152         ywidth = beamSpotHandle->BeamWidthY();
00153         zwidth = beamSpotHandle->sigmaZ();
00154     }
00155 
00156     reco::Vertex::Error err;
00157     for (unsigned i=0; i<3; ++i)
00158         for (unsigned j=0; j<3; ++j)
00159             err[i][j] = 0.0;
00160     err[0][0] = errX*errX;
00161     err[1][1] = errY*errY;
00162     err[2][2] = errZ*errZ;
00163 
00164     for (unsigned iv=0; iv<nVerticesToMake; ++iv)
00165     {
00166         const double x0 = (*rGauss_)(xmean, xwidth);
00167         const double y0 = (*rGauss_)(ymean, ywidth);
00168         const double z0 = (*rGauss_)(zmean, zwidth);
00169         const reco::Vertex::Point position(x0, y0, z0);
00170         pOutput->push_back(reco::Vertex(position, err, chi2, nDof, 0));
00171     }
00172 
00173     if (addExistingVertices)
00174     {
00175         typedef reco::VertexCollection::const_iterator IV;
00176 
00177         edm::Handle<reco::VertexCollection> vertices;
00178         iEvent.getByLabel(existingVerticesLabel, vertices);
00179         if (!vertices.isValid())
00180             throw cms::Exception("FFTJetBadConfig")
00181                 << "ERROR in FFTJetVertexAdder:"
00182                 " could not find existing collection of vertices"
00183                 << std::endl;
00184 
00185         const IV vertend(vertices->end());
00186         for (IV iv=vertices->begin(); iv!=vertend; ++iv)
00187             pOutput->push_back(*iv);
00188     }
00189 
00190     iEvent.put(pOutput, outputLabel);
00191 }
00192 
00193 
00194 // ------------ method called once each job just before starting event loop
00195 void FFTJetVertexAdder::beginJob()
00196 {
00197     if (!rGauss_)
00198     {
00199         edm::Service<edm::RandomNumberGenerator> rng;
00200         if ( !rng.isAvailable() )
00201             throw cms::Exception("FFTJetBadConfig")
00202                 << "ERROR in FFTJetVertexAdder:"
00203                 " failed to initialize the random number generator"
00204                 << std::endl;
00205         rGauss_ = new CLHEP::RandGauss(rng->getEngine());
00206     }
00207 }
00208 
00209 
00210 // ------------ method called once each job just after ending the event loop
00211 void FFTJetVertexAdder::endJob()
00212 {
00213 }
00214 
00215 
00216 //define this as a plug-in
00217 DEFINE_FWK_MODULE(FFTJetVertexAdder);