Belle II Software development
MultivariateNormalGenerator.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11#include <Eigen/Dense>
12#include <Math/Vector3D.h>
13#include <TRandom.h>
14#include <TVectorT.h>
15#include <TVector3.h>
16#include <TMatrixTBase.h>
17
18namespace Belle2 {
54 public:
64 MultivariateNormalGenerator(int n, const double* mean, const double* cov)
65 {
66 setMeanCov(n, mean, cov);
67 }
72 MultivariateNormalGenerator(const Eigen::VectorXd& mean, const Eigen::MatrixXd& cov)
73 {
74 setMeanCov(mean, cov);
75 }
80 Eigen::VectorXd generate() const
81 {
82 //To get the correlated multivariate normal distribution, we multiply
83 //standard normal distributed values (mean=0, sigma=1) with a
84 //transformation matrix we obtained from an LDLT decomposition and add
85 //the mean values.
86 Eigen::VectorXd x(m_mean.rows());
87 for (int i = 0; i < m_mean.rows(); ++i) x(i) = gRandom->Gaus();
88 return m_mean + (m_transform * x);
89 }
93 void reset();
95 size_t size() const { return m_mean.rows(); }
100 void generate(double* output) const
101 {
102 Eigen::VectorXd x = generate();
103 for (int i = 0; i < x.rows(); ++i) { output[i] = x(i); }
104 }
105
112 TVector3 generateVec3() const
113 {
114 Eigen::VectorXd x = generate();
115 TVector3 output(0, 0, 0);
116 for (unsigned int i = 0; i < std::min(3u, (unsigned int)size()); ++i) {
117 output[i] = x(i);
118 }
119 return output;
120 }
121
125 TVectorD generateVecT() const
126 {
127 Eigen::VectorXd x = generate();
128 TVectorD output(x.rows());
129 output.SetElements(x.data());
130 return output;
131 }
132
142 bool setMeanCov(int n, const double* mean, const double* cov);
143
149 bool setMeanCov(const Eigen::VectorXd& mean, const Eigen::MatrixXd& cov);
150
157 template<class value_type> bool setMeanCov(const TVectorT<value_type>& mean,
158 const TMatrixTBase<value_type>& cov);
159
165 template<class value_type> bool setMeanCov(const ROOT::Math::XYZVector& mean,
166 const TMatrixTBase<value_type>& cov);
167
168 private:
170 Eigen::VectorXd m_mean;
173 Eigen::MatrixXd m_transform;
174 };
175
176 template<class value_type> bool MultivariateNormalGenerator::setMeanCov(
177 const TVectorT<value_type>& mean, const TMatrixTBase<value_type>& cov)
178 {
179 Eigen::VectorXd emean(mean.GetNrows());
180 Eigen::MatrixXd ecov(cov.GetNrows(), cov.GetNcols());
181 for (int i = 0; i < mean.GetNrows(); ++i) { emean(i) = mean(i); }
182 for (int i = 0; i < cov.GetNrows(); ++i) {
183 for (int j = 0; j < cov.GetNcols(); ++j) {
184 ecov(i, j) = cov(i, j);
185 }
186 }
187 return setMeanCov(emean, ecov);
188 }
189
190 template<class value_type> bool MultivariateNormalGenerator::setMeanCov(
191 const ROOT::Math::XYZVector& mean, const TMatrixTBase<value_type>& cov)
192 {
193 TVectorT<value_type> tmean(3);
194 tmean[0] = mean.X();
195 tmean[1] = mean.Y();
196 tmean[2] = mean.Z();
197 return setMeanCov(tmean, cov);
198 }
200}
Class to generate normal distributed, correlated random numbers given the mean values and the covaria...
size_t size() const
Return the number of elements to be generated on generate()
MultivariateNormalGenerator(const Eigen::VectorXd &mean, const Eigen::MatrixXd &cov)
constructor with Eigen matrix interface.
void generate(double *output) const
Generate a set of correlated random numbers with the previouly set mean and covariance and store them...
MultivariateNormalGenerator(int n, const double *mean, const double *cov)
constructor with array interface: mean and covariance are passed as double arrays where the covarianc...
Eigen::VectorXd m_mean
Member to store the mean values of the distribution.
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...
Eigen::VectorXd generate() const
Generate a set of correlated random numbers with the previouly set mean and covariance.
MultivariateNormalGenerator()
default constructor to allow later initialization
void reset()
reset the generator setting the size to 0.
TVectorD generateVecT() const
Generate a set of correlated random numbers with the previouly set mean and covariance and return a T...
Eigen::MatrixXd m_transform
Member to store the transformation matrix for standard normal distributed random values.
TVector3 generateVec3() const
Generate a set of correlated random numbers with the previouly set mean and covariance and return a T...
Abstract base class for different kinds of events.