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() ;
00030 virtual void beginRun(edm::Run &run, const edm::EventSetup &es) ;
00031 virtual void produce(edm::Event&, const edm::EventSetup&);
00032 virtual void endJob() ;
00033
00034 edm::InputTag lheTag_;
00035 int _pdfset;
00036 int _pdfmember;
00037 double _origECMS;
00038 double _newECMS;
00039 std::string _label;
00040 };
00041
00042 #include "SimDataFormats/GeneratorProducts/interface/GenEventInfoProduct.h"
00043 namespace LHAPDF {
00044 void initPDFSet(int nset, int setid, int member=0);
00045 int numberPDF(int nset);
00046 void usePDFMember(int nset, int member);
00047 double xfx(int nset, double x, double Q, int fl);
00048 double getXmin(int nset, int member);
00049 double getXmax(int nset, int member);
00050 double getQ2min(int nset, int member);
00051 double getQ2max(int nset, int member);
00052 void extrapolate(bool extrapolate=true);
00053 }
00054
00056 LHECOMWeightProducer::LHECOMWeightProducer(const edm::ParameterSet& pset) :
00057 lheTag_(pset.getParameter<edm::InputTag> ("lheSrc")),
00058 _newECMS(pset.getParameter< double > ("NewECMS"))
00059 {
00060 std::stringstream labelStr;
00061 labelStr << "com" << "To" << _newECMS;
00062 _label = labelStr.str();
00063 produces<GenEventInfoProduct>(_label);
00064 }
00065
00067 LHECOMWeightProducer::~LHECOMWeightProducer(){}
00068
00070 void LHECOMWeightProducer::beginRun(edm::Run &run, const edm::EventSetup &es){
00071 using namespace edm;
00072 Handle<LHERunInfoProduct> lheRun;
00073 run.getByLabel(lheTag_, lheRun);
00074
00075 _pdfset = lheRun->heprup().PDFSUP.first;
00076 _pdfmember = lheRun->heprup().PDFGUP.first;
00077 _origECMS = lheRun->heprup().EBMUP.first + lheRun->heprup().EBMUP.second;
00078 edm::LogInfo("LHECOMWeightProducer") << "PDFSET: " << _pdfset << "; member: " << _pdfmember << "; COM energy: " << _origECMS;
00079 if ( _newECMS > _origECMS )
00080 throw cms::Exception("LHECOMWeightProducer") << "You cannot reweight COM energy to a higher than original energy ";
00081 LHAPDF::initPDFSet(1,_pdfset, _pdfmember);
00082 }
00083
00084
00086 void LHECOMWeightProducer::beginJob() {
00087
00088 }
00089
00091 void LHECOMWeightProducer::endJob(){}
00092
00094 void LHECOMWeightProducer::produce(edm::Event& iEvent, const edm::EventSetup&) {
00095
00096 using namespace std;
00097 bool verbose = false;
00098
00099 if (iEvent.isRealData()) return;
00100
00101 edm::Handle<LHEEventProduct> lheevent;
00102 iEvent.getByLabel(lheTag_, lheevent);
00103
00104 float Q = lheevent->hepeup().SCALUP;
00105
00106 int id1 = lheevent->hepeup().IDUP[0];
00107 double x1 = fabs(lheevent->hepeup().PUP[0][2]/(_origECMS/2));
00108 double x1prime = fabs(lheevent->hepeup().PUP[0][2]/(_newECMS/2));
00109
00110 int id2 = lheevent->hepeup().IDUP[1];
00111 double x2 = fabs(lheevent->hepeup().PUP[1][2]/(_origECMS/2));
00112 double x2prime = fabs(lheevent->hepeup().PUP[1][2]/(_newECMS/2));
00113
00114 LogTrace("LHECOMWeightProducer") << "*******LHECOMWeightProducer*******\n" <<
00115 " Q : " << Q << "\n" <<
00116 " id1: " << id1 << "\n" <<
00117 " x1 : " << x1 << "\n" <<
00118 " x1': " << x1prime << "\n" <<
00119 " id2: " << id2 << "\n" <<
00120 " x2 : " << x2 << "\n" <<
00121 " x2': " << x2prime ;
00122
00123 if (id1 == 21)
00124 id1 = 0;
00125 if (id2 == 21)
00126 id2 = 0;
00127
00128
00129 if (verbose)
00130 cout << " Set : " << _pdfset << " member : " << _pdfmember << endl;
00131
00132 LHAPDF::usePDFMember(1,_pdfmember);
00133 double oldpdf1 = LHAPDF::xfx(1, x1, Q, id1)/x1;
00134 double oldpdf2 = LHAPDF::xfx(1, x2, Q, id2)/x2;
00135 double newpdf1 = LHAPDF::xfx(1, x1prime, Q, id1)/x1prime;
00136 double newpdf2 = LHAPDF::xfx(1, x2prime, Q, id2)/x2prime;
00137 LogTrace("LHECOMWeightProducer") <<
00138 " xfx1 : " << oldpdf1 << "\n" <<
00139 " xfx2 : " << oldpdf2 << "\n" <<
00140 " xfx1': " << newpdf1 << "\n" <<
00141 " xfx2': " << newpdf2 << "\n" <<
00142 " weight:" << (newpdf1/oldpdf1)*(newpdf2/oldpdf2);
00143 double weight = (newpdf1/oldpdf1)*(newpdf2/oldpdf2);
00144 std::vector<double> weights;
00145 weights.push_back(weight);
00146 std::auto_ptr<GenEventInfoProduct> info(new GenEventInfoProduct());
00147 info->setWeights(weights);
00148 iEvent.put(info, _label);
00149 }
00150
00151 DEFINE_FWK_MODULE(LHECOMWeightProducer);