Herwig7HepMC3Hadronizer

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 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250
#include <memory>
#include <sstream>
#include <fstream>

#include <HepMC3/GenEvent.h>
//#include <HepMC/IO_BaseClass.h>
#include "HepMC3/Print.h"

#include <ThePEG/Repository/Repository.h>
#include <ThePEG/EventRecord/Event.h>
#include <ThePEG/Config/ThePEG.h>
#include <ThePEG/LesHouches/LesHouchesReader.h>

#include "FWCore/Framework/interface/LuminosityBlock.h"
#include "FWCore/Framework/interface/EventSetup.h"
#include "FWCore/MessageLogger/interface/MessageLogger.h"

#include "SimDataFormats/GeneratorProducts/interface/HepMC3Product.h"
#include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
#include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct3.h"
#include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"

#include "GeneratorInterface/Core/interface/BaseHadronizer.h"
#include "GeneratorInterface/Core/interface/GeneratorFilter.h"
#include "GeneratorInterface/Core/interface/HadronizerFilter.h"

#include "GeneratorInterface/LHEInterface/interface/LHEProxy.h"

#include "GeneratorInterface/Herwig7Interface/interface/Herwig7HepMC3Interface.h"

#include <Herwig/API/HerwigAPI.h>
#include "CLHEP/Random/RandomEngine.h"

namespace CLHEP {
  class HepRandomEngine;
}

class Herwig7HepMC3Hadronizer : public Herwig7HepMC3Interface, public gen::BaseHadronizer {
public:
  Herwig7HepMC3Hadronizer(const edm::ParameterSet& params);
  ~Herwig7HepMC3Hadronizer() override;

  bool readSettings(int) { return true; }
  bool initializeForInternalPartons();
  bool initializeForExternalPartons();
  bool declareStableParticles(const std::vector<int>& pdgIds);
  bool declareSpecialSettings(const std::vector<std::string>) { return true; }

  void statistics();

  bool generatePartonsAndHadronize();
  bool hadronize();
  bool decay();
  bool residualDecay();
  void finalizeEvent();

  const char* classname() const { return "Herwig7HepMC3Hadronizer"; }
  std::unique_ptr<GenLumiInfoHeader> getGenLumiInfoHeader() const override;
  void beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&);
  void randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine);

private:
  void doSetRandomEngine(CLHEP::HepRandomEngine* v) override { setPEGRandomEngine(v); }

  unsigned int eventsToPrint;

  ThePEG::EventPtr thepegEvent;
  bool haveEvt = false;

  std::shared_ptr<lhef::LHEProxy> proxy_;
  const std::string handlerDirectory_;
  edm::ParameterSet paramSettings;
  const std::string runFileName;

  unsigned int firstLumiBlock = 0;
  unsigned int currentLumiBlock = 0;
};

Herwig7HepMC3Hadronizer::Herwig7HepMC3Hadronizer(const edm::ParameterSet& pset)
    : Herwig7HepMC3Interface(pset),
      BaseHadronizer(pset),
      eventsToPrint(pset.getUntrackedParameter<unsigned int>("eventsToPrint", 0)),
      handlerDirectory_(pset.getParameter<std::string>("eventHandlers")),
      runFileName(pset.getParameter<std::string>("run")) {
  ivhepmc = 3;
  initRepository(pset);
  paramSettings = pset;
}

Herwig7HepMC3Hadronizer::~Herwig7HepMC3Hadronizer() {}

bool Herwig7HepMC3Hadronizer::initializeForInternalPartons() {
  if (currentLumiBlock == firstLumiBlock) {
    std::ifstream runFile(runFileName + ".run");
    if (runFile.fail())  //required for showering of LHE files
    {
      initRepository(paramSettings);
    }
    if (!initGenerator()) {
      edm::LogInfo("Generator|Herwig7HepMC3Hadronizer") << "No run step for Herwig chosen. Program will be aborted.";
      exit(0);
    }
  }
  return true;
}

bool Herwig7HepMC3Hadronizer::initializeForExternalPartons() {
  if (currentLumiBlock == firstLumiBlock) {
    std::ifstream runFile(runFileName + ".run");
    if (runFile.fail())  //required for showering of LHE files
    {
      initRepository(paramSettings);
    }
    if (!initGenerator()) {
      edm::LogInfo("Generator|Herwig7HepMC3Hadronizer") << "No run step for Herwig chosen. Program will be aborted.";
      exit(0);
    }
  }
  return true;
}

bool Herwig7HepMC3Hadronizer::declareStableParticles(const std::vector<int>& pdgIds) { return false; }

void Herwig7HepMC3Hadronizer::statistics() {
  if (eg_) {
    runInfo().setInternalXSec(
        GenRunInfoProduct::XSec(eg_->integratedXSec() / ThePEG::picobarn, eg_->integratedXSecErr() / ThePEG::picobarn));
  }
}

bool Herwig7HepMC3Hadronizer::generatePartonsAndHadronize() {
  edm::LogInfo("Generator|Herwig7HepMC3Hadronizer") << "Start production";

  try {
    thepegEvent = eg_->shoot();
  } catch (std::exception& exc) {
    edm::LogWarning("Generator|Herwig7HepMC3Hadronizer")
        << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
    return false;
  }

  if (!thepegEvent) {
    edm::LogWarning("Generator|Herwig7HepMC3Hadronizer") << "thepegEvent not initialized";
    return false;
  }

  event3() = convert(thepegEvent);
  if (!event3().get()) {
    edm::LogWarning("Generator|Herwig7HepMC3Hadronizer") << "genEvent not initialized";
    return false;
  }

  return true;
}

bool Herwig7HepMC3Hadronizer::hadronize() {
  if (!haveEvt) {
    try {
      thepegEvent = eg_->shoot();
      haveEvt = true;
    } catch (std::exception& exc) {
      edm::LogWarning("Generator|Herwig7HepMC3Hadronizer")
          << "EGPtr::shoot() thrown an exception, event skipped: " << exc.what();
      return false;
    }
  }
  int evtnum = lheEvent()->evtnum();
  if (evtnum == -1) {
    edm::LogError("Generator|Herwig7HepMC3Hadronizer")
        << "Event number not set in lhe file, needed for correctly aligning Herwig and LHE events!";
    return false;
  }
  if (thepegEvent->number() < evtnum) {
    edm::LogError("Herwig7 interface") << "Herwig does not seem to be generating events in order, did you set "
                                          "/Herwig/EventHandlers/FxFxLHReader:AllowedToReOpen Yes?";
    return false;
  } else if (thepegEvent->number() == evtnum) {
    haveEvt = false;
    if (!thepegEvent) {
      edm::LogWarning("Generator|Herwig7HepMC3Hadronizer") << "thepegEvent not initialized";
      return false;
    }

    event3() = convert(thepegEvent);
    if (!event3().get()) {
      edm::LogWarning("Generator|Herwig7HepMC3Hadronizer") << "genEvent not initialized";
      return false;
    }
    return true;
  }
  edm::LogWarning("Generator|Herwig7HepMC3Hadronizer")
      << "Event " << evtnum << " not generated (likely skipped in merging)";
  return false;
}

void Herwig7HepMC3Hadronizer::finalizeEvent() {
  eventInfo3() = std::make_unique<GenEventInfoProduct3>(event3().get());
  eventInfo3()->setBinningValues(std::vector<double>(1, pthat(thepegEvent)));

  if (eventsToPrint) {
    eventsToPrint--;
    //event3()->print();
    HepMC3::Print::listing(*(event3().get()));
  }

  //  if (iobc_.get())
  //    iobc_->write_event(event().get());

  edm::LogInfo("Generator|Herwig7HepMC3Hadronizer") << "Event produced";
}

bool Herwig7HepMC3Hadronizer::decay() { return true; }

bool Herwig7HepMC3Hadronizer::residualDecay() { return true; }

std::unique_ptr<GenLumiInfoHeader> Herwig7HepMC3Hadronizer::getGenLumiInfoHeader() const {
  auto genLumiInfoHeader = BaseHadronizer::getGenLumiInfoHeader();

  if (thepegEvent) {
    int weights_number = thepegEvent->optionalWeights().size();

    if (weights_number > 1) {
      genLumiInfoHeader->weightNames().reserve(weights_number + 1);
      genLumiInfoHeader->weightNames().push_back("nominal");
      std::map<std::string, double> weights_map = thepegEvent->optionalWeights();
      for (std::map<std::string, double>::iterator it = weights_map.begin(); it != weights_map.end(); it++) {
        genLumiInfoHeader->weightNames().push_back(it->first);
      }
    }
  }

  return genLumiInfoHeader;
}

void Herwig7HepMC3Hadronizer::randomizeIndex(edm::LuminosityBlock const& lumi, CLHEP::HepRandomEngine* rengine) {
  BaseHadronizer::randomizeIndex(lumi, rengine);

  if (firstLumiBlock == 0) {
    firstLumiBlock = lumi.id().luminosityBlock();
  }
  currentLumiBlock = lumi.id().luminosityBlock();
}

#include "GeneratorInterface/ExternalDecays/interface/ExternalDecayDriver.h"

typedef edm::GeneratorFilter<Herwig7HepMC3Hadronizer, gen::ExternalDecayDriver> Herwig7HepMC3GeneratorFilter;
DEFINE_FWK_MODULE(Herwig7HepMC3GeneratorFilter);

typedef edm::HadronizerFilter<Herwig7HepMC3Hadronizer, gen::ExternalDecayDriver> Herwig7HepMC3HadronizerFilter;
DEFINE_FWK_MODULE(Herwig7HepMC3HadronizerFilter);