All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
WeightCalc.cxx
Go to the documentation of this file.
2 
3 // art libraries
4 #include "canvas/Utilities/Exception.h"
5 
6 // ROOT libraries
7 #include "TMatrixD.h"
8 #include "TDecompChol.h"
9 
10 // CLHEP libraries
11 #include "CLHEP/Random/RandGaussQ.h"
12 
13 namespace evwgh {
14  std::vector<std::vector<double> > WeightCalc::MultiGaussianSmearing(std::vector<double> const& centralValue,std::vector< std::vector<double> > const& inputCovarianceMatrix,int n_multisims,CLHEP::RandGaussQ& GaussRandom)
15  {
16 
17  std::vector<std::vector<double> > setOfSmearedCentralValues;
18 
19  //Check that covarianceMatrix is of compatible dimension with the central values
20  unsigned int covarianceMatrix_dim = centralValue.size();
21 
22  if (inputCovarianceMatrix.size() != covarianceMatrix_dim)
23  {
24  throw art::Exception(art::errors::Configuration)
25  << "inputCovarianceMatrix rows " << inputCovarianceMatrix.size()
26  << " not equal to entries of centralValue[]: " << covarianceMatrix_dim;
27  }
28 
29  for (auto const & row : inputCovarianceMatrix )
30  { // Range-for (C++11).
31  if(row.size() != covarianceMatrix_dim)
32  {
33  throw art::Exception(art::errors::Configuration)
34  << "inputCovarianceMatrix columns " << row.size()
35  << " not equal to entries of centralValue[]: " << covarianceMatrix_dim;
36  }
37  }
38 
39  //change covariance matrix into a TMartixD object
40  int dim = int(inputCovarianceMatrix.size());
41  TMatrixD covarianceMatrix(dim,dim);
42  int i = 0;
43  for (auto const & row : inputCovarianceMatrix)
44  {
45  int j = 0;
46  for (auto const & element : row)
47  {
48  covarianceMatrix[i][j] = element;
49  ++j;
50  }
51  ++i;
52  }
53 
54  //perform Choleskey Decomposition
55  TDecompChol dc = TDecompChol(covarianceMatrix);
56  if(!dc.Decompose())
57  {
58  throw art::Exception(art::errors::StdException)
59  << "Cannot decompose covariance matrix to begin smearing.";
60  return setOfSmearedCentralValues;
61  }
62 
63  //Get upper triangular matrix. This maintains the relations in the
64  //covariance matrix, but simplifies the structure.
65  TMatrixD U = dc.GetU();
66 
67  //for every multisim
68  for(int n = 0; n < n_multisims; ++n)
69  {
70 
71  //get a gaussian random number for every central value element
72  int dim = centralValue.size();
73 
74  std::vector<double> rands(dim);
75  GaussRandom.fireArray(dim, rands.data());
76 
77  //compute the smeared central values
78  std::vector<double> smearedCentralValues;
79  for(int col = 0; col < dim; ++col)
80  {
81  //find the weight of each col of the upper triangular cov. matrix
82  double weightFromU = 0.;
83  for(int row = 0; row < col+1; ++row)
84  {
85  weightFromU += U(row,col)*rands[row];
86  }
87  //multiply this weight by each col of the central values to obtain
88  //the gaussian smeared and constrained central values
89  // smearedCentralValues.push_back(weightFromU * centralValue[col]);
90  smearedCentralValues.push_back(weightFromU + centralValue[col]);
91  }
92 
93  //collect the smeared central values into a set
94  setOfSmearedCentralValues.push_back(smearedCentralValues);
95  }//loop over multisims
96  return setOfSmearedCentralValues;
97  } // WeightCalc::MultiGaussianSmearing() - gives a vector of sets of smeared parameters
98 
99  /////////////
100  /// Over load of the above function that only returns a single varied parameter set
101  ////////////
102  std::vector<double> WeightCalc::MultiGaussianSmearing(std::vector<double> const& centralValue, TMatrixD* const& inputCovarianceMatrix, std::vector< double > rand)
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
147 
148 
149  /////////
150  //
151  /*
152  static std::vector<double> MultiGaussianSmearing(
153  std::vector<double> const& centralValue,
154  TMatrixD* const& LowerTriangleCovarianceMatrix,
155  bool isDecomposed,
156  std::vector<double> rand);
157  */
158  ////////
159  std::vector<double> WeightCalc::MultiGaussianSmearing(std::vector<double> const& centralValue, TMatrixD* const& LowerTriangleCovarianceMatrix, bool isDecomposed, std::vector< double > rand)
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
190 
191 
192 
193 } // namespace evwgh
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.
auto vector(Vector const &v)
Returns a manipulator which will print the specified array.
Definition: DumpUtils.h:265