CMS 3D CMS Logo

PDFWeightsTest.cc
Go to the documentation of this file.
16 #include "TTree.h"
17 #include <iostream>
18 using namespace std;
19 using namespace edm;
20 using namespace HepMC;
21 
23 private:
28 
29  unsigned int pdfWeightOffset_;
30  unsigned int nPdfWeights_;
31  unsigned int nPdfEigWeights_;
32 
33  TTree* tree_;
34  std::vector<float> pdfweights_;
35  std::vector<float> pdfeigweights_;
36  float weight_;
37 
38 public:
39  explicit PDFWeightsTest(const ParameterSet& cfg)
40  : srcToken_(consumes<LHEEventProduct>(edm::InputTag("externalLHEProducer"))),
41  srcTokenAlt_(consumes<LHEEventProduct>(edm::InputTag("source"))),
42  srcTokenGen_(consumes<GenEventInfoProduct>(edm::InputTag("generator"))),
43  pdfWeightOffset_(cfg.getParameter<unsigned int>("pdfWeightOffset")),
44  nPdfWeights_(cfg.getParameter<unsigned int>("nPdfWeights")),
45  nPdfEigWeights_(cfg.getParameter<unsigned int>("nPdfEigWeights")),
46  pdfweights_(nPdfWeights_),
47  pdfeigweights_(nPdfEigWeights_) {
49  tree_ = fs->make<TTree>("tree", "");
50 
51  tree_->Branch("pdfrep", &pdfweights_);
52  tree_->Branch("pdfeig", &pdfeigweights_);
53  tree_->Branch("weight", &weight_);
54 
55  edm::FileInPath mc2hessianCSV = cfg.getParameter<edm::FileInPath>("mc2hessianCSV");
56  pdfweightshelper_.Init(nPdfWeights_, nPdfEigWeights_, mc2hessianCSV);
57  }
58 
59 private:
60  void analyze(const Event& evt, const EventSetup& es) override {
62  evt.getByToken(srcToken_, lheInfo);
63 
64  if (!lheInfo.isValid()) {
65  evt.getByToken(srcTokenAlt_, lheInfo);
66  }
67 
68  double nomlheweight = lheInfo->weights()[0].wgt;
69 
71  evt.getByToken(srcTokenGen_, genInfo);
72 
73  weight_ = genInfo->weight();
74 
75  //get the original mc replica weights
76  std::vector<double> inpdfweights(nPdfWeights_);
77  for (unsigned int ipdf = 0; ipdf < nPdfWeights_; ++ipdf) {
78  unsigned int iwgt = ipdf + pdfWeightOffset_;
79 
80  //this is the weight to be used for evaluating uncertainties with mc replica weights
81  pdfweights_[ipdf] = lheInfo->weights()[iwgt].wgt * weight_ / nomlheweight;
82 
83  //this is the raw weight to be fed to the mc2hessian convertor
84  inpdfweights[ipdf] = lheInfo->weights()[iwgt].wgt;
85  }
86 
87  std::vector<double> outpdfweights(nPdfEigWeights_);
88  //do the actual conversion, where the nominal lhe weight is needed as the reference point for the linearization
89  pdfweightshelper_.DoMC2Hessian(nomlheweight, inpdfweights.data(), outpdfweights.data());
90 
91  for (unsigned int iwgt = 0; iwgt < nPdfEigWeights_; ++iwgt) {
92  double wgtval = outpdfweights[iwgt];
93 
94  //the is the weight to be used for evaluating uncertainties with hessian weights
95  pdfeigweights_[iwgt] = wgtval * weight_ / nomlheweight;
96  }
97 
98  tree_->Fill();
99  }
100 };
101 
103 
GenEventInfoProduct
Definition: GenEventInfoProduct.h:17
Handle.h
EDAnalyzer.h
PDFWeightsTest::srcToken_
EDGetTokenT< LHEEventProduct > srcToken_
Definition: PDFWeightsTest.cc:25
LHEEventProduct
Definition: LHEEventProduct.h:12
ESHandle.h
edm::EDGetTokenT< LHEEventProduct >
PDFWeightsTest::weight_
float weight_
Definition: PDFWeightsTest.cc:36
edm
HLT enums.
Definition: AlignableModifier.h:19
PDFWeightsTest::nPdfEigWeights_
unsigned int nPdfEigWeights_
Definition: PDFWeightsTest.cc:31
PDFWeightsTest::srcTokenAlt_
EDGetTokenT< LHEEventProduct > srcTokenAlt_
Definition: PDFWeightsTest.cc:26
edm::one::EDAnalyzer
Definition: EDAnalyzer.h:30
edm::Handle
Definition: AssociativeIterator.h:50
PDFWeightsHelper::Init
void Init(unsigned int nreplicas, unsigned int neigenvectors, const edm::FileInPath &incsv)
Definition: PDFWeightsHelper.cc:7
PDFWeightsHelper::DoMC2Hessian
void DoMC2Hessian(double nomweight, const double *inweights, double *outweights) const
Definition: PDFWeightsHelper.cc:28
EDMException.h
edm::FileInPath
Definition: FileInPath.h:64
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Service.h
PDFWeightsTest::srcTokenGen_
EDGetTokenT< GenEventInfoProduct > srcTokenGen_
Definition: PDFWeightsTest.cc:27
PDFWeightsTest::pdfeigweights_
std::vector< float > pdfeigweights_
Definition: PDFWeightsTest.cc:35
edm::Event::getByToken
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
TFileService.h
PDFWeightsTest::pdfWeightOffset_
unsigned int pdfWeightOffset_
Definition: PDFWeightsTest.cc:29
edm::ParameterSet
Definition: ParameterSet.h:36
GenEventInfoProduct.h
Event.h
ParticleDataTable.h
nano_cff.lheInfo
lheInfo
Definition: nano_cff.py:99
PDFWeightsTest::PDFWeightsTest
PDFWeightsTest(const ParameterSet &cfg)
Definition: PDFWeightsTest.cc:39
edm::Service< TFileService >
createfilelist.int
int
Definition: createfilelist.py:10
edm::EventSetup
Definition: EventSetup.h:57
PDFWeightsTest::pdfweightshelper_
PDFWeightsHelper pdfweightshelper_
Definition: PDFWeightsTest.cc:24
InputTag.h
looper.cfg
cfg
Definition: looper.py:297
PDFWeightsHelper
Definition: PDFWeightsHelper.h:10
PDFWeightsTest::pdfweights_
std::vector< float > pdfweights_
Definition: PDFWeightsTest.cc:34
LHEEventProduct.h
std
Definition: JetResolutionObject.h:76
HepMC
Definition: GenParticle.h:15
EventSetup.h
PDFWeightsTest::tree_
TTree * tree_
Definition: PDFWeightsTest.cc:33
PDFWeightsTest
Definition: PDFWeightsTest.cc:22
ParameterSet.h
HepMCProduct.h
PDFWeightsTest::analyze
void analyze(const Event &evt, const EventSetup &es) override
Definition: PDFWeightsTest.cc:60
edm::Event
Definition: Event.h:73
PDFWeightsTest::nPdfWeights_
unsigned int nPdfWeights_
Definition: PDFWeightsTest.cc:30
edm::InputTag
Definition: InputTag.h:15
PDFWeightsHelper.h
TFileService::make
T * make(const Args &... args) const
make new ROOT object
Definition: TFileService.h:64
GenEventInfoProduct::weight
double weight() const
Definition: GenEventInfoProduct.h:35