CMS 3D CMS Logo

LHECOMWeightProducer.cc
Go to the documentation of this file.
1 #include <iostream>
2 #include <memory>
3 
6 
10 
12 
16 
17 #include <sstream>
18 //class that reweights a pure parton level event from the originale COM energy to
19 //an energy that is < than original COM energy.
20 
21 //
22 // class declaration
23 //
24 class LHECOMWeightProducer : public edm::one::EDProducer<edm::one::WatchRuns> {
25 public:
26  explicit LHECOMWeightProducer(const edm::ParameterSet&);
27  ~LHECOMWeightProducer() override = default;
28 
29 private:
30  void beginJob() override;
31  void beginRun(edm::Run const& run, const edm::EventSetup& es) override;
32  void endRun(edm::Run const& run, const edm::EventSetup& es) override {}
33  void produce(edm::Event&, const edm::EventSetup&) override;
34 
36  const double _newECMS;
39  int _pdfset;
41  double _origECMS;
43 };
44 
46 namespace LHAPDF {
47  void initPDFSet(int nset, int setid, int member = 0);
48  int numberPDF(int nset);
49  void usePDFMember(int nset, int member);
50  double xfx(int nset, double x, double Q, int fl);
51  double getXmin(int nset, int member);
52  double getXmax(int nset, int member);
53  double getQ2min(int nset, int member);
54  double getQ2max(int nset, int member);
55  void extrapolate(bool extrapolate = true);
56 } // namespace LHAPDF
57 
60  : lheTag_(pset.getParameter<edm::InputTag>("lheSrc")),
61  _newECMS(pset.getParameter<double>("NewECMS")),
62  tokenLHERun_(consumes<LHERunInfoProduct, edm::InRun>(lheTag_)),
63  tokenLHEEvent_(consumes<LHEEventProduct>(lheTag_)) {
64  std::stringstream labelStr;
65  labelStr << "com"
66  << "To" << _newECMS;
67  _label = labelStr.str();
68  produces<GenEventInfoProduct>(_label);
69 }
70 
73  using namespace edm;
74  const edm::Handle<LHERunInfoProduct>& lheRun = run.getHandle(tokenLHERun_);
75  //assumes the same pdf is used for both beams
76  _pdfset = lheRun->heprup().PDFSUP.first;
77  _pdfmember = lheRun->heprup().PDFGUP.first;
78  _origECMS = lheRun->heprup().EBMUP.first + lheRun->heprup().EBMUP.second;
79  edm::LogInfo("LHECOMWeightProducer") << "PDFSET: " << _pdfset << "; member: " << _pdfmember
80  << "; COM energy: " << _origECMS;
81  if (_newECMS > _origECMS)
82  throw cms::Exception("LHECOMWeightProducer") << "You cannot reweight COM energy to a higher than original energy ";
84 }
85 
88  //LHAPDF::initPDFSet(1,pdfSetName_);
89 }
90 
93  using namespace std;
94  bool verbose = false;
95 
96  if (iEvent.isRealData())
97  return;
98 
99  const edm::Handle<LHEEventProduct>& lheevent = iEvent.getHandle(tokenLHEEvent_);
100 
101  float Q = lheevent->hepeup().SCALUP;
102 
103  int id1 = lheevent->hepeup().IDUP[0];
104  double x1 = fabs(lheevent->hepeup().PUP[0][2] / (_origECMS / 2));
105  double x1prime = fabs(lheevent->hepeup().PUP[0][2] / (_newECMS / 2));
106 
107  int id2 = lheevent->hepeup().IDUP[1];
108  double x2 = fabs(lheevent->hepeup().PUP[1][2] / (_origECMS / 2));
109  double x2prime = fabs(lheevent->hepeup().PUP[1][2] / (_newECMS / 2));
110 
111  LogTrace("LHECOMWeightProducer") << "*******LHECOMWeightProducer*******\n"
112  << " Q : " << Q << "\n"
113  << " id1: " << id1 << "\n"
114  << " x1 : " << x1 << "\n"
115  << " x1': " << x1prime << "\n"
116  << " id2: " << id2 << "\n"
117  << " x2 : " << x2 << "\n"
118  << " x2': " << x2prime;
119  //gluon is 0 in the LHAPDF numbering
120  if (id1 == 21)
121  id1 = 0;
122  if (id2 == 21)
123  id2 = 0;
124 
125  // Put PDF weights in the event
126  if (verbose)
127  cout << " Set : " << _pdfset << " member : " << _pdfmember << endl;
128 
130  double oldpdf1 = LHAPDF::xfx(1, x1, Q, id1) / x1;
131  double oldpdf2 = LHAPDF::xfx(1, x2, Q, id2) / x2;
132  double newpdf1 = LHAPDF::xfx(1, x1prime, Q, id1) / x1prime;
133  double newpdf2 = LHAPDF::xfx(1, x2prime, Q, id2) / x2prime;
134  LogTrace("LHECOMWeightProducer") << " xfx1 : " << oldpdf1 << "\n"
135  << " xfx2 : " << oldpdf2 << "\n"
136  << " xfx1': " << newpdf1 << "\n"
137  << " xfx2': " << newpdf2 << "\n"
138  << " weight:" << (newpdf1 / oldpdf1) * (newpdf2 / oldpdf2);
139  double weight = (newpdf1 / oldpdf1) * (newpdf2 / oldpdf2);
140  std::vector<double> weights;
141  weights.push_back(weight);
142  std::unique_ptr<GenEventInfoProduct> info(new GenEventInfoProduct());
143  info->setWeights(weights);
144  iEvent.put(std::move(info), _label);
145 }
146 
~LHECOMWeightProducer() override=default
static const TGPicture * info(bool iBackgroundIsBlack)
bool verbose
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::pair< double, double > EBMUP
Definition: LesHouches.h:82
Definition: weight.py:1
#define LogTrace(id)
std::pair< int, int > PDFGUP
Definition: LesHouches.h:88
double getXmax(int nset, int member)
double getQ2max(int nset, int member)
int iEvent
Definition: GenABIO.cc:224
const lhef::HEPRUP & heprup() const
void endRun(edm::Run const &run, const edm::EventSetup &es) override
void extrapolate(bool extrapolate=true)
std::vector< FiveVector > PUP
Definition: LesHouches.h:246
void beginRun(edm::Run const &run, const edm::EventSetup &es) override
int numberPDF(int nset)
double getXmin(int nset, int member)
const edm::InputTag lheTag_
std::vector< int > IDUP
Definition: LesHouches.h:223
Log< level::Info, false > LogInfo
void initPDFSet(int nset, int setid, int member=0)
const edm::EDGetTokenT< LHEEventProduct > tokenLHEEvent_
const edm::EDGetTokenT< LHERunInfoProduct > tokenLHERun_
std::pair< int, int > PDFSUP
Definition: LesHouches.h:94
LHECOMWeightProducer(const edm::ParameterSet &)
HLT enums.
double xfx(int nset, double x, double Q, int fl)
float x
void produce(edm::Event &, const edm::EventSetup &) override
double getQ2min(int nset, int member)
def move(src, dest)
Definition: eostools.py:511
Definition: Run.h:45
void usePDFMember(int nset, int member)
double SCALUP
Definition: LesHouches.h:208
const lhef::HEPEUP & hepeup() const