Belle II Software  release-05-02-19
MultivariateNormalGenerator.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2015 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Ritter *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 #include <Eigen/Dense>
14 #include <TRandom.h>
15 #include <TVectorT.h>
16 #include <TVector3.h>
17 #include <TMatrixTBase.h>
18 
19 namespace Belle2 {
54  class MultivariateNormalGenerator {
55  public:
65  MultivariateNormalGenerator(int n, const double* mean, const double* cov)
66  {
67  setMeanCov(n, mean, cov);
68  }
73  MultivariateNormalGenerator(const Eigen::VectorXd& mean, const Eigen::MatrixXd& cov)
74  {
75  setMeanCov(mean, cov);
76  }
81  Eigen::VectorXd generate() const
82  {
83  //To get the correlated multivariate normal distribution, we multiply
84  //standard normal distributed values (mean=0, sigma=1) with a
85  //transformation matrix we obtained from an LDLT decomposition and add
86  //the mean values.
87  Eigen::VectorXd x(m_mean.rows());
88  for (int i = 0; i < m_mean.rows(); ++i) x(i) = gRandom->Gaus();
89  return m_mean + (m_transform * x);
90  }
94  void reset();
96  size_t size() const { return m_mean.rows(); }
101  void generate(double* output) const
102  {
103  Eigen::VectorXd x = generate();
104  for (int i = 0; i < x.rows(); ++i) { output[i] = x(i); }
105  }
106 
113  TVector3 generateVec3() const
114  {
115  Eigen::VectorXd x = generate();
116  TVector3 output(0, 0, 0);
117  for (unsigned int i = 0; i < std::min(3u, (unsigned int)size()); ++i) {
118  output[i] = x(i);
119  }
120  return output;
121  }
122 
126  TVectorD generateVecT() const
127  {
128  Eigen::VectorXd x = generate();
129  TVectorD output(x.rows());
130  output.SetElements(x.data());
131  return output;
132  }
133 
143  bool setMeanCov(int n, const double* mean, const double* cov);
144 
150  bool setMeanCov(const Eigen::VectorXd& mean, const Eigen::MatrixXd& cov);
151 
158  template<class value_type> bool setMeanCov(const TVectorT<value_type>& mean,
159  const TMatrixTBase<value_type>& cov);
160 
167  template<class value_type> bool setMeanCov(const TVector3& mean,
168  const TMatrixTBase<value_type>& cov);
169  private:
171  Eigen::VectorXd m_mean;
174  Eigen::MatrixXd m_transform;
175  };
176 
177  template<class value_type> bool MultivariateNormalGenerator::setMeanCov(
178  const TVectorT<value_type>& mean, const TMatrixTBase<value_type>& cov)
179  {
180  Eigen::VectorXd emean(mean.GetNrows());
181  Eigen::MatrixXd ecov(cov.GetNrows(), cov.GetNcols());
182  for (int i = 0; i < mean.GetNrows(); ++i) { emean(i) = mean(i); }
183  for (int i = 0; i < cov.GetNrows(); ++i) {
184  for (int j = 0; j < cov.GetNcols(); ++j) {
185  ecov(i, j) = cov(i, j);
186  }
187  }
188  return setMeanCov(emean, ecov);
189  }
190 
191  template<class value_type> bool MultivariateNormalGenerator::setMeanCov(
192  const TVector3& mean, const TMatrixTBase<value_type>& cov)
193  {
194  TVectorT<value_type> tmean(3);
195  for (int i = 0; i < 3; ++i) {
196  tmean[i] = mean[i];
197  }
198  return setMeanCov(tmean, cov);
199  }
201 }
Belle2::MultivariateNormalGenerator::generate
Eigen::VectorXd generate() const
Generate a set of correlated random numbers with the previouly set mean and covariance.
Definition: MultivariateNormalGenerator.h:89
Belle2::MultivariateNormalGenerator::MultivariateNormalGenerator
MultivariateNormalGenerator()
default constructor to allow later initialization
Definition: MultivariateNormalGenerator.h:65
Belle2::MultivariateNormalGenerator::generateVecT
TVectorD generateVecT() const
Generate a set of correlated random numbers with the previouly set mean and covariance and return a T...
Definition: MultivariateNormalGenerator.h:134
Belle2::MultivariateNormalGenerator::MultivariateNormalGenerator
MultivariateNormalGenerator(int n, const double *mean, const double *cov)
constructor with array interface: mean and covariance are passed as double arrays where the covarianc...
Definition: MultivariateNormalGenerator.h:73
Belle2::MultivariateNormalGenerator::m_mean
Eigen::VectorXd m_mean
Member to store the mean values of the distribution.
Definition: MultivariateNormalGenerator.h:179
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::MultivariateNormalGenerator::size
size_t size() const
Return the number of elements to be generated on generate()
Definition: MultivariateNormalGenerator.h:104
Belle2::MultivariateNormalGenerator::generateVec3
TVector3 generateVec3() const
Generate a set of correlated random numbers with the previouly set mean and covariance and return a T...
Definition: MultivariateNormalGenerator.h:121
Belle2::MultivariateNormalGenerator::m_transform
Eigen::MatrixXd m_transform
Member to store the transformation matrix for standard normal distributed random values.
Definition: MultivariateNormalGenerator.h:182
Belle2::MultivariateNormalGenerator::setMeanCov
bool setMeanCov(int n, const double *mean, const double *cov)
set the mean and covariance for the distribution with array interface: mean and covariance are passed...
Definition: MultivariateNormalGenerator.cc:22
Belle2::MultivariateNormalGenerator::reset
void reset()
reset the generator setting the size to 0.
Definition: MultivariateNormalGenerator.cc:16