LArSoft  v09_90_00
Liquid Argon Software toolkit - https://larsoft.org/
evwgh::WeightCalc Class Referenceabstract

#include "WeightCalc.h"

Inheritance diagram for evwgh::WeightCalc:
evwgh::GenieWeightCalc evwgh::PPFXCVWeightCalc evwgh::PPFXMIPPKaonWeightCalc evwgh::PPFXMIPPPionWeightCalc evwgh::PPFXOtherWeightCalc evwgh::PPFXTargAttenWeightCalc evwgh::PPFXThinKaonWeightCalc evwgh::PPFXThinMesonWeightCalc evwgh::PPFXThinNeutronPionWeightCalc evwgh::PPFXThinNucAWeightCalc evwgh::PPFXThinNucWeightCalc evwgh::PPFXThinPionWeightCalc evwgh::PPFXTotAbsorpWeightCalc evwgh::PPFXWeightCalc

Public Member Functions

virtual void Configure (fhicl::ParameterSet const &pset, CLHEP::HepRandomEngine &)=0
 
virtual std::vector< std::vector< double > > GetWeight (art::Event &e)=0
 
void SetName (std::string name)
 
std::string GetName ()
 

Static Public Member Functions

static std::vector< std::vector< double > > MultiGaussianSmearing (std::vector< double > const &centralValues, std::vector< std::vector< double >> const &inputCovarianceMatrix, int n_multisims, CLHEP::RandGaussQ &GaussRandom)
 Applies Gaussian smearing to a set of data. More...
 
static std::vector< double > MultiGaussianSmearing (std::vector< double > const &centralValue, TMatrixD *const &inputCovarianceMatrix, std::vector< double > rand)
 Over load of the above function that only returns a single varied parameter set. More...
 
static std::vector< double > MultiGaussianSmearing (std::vector< double > const &centralValue, TMatrixD *const &LowerTriangleCovarianceMatrix, bool isDecomposed, std::vector< double > rand)
 

Private Attributes

std::string fName
 

Detailed Description

Definition at line 21 of file WeightCalc.h.

Member Function Documentation

std::vector< std::vector< double > > evwgh::WeightCalc::MultiGaussianSmearing ( std::vector< double > const &  centralValues,
std::vector< std::vector< double >> const &  inputCovarianceMatrix,
int  n_multisims,
CLHEP::RandGaussQ &  GaussRandom 
)
static

Applies Gaussian smearing to a set of data.

Parameters
centralValuesthe values to be smeared
inputCovarianceMatrixcovariance matrix for smearing
n_multisimsnumber of sets of smeared values to be produced
Returns
a set of n_multisims value sets smeared from the central value

If centralValues is of dimension N, inputCovarianceMatrix needs to be NxN, and each of the returned data sets will be also of dimension N.

Definition at line 14 of file WeightCalc.cxx.

References col, art::errors::Configuration, n, and art::errors::StdException.

19  {
20 
21  std::vector<std::vector<double>> setOfSmearedCentralValues;
22 
23  //Check that covarianceMatrix is of compatible dimension with the central values
24  unsigned int covarianceMatrix_dim = centralValue.size();
25 
26  if (inputCovarianceMatrix.size() != covarianceMatrix_dim) {
28  << "inputCovarianceMatrix rows " << inputCovarianceMatrix.size()
29  << " not equal to entries of centralValue[]: " << covarianceMatrix_dim;
30  }
31 
32  for (auto const& row : inputCovarianceMatrix) { // Range-for (C++11).
33  if (row.size() != covarianceMatrix_dim) {
35  << "inputCovarianceMatrix columns " << row.size()
36  << " not equal to entries of centralValue[]: " << covarianceMatrix_dim;
37  }
38  }
39 
40  //change covariance matrix into a TMartixD object
41  int dim = int(inputCovarianceMatrix.size());
42  TMatrixD covarianceMatrix(dim, dim);
43  int i = 0;
44  for (auto const& row : inputCovarianceMatrix) {
45  int j = 0;
46  for (auto const& element : row) {
47  covarianceMatrix[i][j] = element;
48  ++j;
49  }
50  ++i;
51  }
52 
53  //perform Choleskey Decomposition
54  TDecompChol dc = TDecompChol(covarianceMatrix);
55  if (!dc.Decompose()) {
57  << "Cannot decompose covariance matrix to begin smearing.";
58  return setOfSmearedCentralValues;
59  }
60 
61  //Get upper triangular matrix. This maintains the relations in the
62  //covariance matrix, but simplifies the structure.
63  TMatrixD U = dc.GetU();
64 
65  //for every multisim
66  for (int n = 0; n < n_multisims; ++n) {
67 
68  //get a gaussian random number for every central value element
69  int dim = centralValue.size();
70 
71  std::vector<double> rands(dim);
72  GaussRandom.fireArray(dim, rands.data());
73 
74  //compute the smeared central values
75  std::vector<double> smearedCentralValues;
76  for (int col = 0; col < dim; ++col) {
77  //find the weight of each col of the upper triangular cov. matrix
78  double weightFromU = 0.;
79  for (int row = 0; row < col + 1; ++row) {
80  weightFromU += U(row, col) * rands[row];
81  }
82  //multiply this weight by each col of the central values to obtain
83  //the gaussian smeared and constrained central values
84  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
85  smearedCentralValues.push_back(weightFromU + centralValue[col]);
86  }
87 
88  //collect the smeared central values into a set
89  setOfSmearedCentralValues.push_back(smearedCentralValues);
90  } //loop over multisims
91  return setOfSmearedCentralValues;
92  } // WeightCalc::MultiGaussianSmearing() - gives a vector of sets of smeared parameters
Int_t col[ntarg]
Definition: Style.C:29
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
Char_t n[5]
std::vector< double > evwgh::WeightCalc::MultiGaussianSmearing ( std::vector< double > const &  centralValue,
TMatrixD *const &  inputCovarianceMatrix,
std::vector< double >  rand 
)
static

Over load of the above function that only returns a single varied parameter set.

Definition at line 97 of file WeightCalc.cxx.

References col, and art::errors::StdException.

100  {
101 
102  //compute the smeared central values
103  std::vector<double> smearedCentralValues;
104 
105  //
106  //perform Choleskey Decomposition
107  //
108  // Best description I have found
109  // is in the PDG (Monte Carlo techniques, Algorithms, Gaussian distribution)
110  //
111  // http://pdg.lbl.gov/2016/reviews/rpp2016-rev-monte-carlo-techniques.pdf (Page 5)
112  //
113  //
114  TDecompChol dc = TDecompChol(*(inputCovarianceMatrix));
115  if (!dc.Decompose()) {
117  << "Cannot decompose covariance matrix to begin smearing.";
118  return smearedCentralValues;
119  }
120 
121  //Get upper triangular matrix. This maintains the relations in the
122  //covariance matrix, but simplifies the structure.
123  TMatrixD U = dc.GetU();
124 
125  for (unsigned int col = 0; col < centralValue.size(); ++col) {
126  //find the weight of each col of the upper triangular cov. matrix
127  double weightFromU = 0.;
128 
129  for (unsigned int row = 0; row < col + 1; ++row) {
130  weightFromU += U(row, col) * rand[row];
131  }
132 
133  //multiply this weight by each col of the central values to obtain
134  //the gaussian smeared and constrained central values
135  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
136  smearedCentralValues.push_back(weightFromU + centralValue[col]);
137  }
138  return smearedCentralValues;
139  } // WeightCalc::MultiGaussianSmearing() - Gives a single set of smeared parameters
Int_t col[ntarg]
Definition: Style.C:29
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
std::vector< double > evwgh::WeightCalc::MultiGaussianSmearing ( std::vector< double > const &  centralValue,
TMatrixD *const &  LowerTriangleCovarianceMatrix,
bool  isDecomposed,
std::vector< double >  rand 
)
static

Definition at line 151 of file WeightCalc.cxx.

References col, and art::errors::StdException.

156  {
157 
158  //compute the smeared central values
159  std::vector<double> smearedCentralValues;
160 
161  // This guards against accidentally
162  if (!isDecomposed) {
164  << "Must supply the decomposed lower triangular covariance matrix.";
165  return smearedCentralValues;
166  }
167 
168  for (unsigned int col = 0; col < centralValue.size(); ++col) {
169  //find the weight of each col of the upper triangular cov. matrix
170  double weightFromU = 0.;
171 
172  for (unsigned int row = 0; row < col + 1; ++row) {
173  weightFromU += LowerTriangleCovarianceMatrix[0][row][col] * rand[row];
174  }
175 
176  //multiply this weight by each col of the central values to obtain
177  //the gaussian smeared and constrained central values
178  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
179  smearedCentralValues.push_back(weightFromU + centralValue[col]);
180  }
181  return smearedCentralValues;
182  } // WeightCalc::MultiGaussianSmearing() - Gives a single set of smeared parameters
Int_t col[ntarg]
Definition: Style.C:29
cet::coded_exception< errors::ErrorCodes, ExceptionDetail::translate > Exception
Definition: Exception.h:66
void evwgh::WeightCalc::SetName ( std::string  name)
inline

Definition at line 25 of file WeightCalc.h.

25 { fName = name; }
std::string fName
Definition: WeightCalc.h:54

Member Data Documentation

std::string evwgh::WeightCalc::fName
private

Definition at line 54 of file WeightCalc.h.


The documentation for this class was generated from the following files: