LHE2HepMCConverter

Line Code
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157
// -*- C++ -*-
//
// Package:    LHE2HepMCConverter
// Class:      LHE2HepMCConverter
//
/**\class LHE2HepMCConverter LHE2HepMCConverter.cc GeneratorInterface/LHE2HepMCConverter/src/LHE2HepMCConverter.cc

 Description: [one line class summary]

 Implementation:
     [Notes on implementation]
*/
//
// Original Author:  Piergiulio Lenzi,40 1-B01,+41227671638,
//         Created:  Wed Aug 31 19:02:24 CEST 2011
//
//

// system include files
#include <memory>

// user include files
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/one/EDProducer.h"

#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/Run.h"
#include "FWCore/Framework/interface/MakerMacros.h"

#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "FWCore/Utilities/interface/InputTag.h"
#include "DataFormats/Common/interface/Handle.h"

#include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
#include "SimDataFormats/GeneratorProducts/interface/LHERunInfoProduct.h"
#include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"

//
// class declaration
//

class LHE2HepMCConverter : public edm::one::EDProducer<edm::one::WatchRuns> {
public:
  explicit LHE2HepMCConverter(const edm::ParameterSet&);
  ~LHE2HepMCConverter() override = default;

  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);

protected:
  //      lhef::LHERunInfo *lheRunInfo() { return lheRunInfo_.get(); }

private:
  void produce(edm::Event&, const edm::EventSetup&) override;
  void beginRun(edm::Run const&, edm::EventSetup const&) override;
  void endRun(edm::Run const&, edm::EventSetup const&) override {}

  // ----------member data ---------------------------
  const LHERunInfoProduct* _lheRunSrc;
  const edm::InputTag _lheEventSrcTag;
  const edm::InputTag _lheRunSrcTag;
  const edm::EDGetTokenT<LHEEventProduct> eventProductToken_;
  const edm::EDGetTokenT<LHERunInfoProduct> runInfoProductToken_;
  //      std::shared_ptr<lhef::LHERunInfo> lheRunInfo_;
};

//
// constants, enums and typedefs
//

//
// static data member definitions
//

//
// constructors and destructor
//
LHE2HepMCConverter::LHE2HepMCConverter(const edm::ParameterSet& iConfig)
    : _lheRunSrc(nullptr),
      _lheEventSrcTag(iConfig.getParameter<edm::InputTag>("LHEEventProduct")),
      _lheRunSrcTag(iConfig.getParameter<edm::InputTag>("LHERunInfoProduct")),
      eventProductToken_(consumes<LHEEventProduct>(_lheEventSrcTag)),
      runInfoProductToken_(consumes<LHERunInfoProduct, edm::InRun>(_lheRunSrcTag)) {
  //register your products
  produces<edm::HepMCProduct>("unsmeared");
}

//
// member functions
//

// ------------ method called to produce the data  ------------
void LHE2HepMCConverter::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
  using namespace edm;
  const edm::Handle<LHEEventProduct>& lheEventSrc = iEvent.getHandle(eventProductToken_);

  HepMC::GenEvent* evt = new HepMC::GenEvent();
  HepMC::GenVertex* v = new HepMC::GenVertex();
  evt->add_vertex(v);
  if (_lheRunSrc) {
    HepMC::FourVector beam1(0, 0, _lheRunSrc->heprup().EBMUP.first, _lheRunSrc->heprup().EBMUP.first);
    HepMC::GenParticle* gp1 = new HepMC::GenParticle(beam1, _lheRunSrc->heprup().IDBMUP.first, 4);
    v->add_particle_in(gp1);
    HepMC::FourVector beam2(0, 0, _lheRunSrc->heprup().EBMUP.second, _lheRunSrc->heprup().EBMUP.second);
    HepMC::GenParticle* gp2 = new HepMC::GenParticle(beam2, _lheRunSrc->heprup().IDBMUP.second, 4);
    v->add_particle_in(gp2);
    evt->set_beam_particles(gp1, gp2);
  } else {
    LogWarning("LHE2HepMCConverter") << "Could not retrieve the LHERunInfoProduct for this event. You'll miss the beam "
                                        "particles in your HepMC product.";
  }

  for (int i = 0; i < lheEventSrc->hepeup().NUP; ++i) {
    if (lheEventSrc->hepeup().ISTUP[i] != 1) {
      //cout << reader->hepeup.ISTUP[i] << ", " << reader->hepeup.IDUP[i] << endl;
      continue;
    }
    HepMC::FourVector p(lheEventSrc->hepeup().PUP[i][0],
                        lheEventSrc->hepeup().PUP[i][1],
                        lheEventSrc->hepeup().PUP[i][2],
                        lheEventSrc->hepeup().PUP[i][3]);
    HepMC::GenParticle* gp = new HepMC::GenParticle(p, lheEventSrc->hepeup().IDUP[i], 1);
    gp->set_generated_mass(lheEventSrc->hepeup().PUP[i][4]);
    v->add_particle_out(gp);
  }

  std::unique_ptr<HepMCProduct> pOut(new HepMCProduct(evt));
  iEvent.put(std::move(pOut), "unsmeared");
}

// ------------ method called when starting to processes a run  ------------
void LHE2HepMCConverter::beginRun(edm::Run const& iRun, edm::EventSetup const&) {
  edm::Handle<LHERunInfoProduct> lheRunSrcHandle = iRun.getHandle(runInfoProductToken_);
  if (lheRunSrcHandle.isValid()) {
    _lheRunSrc = lheRunSrcHandle.product();
  } else {
    if (_lheRunSrcTag.label() != "source") {
      iRun.getByLabel("source", lheRunSrcHandle);
      if (lheRunSrcHandle.isValid()) {
        _lheRunSrc = lheRunSrcHandle.product();
        edm::LogInfo("LHE2HepMCConverter") << "Taking LHERunInfoproduct from source";
      } else
        edm::LogWarning("LHE2HepMCConverter") << "No LHERunInfoProduct from source";
    }
  }
}

// ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
void LHE2HepMCConverter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
  //The following says we do not know what parameters are allowed so do no validation
  // Please change this to state exactly what you do use, even if it is no parameters
  edm::ParameterSetDescription desc;
  desc.setUnknown();
  descriptions.addDefault(desc);
}

//define this as a plug-in
DEFINE_FWK_MODULE(LHE2HepMCConverter);