Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:23

0001 // -*- C++ -*-
0002 //
0003 // Package:    FFTJetProducers
0004 // Class:      FFTJetVertexAdder
0005 //
0006 /**\class FFTJetVertexAdder FFTJetVertexAdder.cc RecoJets/FFTJetProducers/plugins/FFTJetVertexAdder.cc
0007 
0008  Description: adds a collection of fake vertices to the event record
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  Igor Volobouev
0015 //         Created:  Thu Jun 21 19:19:40 CDT 2012
0016 //
0017 //
0018 
0019 #include <iostream>
0020 #include "CLHEP/Random/RandGauss.h"
0021 
0022 // framework include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDProducer.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/MakerMacros.h"
0027 #include "FWCore/ServiceRegistry/interface/Service.h"
0028 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
0029 #include "FWCore/Utilities/interface/Exception.h"
0030 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0031 
0032 #include "DataFormats/Common/interface/View.h"
0033 #include "DataFormats/Common/interface/Handle.h"
0034 
0035 #include "DataFormats/VertexReco/interface/Vertex.h"
0036 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0037 
0038 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0039 
0040 #define init_param(type, varname) varname(ps.getParameter<type>(#varname))
0041 
0042 //
0043 // class declaration
0044 //
0045 class FFTJetVertexAdder : public edm::stream::EDProducer<> {
0046 public:
0047   explicit FFTJetVertexAdder(const edm::ParameterSet&);
0048   FFTJetVertexAdder() = delete;
0049   FFTJetVertexAdder(const FFTJetVertexAdder&) = delete;
0050   FFTJetVertexAdder& operator=(const FFTJetVertexAdder&) = delete;
0051   ~FFTJetVertexAdder() override;
0052 
0053 protected:
0054   // methods
0055   void beginStream(edm::StreamID) override;
0056   void produce(edm::Event&, const edm::EventSetup&) override;
0057 
0058 private:
0059   const edm::InputTag beamSpotLabel;
0060   const edm::InputTag existingVerticesLabel;
0061 
0062   edm::EDGetTokenT<reco::BeamSpot> beamSpotToken;
0063   edm::EDGetTokenT<reco::VertexCollection> existingVerticesToken;
0064 
0065   const std::string outputLabel;
0066 
0067   const bool useBeamSpot;
0068   const bool addExistingVertices;
0069 
0070   const double fixedX;
0071   const double fixedY;
0072   const double fixedZ;
0073 
0074   const double sigmaX;
0075   const double sigmaY;
0076   const double sigmaZ;
0077 
0078   const double nDof;
0079   const double chi2;
0080   const double errX;
0081   const double errY;
0082   const double errZ;
0083 
0084   const unsigned nVerticesToMake;
0085 };
0086 
0087 //
0088 // constructors and destructor
0089 //
0090 FFTJetVertexAdder::FFTJetVertexAdder(const edm::ParameterSet& ps)
0091     : init_param(edm::InputTag, beamSpotLabel),
0092       init_param(edm::InputTag, existingVerticesLabel),
0093       init_param(std::string, outputLabel),
0094       init_param(bool, useBeamSpot),
0095       init_param(bool, addExistingVertices),
0096       init_param(double, fixedX),
0097       init_param(double, fixedY),
0098       init_param(double, fixedZ),
0099       init_param(double, sigmaX),
0100       init_param(double, sigmaY),
0101       init_param(double, sigmaZ),
0102       init_param(double, nDof),
0103       init_param(double, chi2),
0104       init_param(double, errX),
0105       init_param(double, errY),
0106       init_param(double, errZ),
0107       init_param(unsigned, nVerticesToMake) {
0108   if (useBeamSpot)
0109     beamSpotToken = consumes<reco::BeamSpot>(beamSpotLabel);
0110   if (addExistingVertices)
0111     existingVerticesToken = consumes<reco::VertexCollection>(existingVerticesLabel);
0112   produces<reco::VertexCollection>(outputLabel);
0113 }
0114 
0115 FFTJetVertexAdder::~FFTJetVertexAdder() {}
0116 
0117 // ------------ method called to produce the data  ------------
0118 void FFTJetVertexAdder::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0119   edm::Service<edm::RandomNumberGenerator> rng;
0120   CLHEP::RandGauss rGauss(rng->getEngine(iEvent.streamID()));
0121 
0122   // get PFCandidates
0123   auto pOutput = std::make_unique<reco::VertexCollection>();
0124 
0125   double xmean = fixedX;
0126   double ymean = fixedY;
0127   double zmean = fixedZ;
0128 
0129   double xwidth = sigmaX;
0130   double ywidth = sigmaY;
0131   double zwidth = sigmaZ;
0132 
0133   if (useBeamSpot) {
0134     edm::Handle<reco::BeamSpot> beamSpotHandle;
0135     iEvent.getByToken(beamSpotToken, beamSpotHandle);
0136     if (!beamSpotHandle.isValid())
0137       throw cms::Exception("FFTJetBadConfig") << "ERROR in FFTJetVertexAdder:"
0138                                                  " could not find beam spot information"
0139                                               << std::endl;
0140 
0141     xmean = beamSpotHandle->x0();
0142     ymean = beamSpotHandle->y0();
0143     zmean = beamSpotHandle->z0();
0144 
0145     xwidth = beamSpotHandle->BeamWidthX();
0146     ywidth = beamSpotHandle->BeamWidthY();
0147     zwidth = beamSpotHandle->sigmaZ();
0148   }
0149 
0150   reco::Vertex::Error err;
0151   for (unsigned i = 0; i < 3; ++i)
0152     for (unsigned j = 0; j < 3; ++j)
0153       err[i][j] = 0.0;
0154   err[0][0] = errX * errX;
0155   err[1][1] = errY * errY;
0156   err[2][2] = errZ * errZ;
0157 
0158   for (unsigned iv = 0; iv < nVerticesToMake; ++iv) {
0159     const double x0 = rGauss(xmean, xwidth);
0160     const double y0 = rGauss(ymean, ywidth);
0161     const double z0 = rGauss(zmean, zwidth);
0162     const reco::Vertex::Point position(x0, y0, z0);
0163     pOutput->push_back(reco::Vertex(position, err, chi2, nDof, 0));
0164   }
0165 
0166   if (addExistingVertices) {
0167     typedef reco::VertexCollection::const_iterator IV;
0168 
0169     edm::Handle<reco::VertexCollection> vertices;
0170     iEvent.getByToken(existingVerticesToken, vertices);
0171     if (!vertices.isValid())
0172       throw cms::Exception("FFTJetBadConfig") << "ERROR in FFTJetVertexAdder:"
0173                                                  " could not find existing collection of vertices"
0174                                               << std::endl;
0175 
0176     const IV vertend(vertices->end());
0177     for (IV iv = vertices->begin(); iv != vertend; ++iv)
0178       pOutput->push_back(*iv);
0179   }
0180 
0181   iEvent.put(std::move(pOutput), outputLabel);
0182 }
0183 
0184 // ------------ method called once each job just before starting event loop
0185 void FFTJetVertexAdder::beginStream(edm::StreamID) {
0186   edm::Service<edm::RandomNumberGenerator> rng;
0187   if (!rng.isAvailable()) {
0188     throw cms::Exception("FFTJetBadConfig") << "ERROR in FFTJetVertexAdder:"
0189                                                " failed to initialize the random number generator"
0190                                             << std::endl;
0191   }
0192 }
0193 
0194 //define this as a plug-in
0195 DEFINE_FWK_MODULE(FFTJetVertexAdder);