All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
evwgh::WeightCalc Class Referenceabstract

#include <WeightCalc.h>

Inheritance diagram for evwgh::WeightCalc:
evwgh::GenieWeightCalc

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 19 of file larsim/larsim/EventWeight/Base/WeightCalc.h.

Member Function Documentation

virtual void evwgh::WeightCalc::Configure ( fhicl::ParameterSet const &  pset,
CLHEP::HepRandomEngine &   
)
pure virtual

Implemented in evwgh::GenieWeightCalc.

std::string evwgh::WeightCalc::GetName ( )
inline
virtual std::vector<std::vector<double> > evwgh::WeightCalc::GetWeight ( art::Event &  e)
pure virtual

Implemented in evwgh::GenieWeightCalc.

static 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.

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 102 of file WeightCalc.cxx.

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

Definition at line 159 of file WeightCalc.cxx.

160  {
161 
162  //compute the smeared central values
163  std::vector<double> smearedCentralValues;
164 
165  // This guards against accidentally
166  if(!isDecomposed){
167  throw art::Exception(art::errors::StdException)
168  << "Must supply the decomposed lower triangular covariance matrix.";
169  return smearedCentralValues;
170  }
171 
172 
173  for(unsigned int col = 0; col < centralValue.size(); ++col)
174  {
175  //find the weight of each col of the upper triangular cov. matrix
176  double weightFromU = 0.;
177 
178  for(unsigned int row = 0; row < col+1; ++row)
179  {
180  weightFromU += LowerTriangleCovarianceMatrix[0][row][col]*rand[row];
181  }
182 
183  //multiply this weight by each col of the central values to obtain
184  //the gaussian smeared and constrained central values
185  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
186  smearedCentralValues.push_back(weightFromU + centralValue[col]);
187  }
188  return smearedCentralValues;
189  } // WeightCalc::MultiGaussianSmearing() - Gives a single set of smeared parameters
void evwgh::WeightCalc::SetName ( std::string  name)
inline

Member Data Documentation

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

Definition at line 56 of file larsim/larsim/EventWeight/Base/WeightCalc.h.


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