CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
LHECOMWeightProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
9 
11 
15 
16 #include <sstream>
17 //class that reweights a pure parton level event from the originale COM energy to
18 //an energy that is < than original COM energy.
19 
20 //
21 // class declaration
22 //
24  public:
25  explicit LHECOMWeightProducer(const edm::ParameterSet&);
27 
28  private:
29  virtual void beginJob() override;
30  virtual void beginRun(edm::Run const& run, const edm::EventSetup &es) override;
31  virtual void produce(edm::Event&, const edm::EventSetup&) override;
32 
34  int _pdfset;
36  double _origECMS;
37  double _newECMS;
39 };
40 
42 namespace LHAPDF {
43  void initPDFSet(int nset, int setid, int member=0);
44  int numberPDF(int nset);
45  void usePDFMember(int nset, int member);
46  double xfx(int nset, double x, double Q, int fl);
47  double getXmin(int nset, int member);
48  double getXmax(int nset, int member);
49  double getQ2min(int nset, int member);
50  double getQ2max(int nset, int member);
51  void extrapolate(bool extrapolate=true);
52 }
53 
56  lheTag_(pset.getParameter<edm::InputTag> ("lheSrc")),
57  _newECMS(pset.getParameter< double > ("NewECMS"))
58 {
59  std::stringstream labelStr;
60  labelStr << "com" << "To" << _newECMS;
61  _label = labelStr.str();
62  produces<GenEventInfoProduct>(_label);
63 }
64 
67 
70  using namespace edm;
72  run.getByLabel(lheTag_, lheRun);
73  //assumes the same pdf is used for both beams
74  _pdfset = lheRun->heprup().PDFSUP.first;
75  _pdfmember = lheRun->heprup().PDFGUP.first;
76  _origECMS = lheRun->heprup().EBMUP.first + lheRun->heprup().EBMUP.second;
77  edm::LogInfo("LHECOMWeightProducer") << "PDFSET: " << _pdfset << "; member: " << _pdfmember << "; COM energy: " << _origECMS;
78  if ( _newECMS > _origECMS )
79  throw cms::Exception("LHECOMWeightProducer") << "You cannot reweight COM energy to a higher than original energy ";
81 }
82 
83 
86  //LHAPDF::initPDFSet(1,pdfSetName_);
87 }
88 
91 
92  using namespace std;
93  bool verbose = false;
94 
95  if (iEvent.isRealData()) return;
96 
98  iEvent.getByLabel(lheTag_, lheevent);
99 
100  float Q = lheevent->hepeup().SCALUP;
101 
102  int id1 = lheevent->hepeup().IDUP[0];
103  double x1 = fabs(lheevent->hepeup().PUP[0][2]/(_origECMS/2));
104  double x1prime = fabs(lheevent->hepeup().PUP[0][2]/(_newECMS/2));
105 
106  int id2 = lheevent->hepeup().IDUP[1];
107  double x2 = fabs(lheevent->hepeup().PUP[1][2]/(_origECMS/2));
108  double x2prime = fabs(lheevent->hepeup().PUP[1][2]/(_newECMS/2));
109 
110  LogTrace("LHECOMWeightProducer") << "*******LHECOMWeightProducer*******\n" <<
111  " Q : " << Q << "\n" <<
112  " id1: " << id1 << "\n" <<
113  " x1 : " << x1 << "\n" <<
114  " x1': " << x1prime << "\n" <<
115  " id2: " << id2 << "\n" <<
116  " x2 : " << x2 << "\n" <<
117  " x2': " << x2prime ;
118  //gluon is 0 in the LHAPDF numbering
119  if (id1 == 21)
120  id1 = 0;
121  if (id2 == 21)
122  id2 = 0;
123 
124  // Put PDF weights in the event
125  if (verbose)
126  cout << " Set : " << _pdfset << " member : " << _pdfmember << endl;
127 
129  double oldpdf1 = LHAPDF::xfx(1, x1, Q, id1)/x1;
130  double oldpdf2 = LHAPDF::xfx(1, x2, Q, id2)/x2;
131  double newpdf1 = LHAPDF::xfx(1, x1prime, Q, id1)/x1prime;
132  double newpdf2 = LHAPDF::xfx(1, x2prime, Q, id2)/x2prime;
133  LogTrace("LHECOMWeightProducer") <<
134  " 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::auto_ptr<GenEventInfoProduct> info(new GenEventInfoProduct());
143  info->setWeights(weights);
144  iEvent.put(info, _label);
145 }
146 
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:176
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool isRealData() const
Definition: EventBase.h:60
void initPDFSet(int nset, const std::string &filename, int member=0)
double getXmax(int nset, int member)
double getQ2max(int nset, int member)
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
void extrapolate(bool extrapolate=true)
virtual void beginRun(edm::Run const &run, const edm::EventSetup &es) override
int numberPDF(int nset)
double getXmin(int nset, int member)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
#define LogTrace(id)
LHECOMWeightProducer(const edm::ParameterSet &)
double xfx(int nset, double x, double Q, int fl)
virtual void produce(edm::Event &, const edm::EventSetup &) override
tuple cout
Definition: gather_cfg.py:121
double getQ2min(int nset, int member)
Definition: DDAxes.h:10
virtual void beginJob() override
int weight
Definition: histoStyle.py:50
Definition: Run.h:36
void usePDFMember(int nset, int member)