All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Typedefs | Functions
sbn::evwgh Namespace Reference

Classes

class  SBNEventWeight
 
class  WeightCalc
 
class  WeightCalcCreator
 
class  WeightCalcImpl
 
class  WeightCalcFactory
 
class  WeightManager
 
class  FluxWeightCalc
 
class  GenieWeightCalc
 
struct  EventWeightParameter
 A single parameter to be reweighted. More...
 
class  EventWeightParameterSet
 Container for a set of reweightable parameters. More...
 

Typedefs

typedef std::map< std::string,
std::vector< float > > 
EventWeightMap
 Container for event-level weights. More...
 

Functions

std::ostream & operator<< (std::ostream &os, const sbn::evwgh::EventWeightParameterSet &p)
 
std::vector< std::vector
< double > > 
MultiGaussianSmearing (std::vector< double > const &centralValue, std::vector< std::vector< double > > const &inputCovarianceMatrix, int n_multisims, CLHEP::RandGaussQ &GaussRandom)
 
std::vector< double > MultiGaussianSmearing (std::vector< double > const &centralValue, TMatrixD *const &inputCovarianceMatrix, std::vector< float > rand)
 
std::vector< double > MultiGaussianSmearing (std::vector< double > const &centralValue, TMatrixD *const &LowerTriangleCovarianceMatrix, bool isDecomposed, std::vector< float > rand)
 

Typedef Documentation

Container for event-level weights.

Provides a mapping from a string identifier for a particular weight calculator to the corresponding set of weights for each universe.

Definition at line 18 of file EventWeightMap.h.

Function Documentation

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

Apply Gaussian smearing to a set of data.

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

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

Definition at line 10 of file SmearingUtils.cxx.

14  {
15  std::vector<std::vector<double> > setOfSmearedCentralValues;
16 
17  // Check that covarianceMatrix is of compatible dimension with the central values
18  unsigned int covarianceMatrix_dim = centralValue.size();
19 
20  if (inputCovarianceMatrix.size() != covarianceMatrix_dim) {
21  throw art::Exception(art::errors::Configuration)
22  << "inputCovarianceMatrix rows " << inputCovarianceMatrix.size()
23  << " not equal to entries of centralValue[]: " << covarianceMatrix_dim;
24  }
25 
26  for (auto const & row : inputCovarianceMatrix) {
27  if (row.size() != covarianceMatrix_dim) {
28  throw art::Exception(art::errors::Configuration)
29  << "inputCovarianceMatrix columns " << row.size()
30  << " not equal to entries of centralValue[]: " << covarianceMatrix_dim;
31  }
32  }
33 
34  // Change covariance matrix into a TMatrixD object
35  int dim = int(inputCovarianceMatrix.size());
36  TMatrixD covarianceMatrix(dim, dim);
37  for (size_t i=0; i<inputCovarianceMatrix.size(); i++) {
38  for (size_t j=0; j<inputCovarianceMatrix[i].size(); j++) {
39  covarianceMatrix[i][j] = inputCovarianceMatrix[i][j];
40  }
41  }
42 
43  // Perform Choleskey Decomposition
44  TDecompChol dc = TDecompChol(covarianceMatrix);
45  if(!dc.Decompose()) {
46  throw art::Exception(art::errors::StdException)
47  << "Cannot decompose covariance matrix to begin smearing.";
48  return setOfSmearedCentralValues;
49  }
50 
51  // Get upper triangular matrix. This maintains the relations in the
52  // covariance matrix, but simplifies the structure.
53  TMatrixD U = dc.GetU();
54 
55  // for every multisim
56  for (int n=0; n<n_multisims; n++) {
57  // Get a gaussian random number for every central value element
58  int dim = centralValue.size();
59 
60  std::vector<double> rands(dim);
61  GaussRandom.fireArray(dim, rands.data());
62 
63  // Compute the smeared central values
64  std::vector<double> smearedCentralValues;
65  for(int col=0; col<dim; col++) {
66  // Find the weight of each col of the upper triangular cov. matrix
67  double weightFromU = 0.;
68  for (int row=0; row<col+1; row++) {
69  weightFromU += U(row,col) * rands[row];
70  }
71 
72  // multiply this weight by each col of the central values to obtain
73  // the gaussian smeared and constrained central values
74  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
75  smearedCentralValues.push_back(weightFromU + centralValue[col]);
76  }
77 
78  // collect the smeared central values into a set
79  setOfSmearedCentralValues.push_back(smearedCentralValues);
80  }
81 
82  return setOfSmearedCentralValues;
83 }
std::vector< double > sbn::evwgh::MultiGaussianSmearing ( std::vector< double > const &  centralValue,
TMatrixD *const &  inputCovarianceMatrix,
std::vector< float >  rand 
)

Definition at line 86 of file SmearingUtils.cxx.

89  {
90  std::vector<double> smearedCentralValues;
91 
92  // Perform Choleskey Decomposition
93  // see http://pdg.lbl.gov/2016/reviews/rpp2016-rev-monte-carlo-techniques.pdf (Page 5)
94  TDecompChol dc = TDecompChol(*(inputCovarianceMatrix));
95  if (!dc.Decompose()) {
96  throw art::Exception(art::errors::StdException)
97  << "Cannot decompose covariance matrix to begin smearing.";
98  return smearedCentralValues;
99  }
100 
101  // Get upper triangular matrix. This maintains the relations in the
102  // covariance matrix, but simplifies the structure.
103  TMatrixD U = dc.GetU();
104 
105  for (size_t col=0; col<centralValue.size(); col++) {
106  // Find the weight of each col of the upper triangular cov. matrix
107  double weightFromU = 0;
108  for (unsigned int row=0; row<col+1; row++) {
109  weightFromU += U(row,col) * rand[row];
110  }
111 
112  // multiply this weight by each col of the central values to obtain
113  // the gaussian smeared and constrained central values
114  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
115  smearedCentralValues.push_back(weightFromU + centralValue[col]);
116  }
117 
118  return smearedCentralValues;
119 }
std::vector< double > sbn::evwgh::MultiGaussianSmearing ( std::vector< double > const &  centralValue,
TMatrixD *const &  LowerTriangleCovarianceMatrix,
bool  isDecomposed,
std::vector< float >  rand 
)

Definition at line 122 of file SmearingUtils.cxx.

126  {
127  std::vector<double> smearedCentralValues;
128 
129  if (!isDecomposed) {
130  throw art::Exception(art::errors::StdException)
131  << "Must supply the decomposed lower triangular covariance matrix.";
132  return smearedCentralValues;
133  }
134 
135  for(unsigned int col = 0; col < centralValue.size(); ++col) {
136  // Find the weight of each col of the upper triangular cov. matrix
137  double weightFromU = 0;
138  for(size_t row=0; row<col+1; row++) {
139  weightFromU += LowerTriangleCovarianceMatrix[0][row][col] * rand[row];
140  }
141 
142  // multiply this weight by each col of the central values to obtain
143  // the gaussian smeared and constrained central values
144  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
145  smearedCentralValues.push_back(weightFromU + centralValue[col]);
146  }
147 
148  return smearedCentralValues;
149 }
std::ostream& sbn::evwgh::operator<< ( std::ostream &  os,
const sbn::evwgh::EventWeightParameterSet p 
)

Definition at line 132 of file CAFMaker_module.cc.

133  {
134  // TODO proper implementation of this should be added in sbnobj
135  os << p.fName << " " << p.fRWType << std::endl;
136  for(const auto& it: p.fParameterMap){
137  os << it.first.fName << " " << it.first.fMean << " " << it.first.fWidth << std::endl << " ";
138  for(float v: it.second) os << " " << v;
139  os << std::endl;
140  }
141  return os;
142  }
ReweightType fRWType
Type of throws (the same for all parameters in a set)
std::map< EventWeightParameter, std::vector< float > > fParameterMap
Mapping of definitions to the set of values.
std::string fName
Name of the parameter set.