Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 10:05:17

0001 // -*- C++ -*-
0002 //
0003 // Package:    TauAnalysis/EmbeddingProducer
0004 // Class:      EmbeddingVertexCorrector
0005 //
0006 /**\class EmbeddingVertexCorrector EmbeddingVertexCorrector.cc TauAnalysis/EmbeddingProducer/plugins/EmbeddingVertexCorrector.cc
0007 
0008  Description: [one line class summary]
0009 
0010  Implementation:
0011      [Notes on implementation]
0012 */
0013 //
0014 // Original Author:  Artur Akhmetshin
0015 //         Created:  Sat, 23 Apr 2016 21:47:13 GMT
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 
0022 // user include files
0023 #include "FWCore/Framework/interface/Frameworkfwd.h"
0024 #include "FWCore/Framework/interface/stream/EDProducer.h"
0025 
0026 #include "FWCore/Framework/interface/Event.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 
0029 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0030 #include "FWCore/Utilities/interface/StreamID.h"
0031 
0032 #include "CLHEP/Units/GlobalSystemOfUnits.h"
0033 #include "CLHEP/Units/GlobalPhysicalConstants.h"
0034 
0035 #include "DataFormats/Math/interface/LorentzVector.h"
0036 #include "DataFormats/Math/interface/LorentzVectorFwd.h"
0037 #include "DataFormats/VertexReco/interface/Vertex.h"
0038 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0039 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0040 
0041 #include "CondFormats/DataRecord/interface/SimBeamSpotObjectsRcd.h"
0042 #include "CondFormats/BeamSpotObjects/interface/SimBeamSpotObjects.h"
0043 #include "CondFormats/DataRecord/interface/BeamSpotObjectsRcd.h"
0044 #include "CondFormats/BeamSpotObjects/interface/BeamSpotObjects.h"
0045 
0046 #include "FWCore/Framework/interface/EventSetup.h"
0047 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0048 
0049 namespace HepMC {
0050   class FourVector;
0051 }
0052 
0053 //
0054 // class declaration
0055 //
0056 
0057 class EmbeddingVertexCorrector : public edm::stream::EDProducer<> {
0058 public:
0059   explicit EmbeddingVertexCorrector(const edm::ParameterSet&);
0060   ~EmbeddingVertexCorrector() override;
0061 
0062   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0063 
0064 private:
0065   void produce(edm::Event&, const edm::EventSetup&) override;
0066 
0067   // ----------member data ---------------------------
0068   edm::InputTag sourceLabel;
0069   edm::InputTag vertexPositionLabel;
0070 };
0071 
0072 //
0073 // constructors and destructor
0074 //
0075 EmbeddingVertexCorrector::EmbeddingVertexCorrector(const edm::ParameterSet& iConfig) {
0076   produces<edm::HepMCProduct>();
0077 
0078   sourceLabel = iConfig.getParameter<edm::InputTag>("src");
0079   consumes<edm::HepMCProduct>(sourceLabel);
0080   vertexPositionLabel = edm::InputTag("externalLHEProducer", "vertexPosition");
0081   consumes<math::XYZTLorentzVectorD>(vertexPositionLabel);
0082 }
0083 
0084 EmbeddingVertexCorrector::~EmbeddingVertexCorrector() {}
0085 
0086 //
0087 // member functions
0088 //
0089 
0090 // ------------ method called to produce the data  ------------
0091 void EmbeddingVertexCorrector::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0092   using namespace edm;
0093 
0094   // Retrieving generated Z to TauTau Event
0095   Handle<edm::HepMCProduct> InputGenEvent;
0096   iEvent.getByLabel(sourceLabel, InputGenEvent);
0097   HepMC::GenEvent* genevent = new HepMC::GenEvent(*InputGenEvent->GetEvent());
0098   std::unique_ptr<edm::HepMCProduct> CorrectedGenEvent(new edm::HepMCProduct(genevent));
0099 
0100   //Retrieving vertex position from input and creating vertex shift
0101   Handle<math::XYZTLorentzVectorD> vertex_position;
0102   iEvent.getByLabel(vertexPositionLabel, vertex_position);
0103   HepMC::FourVector vertex_shift(
0104       vertex_position.product()->x() * cm, vertex_position.product()->y() * cm, vertex_position.product()->z() * cm);
0105 
0106   // Apply vertex shift to all production vertices of the event
0107   CorrectedGenEvent->applyVtxGen(&vertex_shift);
0108   iEvent.put(std::move(CorrectedGenEvent));
0109 }
0110 
0111 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
0112 void EmbeddingVertexCorrector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0113   //The following says we do not know what parameters are allowed so do no validation
0114   // Please change this to state exactly what you do use, even if it is no parameters
0115   edm::ParameterSetDescription desc;
0116   desc.setUnknown();
0117   descriptions.addDefault(desc);
0118 }
0119 
0120 //define this as a plug-in
0121 DEFINE_FWK_MODULE(EmbeddingVertexCorrector);