Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:25

0001 #include "PhysicsTools/HepMCCandAlgos/interface/PDFWeightsHelper.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include <fstream>
0004 
0005 PDFWeightsHelper::PDFWeightsHelper() {}
0006 
0007 void PDFWeightsHelper::Init(unsigned int nreplicas, unsigned int neigenvectors, const edm::FileInPath &incsv) {
0008   transformation_.resize(nreplicas, neigenvectors);
0009 
0010   std::ifstream instream(incsv.fullPath());
0011   if (!instream.is_open()) {
0012     throw cms::Exception("PDFWeightsHelper") << "Could not open csv file" << incsv.relativePath();
0013   }
0014 
0015   for (unsigned int ireplica = 0; ireplica < nreplicas; ++ireplica) {
0016     std::string linestr;
0017     getline(instream, linestr);
0018     std::istringstream line(linestr);
0019     for (unsigned int ieigen = 0; ieigen < neigenvectors; ++ieigen) {
0020       std::string valstr;
0021       getline(line, valstr, ',');
0022       std::istringstream val(valstr);
0023       val >> transformation_(ireplica, ieigen);
0024     }
0025   }
0026 }
0027 
0028 void PDFWeightsHelper::DoMC2Hessian(double nomweight, const double *inweights, double *outweights) const {
0029   const unsigned int nreplicas = transformation_.rows();
0030   const unsigned int neigenvectors = transformation_.cols();
0031 
0032   Eigen::VectorXd inweightv(nreplicas);
0033   for (unsigned int irep = 0; irep < nreplicas; ++irep) {
0034     inweightv[irep] = inweights[irep] - nomweight;
0035   }
0036 
0037   Eigen::VectorXd outweightv = transformation_.transpose() * inweightv;
0038 
0039   for (unsigned int ieig = 0; ieig < neigenvectors; ++ieig) {
0040     outweights[ieig] = outweightv[ieig] + nomweight;
0041   }
0042 }