LHECOMWeightProducer

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
#include <iostream>
#include <memory>

#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 "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
#include "SimDataFormats/GeneratorProducts/interface/LHERunInfoProduct.h"
#include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"

#include <sstream>
//class that reweights a pure parton level event from the originale COM energy to
//an energy that is < than original COM energy.

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

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

  const edm::InputTag lheTag_;
  const double _newECMS;
  const edm::EDGetTokenT<LHERunInfoProduct> tokenLHERun_;
  const edm::EDGetTokenT<LHEEventProduct> tokenLHEEvent_;
  int _pdfset;
  int _pdfmember;
  double _origECMS;
  std::string _label;
};

#include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
namespace LHAPDF {
  void initPDFSet(int nset, int setid, int member = 0);
  int numberPDF(int nset);
  void usePDFMember(int nset, int member);
  double xfx(int nset, double x, double Q, int fl);
  double getXmin(int nset, int member);
  double getXmax(int nset, int member);
  double getQ2min(int nset, int member);
  double getQ2max(int nset, int member);
  void extrapolate(bool extrapolate = true);
}  // namespace LHAPDF

/////////////////////////////////////////////////////////////////////////////////////
LHECOMWeightProducer::LHECOMWeightProducer(const edm::ParameterSet& pset)
    : lheTag_(pset.getParameter<edm::InputTag>("lheSrc")),
      _newECMS(pset.getParameter<double>("NewECMS")),
      tokenLHERun_(consumes<LHERunInfoProduct, edm::InRun>(lheTag_)),
      tokenLHEEvent_(consumes<LHEEventProduct>(lheTag_)) {
  std::stringstream labelStr;
  labelStr << "com"
           << "To" << _newECMS;
  _label = labelStr.str();
  produces<GenEventInfoProduct>(_label);
}

/////////////////////////////////////////////////////////////////////////////////////
void LHECOMWeightProducer::beginRun(edm::Run const& run, const edm::EventSetup& es) {
  using namespace edm;
  const edm::Handle<LHERunInfoProduct>& lheRun = run.getHandle(tokenLHERun_);
  //assumes the same pdf is used for both beams
  _pdfset = lheRun->heprup().PDFSUP.first;
  _pdfmember = lheRun->heprup().PDFGUP.first;
  _origECMS = lheRun->heprup().EBMUP.first + lheRun->heprup().EBMUP.second;
  edm::LogInfo("LHECOMWeightProducer") << "PDFSET: " << _pdfset << "; member: " << _pdfmember
                                       << "; COM energy: " << _origECMS;
  if (_newECMS > _origECMS)
    throw cms::Exception("LHECOMWeightProducer") << "You cannot reweight COM energy to a higher than original energy ";
  LHAPDF::initPDFSet(1, _pdfset, _pdfmember);
}

/////////////////////////////////////////////////////////////////////////////////////
void LHECOMWeightProducer::beginJob() {
  //LHAPDF::initPDFSet(1,pdfSetName_);
}

/////////////////////////////////////////////////////////////////////////////////////
void LHECOMWeightProducer::produce(edm::Event& iEvent, const edm::EventSetup&) {
  using namespace std;
  bool verbose = false;

  if (iEvent.isRealData())
    return;

  const edm::Handle<LHEEventProduct>& lheevent = iEvent.getHandle(tokenLHEEvent_);

  float Q = lheevent->hepeup().SCALUP;

  int id1 = lheevent->hepeup().IDUP[0];
  double x1 = fabs(lheevent->hepeup().PUP[0][2] / (_origECMS / 2));
  double x1prime = fabs(lheevent->hepeup().PUP[0][2] / (_newECMS / 2));

  int id2 = lheevent->hepeup().IDUP[1];
  double x2 = fabs(lheevent->hepeup().PUP[1][2] / (_origECMS / 2));
  double x2prime = fabs(lheevent->hepeup().PUP[1][2] / (_newECMS / 2));

  LogTrace("LHECOMWeightProducer") << "*******LHECOMWeightProducer*******\n"
                                   << " Q  : " << Q << "\n"
                                   << " id1: " << id1 << "\n"
                                   << " x1 : " << x1 << "\n"
                                   << " x1': " << x1prime << "\n"
                                   << " id2: " << id2 << "\n"
                                   << " x2 : " << x2 << "\n"
                                   << " x2': " << x2prime;
  //gluon is 0 in the LHAPDF numbering
  if (id1 == 21)
    id1 = 0;
  if (id2 == 21)
    id2 = 0;

  // Put PDF weights in the event
  if (verbose)
    cout << " Set : " << _pdfset << "  member : " << _pdfmember << endl;

  LHAPDF::usePDFMember(1, _pdfmember);
  double oldpdf1 = LHAPDF::xfx(1, x1, Q, id1) / x1;
  double oldpdf2 = LHAPDF::xfx(1, x2, Q, id2) / x2;
  double newpdf1 = LHAPDF::xfx(1, x1prime, Q, id1) / x1prime;
  double newpdf2 = LHAPDF::xfx(1, x2prime, Q, id2) / x2prime;
  LogTrace("LHECOMWeightProducer") << "     xfx1 : " << oldpdf1 << "\n"
                                   << "     xfx2 : " << oldpdf2 << "\n"
                                   << "     xfx1': " << newpdf1 << "\n"
                                   << "     xfx2': " << newpdf2 << "\n"
                                   << "     weight:" << (newpdf1 / oldpdf1) * (newpdf2 / oldpdf2);
  double weight = (newpdf1 / oldpdf1) * (newpdf2 / oldpdf2);
  std::vector<double> weights;
  weights.push_back(weight);
  std::unique_ptr<GenEventInfoProduct> info(new GenEventInfoProduct());
  info->setWeights(weights);
  iEvent.put(std::move(info), _label);
}

DEFINE_FWK_MODULE(LHECOMWeightProducer);