CMS 3D CMS Logo

PdfWeightProducer.cc
Go to the documentation of this file.
1 #include <memory>
2 
5 
8 
10 
16 #include <Math/VectorUtil.h>
17 
18 //
19 // class declaration
20 //
22 public:
23  explicit PdfWeightProducer(const edm::ParameterSet&);
24  ~PdfWeightProducer() override;
25 
26 private:
27  void beginJob() override;
28  void produce(edm::Event&, const edm::EventSetup&) override;
29  void endJob() override;
30 
37  std::vector<std::string> pdfSetNames_;
38  std::vector<std::string> pdfShortNames_;
39 };
40 
41 namespace LHAPDF {
42  void initPDFSet(int nset, const std::string& filename, int member = 0);
43  int numberPDF(int nset);
44  void usePDFMember(int nset, int member);
45  double xfx(int nset, double x, double Q, int fl);
46  double getXmin(int nset, int member);
47  double getXmax(int nset, int member);
48  double getQ2min(int nset, int member);
49  double getQ2max(int nset, int member);
50  void extrapolate(bool extrapolate = true);
51 } // namespace LHAPDF
52 
55  : fixPOWHEG_(pset.getUntrackedParameter<std::string>("FixPOWHEG", "")),
56  useFirstAsDefault_(pset.getUntrackedParameter<bool>("useFirstAsDefault", false)),
57  genTag_(pset.getUntrackedParameter<edm::InputTag>("GenTag", edm::InputTag("genParticles"))),
58  genToken_(mayConsume<reco::GenParticleCollection>(genTag_)),
59  pdfInfoTag_(pset.getUntrackedParameter<edm::InputTag>("PdfInfoTag", edm::InputTag("generator"))),
60  pdfInfoToken_(consumes<GenEventInfoProduct>(pdfInfoTag_)),
61  pdfSetNames_(pset.getUntrackedParameter<std::vector<std::string> >("PdfSetNames")) {
62  if (!fixPOWHEG_.empty())
63  pdfSetNames_.insert(pdfSetNames_.begin(), fixPOWHEG_);
64 
65  if (pdfSetNames_.size() > 3) {
66  edm::LogWarning("") << pdfSetNames_.size()
67  << " PDF sets requested on input. Using only the first 3 sets and ignoring the rest!!";
68  pdfSetNames_.erase(pdfSetNames_.begin() + 3, pdfSetNames_.end());
69  }
70 
71  for (unsigned int k = 0; k < pdfSetNames_.size(); k++) {
72  size_t dot = pdfSetNames_[k].find_first_of('.');
73  size_t underscore = pdfSetNames_[k].find_first_of('_');
74  if (underscore < dot) {
75  pdfShortNames_.push_back(pdfSetNames_[k].substr(0, underscore));
76  } else {
77  pdfShortNames_.push_back(pdfSetNames_[k].substr(0, dot));
78  }
79  produces<std::vector<double> >(pdfShortNames_[k]);
80  }
81 }
82 
85 
88  for (unsigned int k = 1; k <= pdfSetNames_.size(); k++) {
90  }
91 }
92 
95 
98  if (iEvent.isRealData())
99  return;
100 
102  if (!iEvent.getByToken(pdfInfoToken_, pdfstuff)) {
103  edm::LogError("PDFWeightProducer") << ">>> PdfInfo not found: " << pdfInfoTag_.encode() << " !!!";
104  return;
105  }
106 
107  float Q = pdfstuff->pdf()->scalePDF;
108 
109  int id1 = pdfstuff->pdf()->id.first;
110  double x1 = pdfstuff->pdf()->x.first;
111  double pdf1 = pdfstuff->pdf()->xPDF.first;
112 
113  int id2 = pdfstuff->pdf()->id.second;
114  double x2 = pdfstuff->pdf()->x.second;
115  double pdf2 = pdfstuff->pdf()->xPDF.second;
116  if (useFirstAsDefault_ && pdf1 == -1. && pdf2 == -1.) {
117  LHAPDF::usePDFMember(1, 0);
118  pdf1 = LHAPDF::xfx(1, x1, Q, id1) / x1;
119  pdf2 = LHAPDF::xfx(1, x2, Q, id2) / x2;
120  }
121 
122  // Ad-hoc fix for POWHEG
123  if (!fixPOWHEG_.empty()) {
125  if (!iEvent.getByToken(genToken_, genParticles)) {
126  edm::LogError("PDFWeightProducer") << ">>> genParticles not found: " << genTag_.encode() << " !!!";
127  return;
128  }
129  unsigned int gensize = genParticles->size();
130  double mboson = 0.;
131  for (unsigned int i = 0; i < gensize; ++i) {
132  const reco::GenParticle& part = (*genParticles)[i];
133  int status = part.status();
134  if (status != 3)
135  continue;
136  int id = part.pdgId();
137  if (id != 23 && abs(id) != 24)
138  continue;
139  mboson = part.mass();
140  break;
141  }
142  Q = sqrt(mboson * mboson + Q * Q);
143  LHAPDF::usePDFMember(1, 0);
144  pdf1 = LHAPDF::xfx(1, x1, Q, id1) / x1;
145  pdf2 = LHAPDF::xfx(1, x2, Q, id2) / x2;
146  }
147 
148  // Put PDF weights in the event
149  for (unsigned int k = 1; k <= pdfSetNames_.size(); ++k) {
150  std::unique_ptr<std::vector<double> > weights(new std::vector<double>);
151  unsigned int nweights = 1;
152  if (LHAPDF::numberPDF(k) > 1)
153  nweights += LHAPDF::numberPDF(k);
154  weights->reserve(nweights);
155 
156  for (unsigned int i = 0; i < nweights; ++i) {
158  double newpdf1 = LHAPDF::xfx(k, x1, Q, id1) / x1;
159  double newpdf2 = LHAPDF::xfx(k, x2, Q, id2) / x2;
160  weights->push_back(newpdf1 / pdf1 * newpdf2 / pdf2);
161  }
163  }
164 }
165 
dot
T dot(const Basic3DVector &v) const
Scalar product, or "dot" product, with a vector of same type.
Definition: Basic3DVectorLD.h:212
GenEventInfoProduct
Definition: GenEventInfoProduct.h:17
electrons_cff.bool
bool
Definition: electrons_cff.py:393
mps_fire.i
i
Definition: mps_fire.py:428
genParticles2HepMC_cfi.genParticles
genParticles
Definition: genParticles2HepMC_cfi.py:4
GenEventInfoProduct::pdf
const PDF * pdf() const
Definition: GenEventInfoProduct.h:45
CompositeCandidate.h
funct::false
false
Definition: Factorize.h:29
EDProducer.h
PdfWeightProducer::genToken_
edm::EDGetTokenT< reco::GenParticleCollection > genToken_
Definition: PdfWeightProducer.cc:34
gen::PdfInfo::x
std::pair< double, double > x
Definition: PdfInfo.h:13
reco::GenParticle
Definition: GenParticle.h:21
mps_update.status
status
Definition: mps_update.py:69
edm::EDGetTokenT< reco::GenParticleCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
testProducerWithPsetDescEmpty_cfi.x2
x2
Definition: testProducerWithPsetDescEmpty_cfi.py:28
globals_cff.id1
id1
Definition: globals_cff.py:33
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
reco::GenParticleCollection
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
Definition: GenParticleFwd.h:13
PdfWeightProducer::genTag_
edm::InputTag genTag_
Definition: PdfWeightProducer.cc:33
Booster.h
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle
Definition: AssociativeIterator.h:50
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
class-composition.Q
Q
Definition: class-composition.py:82
GenParticle.h
MakerMacros.h
LHAPDF::usePDFMember
void usePDFMember(int nset, int member)
testProducerWithPsetDescEmpty_cfi.x1
x1
Definition: testProducerWithPsetDescEmpty_cfi.py:33
part
part
Definition: HCALResponse.h:20
LHAPDF::initPDFSet
void initPDFSet(int nset, const std::string &filename, int member=0)
PdfWeightProducer::beginJob
void beginJob() override
Definition: PdfWeightProducer.cc:87
LHAPDF::getQ2max
double getQ2max(int nset, int member)
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
HLT_FULL_cff.weights
weights
Definition: HLT_FULL_cff.py:99178
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
PdfWeightProducer::pdfSetNames_
std::vector< std::string > pdfSetNames_
Definition: PdfWeightProducer.cc:37
gen::PdfInfo::xPDF
std::pair< double, double > xPDF
Definition: PdfInfo.h:14
LHAPDF::extrapolate
void extrapolate(bool extrapolate=true)
corrVsCorr.filename
filename
Definition: corrVsCorr.py:123
LHAPDF::xfx
double xfx(int nset, double x, double Q, int fl)
PdfWeightProducer::PdfWeightProducer
PdfWeightProducer(const edm::ParameterSet &)
Definition: PdfWeightProducer.cc:54
dqmdumpme.k
k
Definition: dqmdumpme.py:60
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
PdfWeightProducer::pdfInfoTag_
edm::InputTag pdfInfoTag_
Definition: PdfWeightProducer.cc:35
LHAPDF::getXmin
double getXmin(int nset, int member)
edm::ParameterSet
Definition: ParameterSet.h:47
GenEventInfoProduct.h
Event.h
PdfWeightProducer::pdfInfoToken_
edm::EDGetTokenT< GenEventInfoProduct > pdfInfoToken_
Definition: PdfWeightProducer.cc:36
AddFourMomenta.h
iEvent
int iEvent
Definition: GenABIO.cc:224
LHAPDF::getQ2min
double getQ2min(int nset, int member)
trackerHitRTTI::vector
Definition: trackerHitRTTI.h:21
edm::InputTag::encode
std::string encode() const
Definition: InputTag.cc:159
edm::EventSetup
Definition: EventSetup.h:57
edm::LogError
Log< level::Error, false > LogError
Definition: MessageLogger.h:123
gen::PdfInfo::scalePDF
double scalePDF
Definition: PdfInfo.h:15
PdfWeightProducer::useFirstAsDefault_
bool useFirstAsDefault_
Definition: PdfWeightProducer.cc:32
PdfWeightProducer::~PdfWeightProducer
~PdfWeightProducer() override
Definition: PdfWeightProducer.cc:84
LHAPDF::getXmax
double getXmax(int nset, int member)
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
Frameworkfwd.h
LHAPDF::numberPDF
int numberPDF(int nset)
gen::PdfInfo::id
std::pair< int, int > id
Definition: PdfInfo.h:12
PdfWeightProducer::fixPOWHEG_
std::string fixPOWHEG_
Definition: PdfWeightProducer.cc:31
edm::EDProducer
Definition: EDProducer.h:35
PdfWeightProducer::endJob
void endJob() override
Definition: PdfWeightProducer.cc:94
PdfWeightProducer::produce
void produce(edm::Event &, const edm::EventSetup &) override
Definition: PdfWeightProducer.cc:97
PdfWeightProducer
Definition: PdfWeightProducer.cc:21
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
LHAPDF
Definition: PdfWeightProducer.cc:41
ParameterSet.h
PdfWeightProducer::pdfShortNames_
std::vector< std::string > pdfShortNames_
Definition: PdfWeightProducer.cc:38
globals_cff.id2
id2
Definition: globals_cff.py:34
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
muonDTDigis_cfi.pset
pset
Definition: muonDTDigis_cfi.py:27