CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Attributes
PDFWeightsHelper Class Reference

#include <PDFWeightsHelper.h>

Public Member Functions

void DoMC2Hessian (double nomweight, const double *inweights, double *outweights) const
 
void Init (unsigned int nreplicas, unsigned int neigenvectors, const edm::FileInPath &incsv)
 
unsigned int neigenvectors () const
 
 PDFWeightsHelper ()
 

Protected Attributes

Eigen::MatrixXd transformation_
 

Detailed Description

Definition at line 10 of file PDFWeightsHelper.h.

Constructor & Destructor Documentation

◆ PDFWeightsHelper()

PDFWeightsHelper::PDFWeightsHelper ( )

Definition at line 5 of file PDFWeightsHelper.cc.

5 {}

Member Function Documentation

◆ DoMC2Hessian()

void PDFWeightsHelper::DoMC2Hessian ( double  nomweight,
const double *  inweights,
double *  outweights 
) const

Definition at line 28 of file PDFWeightsHelper.cc.

References neigenvectors(), and transformation_.

Referenced by PDFWeightsTest::analyze().

28  {
29  const unsigned int nreplicas = transformation_.rows();
30  const unsigned int neigenvectors = transformation_.cols();
31 
32  Eigen::VectorXd inweightv(nreplicas);
33  for (unsigned int irep = 0; irep < nreplicas; ++irep) {
34  inweightv[irep] = inweights[irep] - nomweight;
35  }
36 
37  Eigen::VectorXd outweightv = transformation_.transpose() * inweightv;
38 
39  for (unsigned int ieig = 0; ieig < neigenvectors; ++ieig) {
40  outweights[ieig] = outweightv[ieig] + nomweight;
41  }
42 }
Eigen::VectorXd VectorXd
Definition: FitUtils.h:16
unsigned int neigenvectors() const
Eigen::MatrixXd transformation_

◆ Init()

void PDFWeightsHelper::Init ( unsigned int  nreplicas,
unsigned int  neigenvectors,
const edm::FileInPath incsv 
)

Definition at line 7 of file PDFWeightsHelper.cc.

References Exception, edm::FileInPath::fullPath(), mps_splice::line, neigenvectors(), edm::FileInPath::relativePath(), AlCaHLTBitMon_QueryRunRegistry::string, transformation_, and heppy_batch::val.

Referenced by PDFWeightsTest::PDFWeightsTest().

7  {
8  transformation_.resize(nreplicas, neigenvectors);
9 
10  std::ifstream instream(incsv.fullPath());
11  if (!instream.is_open()) {
12  throw cms::Exception("PDFWeightsHelper") << "Could not open csv file" << incsv.relativePath();
13  }
14 
15  for (unsigned int ireplica = 0; ireplica < nreplicas; ++ireplica) {
16  std::string linestr;
17  getline(instream, linestr);
18  std::istringstream line(linestr);
19  for (unsigned int ieigen = 0; ieigen < neigenvectors; ++ieigen) {
20  std::string valstr;
21  getline(line, valstr, ',');
22  std::istringstream val(valstr);
23  val >> transformation_(ireplica, ieigen);
24  }
25  }
26 }
std::string fullPath() const
Definition: FileInPath.cc:161
unsigned int neigenvectors() const
Eigen::MatrixXd transformation_
std::string relativePath() const
Definition: FileInPath.cc:157

◆ neigenvectors()

unsigned int PDFWeightsHelper::neigenvectors ( ) const
inline

Definition at line 17 of file PDFWeightsHelper.h.

References transformation_.

Referenced by DoMC2Hessian(), and Init().

17 { return transformation_.cols(); }
Eigen::MatrixXd transformation_

Member Data Documentation

◆ transformation_

Eigen::MatrixXd PDFWeightsHelper::transformation_
protected

Definition at line 20 of file PDFWeightsHelper.h.

Referenced by DoMC2Hessian(), Init(), and neigenvectors().