Go to the documentation of this file.00001 #include <memory>
00002
00003 #include "FWCore/Framework/interface/Frameworkfwd.h"
00004 #include "FWCore/Framework/interface/EDProducer.h"
00005
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/Framework/interface/Run.h"
00008 #include "FWCore/Framework/interface/MakerMacros.h"
00009
00010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00011
00012 #include "SimDataFormats/GeneratorProducts/interface/LHEEventProduct.h"
00013 #include "SimDataFormats/GeneratorProducts/interface/LHERunInfoProduct.h"
00014 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00015
00016 #include <sstream>
00017
00018
00019
00020
00021
00022
00023 class LHECOMWeightProducer : public edm::EDProducer {
00024 public:
00025 explicit LHECOMWeightProducer(const edm::ParameterSet&);
00026 ~LHECOMWeightProducer();
00027
00028 private:
00029 virtual void beginJob() override;
00030 virtual void beginRun(edm::Run const& run, const edm::EventSetup &es) override;
00031 virtual void produce(edm::Event&, const edm::EventSetup&) override;
00032
00033 edm::InputTag lheTag_;
00034 int _pdfset;
00035 int _pdfmember;
00036 double _origECMS;
00037 double _newECMS;
00038 std::string _label;
00039 };
00040
00041 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00042 namespace LHAPDF {
00043 void initPDFSet(int nset, int setid, int member=0);
00044 int numberPDF(int nset);
00045 void usePDFMember(int nset, int member);
00046 double xfx(int nset, double x, double Q, int fl);
00047 double getXmin(int nset, int member);
00048 double getXmax(int nset, int member);
00049 double getQ2min(int nset, int member);
00050 double getQ2max(int nset, int member);
00051 void extrapolate(bool extrapolate=true);
00052 }
00053
00055 LHECOMWeightProducer::LHECOMWeightProducer(const edm::ParameterSet& pset) :
00056 lheTag_(pset.getParameter<edm::InputTag> ("lheSrc")),
00057 _newECMS(pset.getParameter< double > ("NewECMS"))
00058 {
00059 std::stringstream labelStr;
00060 labelStr << "com" << "To" << _newECMS;
00061 _label = labelStr.str();
00062 produces<GenEventInfoProduct>(_label);
00063 }
00064
00066 LHECOMWeightProducer::~LHECOMWeightProducer(){}
00067
00069 void LHECOMWeightProducer::beginRun(edm::Run const& run, const edm::EventSetup &es){
00070 using namespace edm;
00071 Handle<LHERunInfoProduct> lheRun;
00072 run.getByLabel(lheTag_, lheRun);
00073
00074 _pdfset = lheRun->heprup().PDFSUP.first;
00075 _pdfmember = lheRun->heprup().PDFGUP.first;
00076 _origECMS = lheRun->heprup().EBMUP.first + lheRun->heprup().EBMUP.second;
00077 edm::LogInfo("LHECOMWeightProducer") << "PDFSET: " << _pdfset << "; member: " << _pdfmember << "; COM energy: " << _origECMS;
00078 if ( _newECMS > _origECMS )
00079 throw cms::Exception("LHECOMWeightProducer") << "You cannot reweight COM energy to a higher than original energy ";
00080 LHAPDF::initPDFSet(1,_pdfset, _pdfmember);
00081 }
00082
00083
00085 void LHECOMWeightProducer::beginJob() {
00086
00087 }
00088
00090 void LHECOMWeightProducer::produce(edm::Event& iEvent, const edm::EventSetup&) {
00091
00092 using namespace std;
00093 bool verbose = false;
00094
00095 if (iEvent.isRealData()) return;
00096
00097 edm::Handle<LHEEventProduct> lheevent;
00098 iEvent.getByLabel(lheTag_, lheevent);
00099
00100 float Q = lheevent->hepeup().SCALUP;
00101
00102 int id1 = lheevent->hepeup().IDUP[0];
00103 double x1 = fabs(lheevent->hepeup().PUP[0][2]/(_origECMS/2));
00104 double x1prime = fabs(lheevent->hepeup().PUP[0][2]/(_newECMS/2));
00105
00106 int id2 = lheevent->hepeup().IDUP[1];
00107 double x2 = fabs(lheevent->hepeup().PUP[1][2]/(_origECMS/2));
00108 double x2prime = fabs(lheevent->hepeup().PUP[1][2]/(_newECMS/2));
00109
00110 LogTrace("LHECOMWeightProducer") << "*******LHECOMWeightProducer*******\n" <<
00111 " Q : " << Q << "\n" <<
00112 " id1: " << id1 << "\n" <<
00113 " x1 : " << x1 << "\n" <<
00114 " x1': " << x1prime << "\n" <<
00115 " id2: " << id2 << "\n" <<
00116 " x2 : " << x2 << "\n" <<
00117 " x2': " << x2prime ;
00118
00119 if (id1 == 21)
00120 id1 = 0;
00121 if (id2 == 21)
00122 id2 = 0;
00123
00124
00125 if (verbose)
00126 cout << " Set : " << _pdfset << " member : " << _pdfmember << endl;
00127
00128 LHAPDF::usePDFMember(1,_pdfmember);
00129 double oldpdf1 = LHAPDF::xfx(1, x1, Q, id1)/x1;
00130 double oldpdf2 = LHAPDF::xfx(1, x2, Q, id2)/x2;
00131 double newpdf1 = LHAPDF::xfx(1, x1prime, Q, id1)/x1prime;
00132 double newpdf2 = LHAPDF::xfx(1, x2prime, Q, id2)/x2prime;
00133 LogTrace("LHECOMWeightProducer") <<
00134 " xfx1 : " << oldpdf1 << "\n" <<
00135 " xfx2 : " << oldpdf2 << "\n" <<
00136 " xfx1': " << newpdf1 << "\n" <<
00137 " xfx2': " << newpdf2 << "\n" <<
00138 " weight:" << (newpdf1/oldpdf1)*(newpdf2/oldpdf2);
00139 double weight = (newpdf1/oldpdf1)*(newpdf2/oldpdf2);
00140 std::vector<double> weights;
00141 weights.push_back(weight);
00142 std::auto_ptr<GenEventInfoProduct> info(new GenEventInfoProduct());
00143 info->setWeights(weights);
00144 iEvent.put(info, _label);
00145 }
00146
00147 DEFINE_FWK_MODULE(LHECOMWeightProducer);