17 #include <type_traits> 23 template <
typename _MatrixType>
54 std::uniform_int_distribution<uint32_t>
m_uint32;
55 std::uniform_int_distribution<uint64_t>
m_uint64;
99 template <
typename T,
typename U,
typename V>
101 T& ret_number,
const U min_val,
const V max_val)
103 const T
range = max_val - min_val + 1;
106 ret_number = min_val + (rnd %
range);
115 2.3283064370807973754314699618685e-10;
125 MAT&
matrix,
const double unif_min = 0,
const double unif_max = 1)
127 for (
size_t r = 0;
r <
matrix.rows();
r++)
128 for (
size_t c = 0;
c <
matrix.cols();
c++)
138 VEC&
v,
const double unif_min = 0,
const double unif_max = 1)
140 const size_t N =
v.size();
141 for (
size_t c = 0;
c < N;
c++)
187 template <
class MATRIX,
class AUXVECTOR_T = MATRIX>
189 const size_t dim,
const double std_scale = 1.0,
190 const double diagonal_epsilon = 1e-8)
192 AUXVECTOR_T
r(dim, 1);
195 cov.resize(dim, dim);
197 for (
size_t i = 0; i < dim; i++)
198 cov(i, i) += diagonal_epsilon;
208 VEC&
v,
const double mean = 0,
const double std = 1)
210 const size_t N =
v.size();
211 for (
size_t c = 0;
c < N;
c++)
212 v[
c] =
static_cast<std::remove_reference_t<decltype(v[c])>
>(
222 template <
typename T,
typename MATRIX>
224 std::vector<T>& out_result,
const MATRIX&
cov,
225 const std::vector<T>*
mean =
nullptr)
227 const size_t dim =
cov.cols();
228 if (
cov.rows() !=
cov.cols())
229 throw std::runtime_error(
230 "drawGaussianMultivariate(): cov is not square.");
232 throw std::runtime_error(
233 "drawGaussianMultivariate(): mean and cov sizes ");
237 out_result.resize(dim, 0);
243 cov.eigenVectors(Z, D);
245 D = D.array().sqrt().matrix();
247 for (
size_t i = 0; i < dim; i++)
250 for (
size_t d = 0; d < dim; d++)
251 out_result[d] += (Z.get_unsafe(d, i) * rnd);
254 for (
size_t d = 0; d < dim; d++) out_result[d] += (*
mean)[d];
263 template <
class VECTORLIKE,
class COVMATRIX>
265 VECTORLIKE& out_result,
const COVMATRIX&
cov,
266 const VECTORLIKE*
mean =
nullptr)
268 const size_t N =
cov.rows();
269 if (
cov.rows() !=
cov.cols())
270 throw std::runtime_error(
271 "drawGaussianMultivariate(): cov is not square.");
272 if (
mean &&
size_t(
mean->size()) != N)
273 throw std::runtime_error(
274 "drawGaussianMultivariate(): mean and cov sizes ");
280 auto eigVecs = eigensolver.eigenvectors();
281 auto eigVals = eigensolver.eigenvalues();
285 eigVals = eigVals.array().sqrt();
286 for (
typename COVMATRIX::Index i = 0; i < eigVecs.cols(); i++)
287 eigVecs.col(i) *= eigVals[i];
290 out_result.assign(N, 0);
292 for (
size_t i = 0; i < N; i++)
295 for (
size_t d = 0; d < N; d++)
296 out_result[d] += eigVecs.coeff(d, i) * rnd;
299 for (
size_t d = 0; d < N; d++) out_result[d] += (*
mean)[d];
309 template <
typename VECTOR_OF_VECTORS,
typename COVMATRIX>
311 VECTOR_OF_VECTORS& ret,
size_t desiredSamples,
const COVMATRIX&
cov,
312 const typename VECTOR_OF_VECTORS::value_type*
mean =
nullptr)
314 const size_t N =
cov.rows();
315 if (
cov.rows() !=
cov.cols())
316 throw std::runtime_error(
317 "drawGaussianMultivariateMany(): cov is not square.");
318 if (
mean &&
size_t(
mean->size()) != N)
319 throw std::runtime_error(
320 "drawGaussianMultivariateMany(): mean and cov sizes ");
327 typename COVMATRIX::PlainObject>::MatrixType eigVecs =
328 eigensolver.eigenvectors();
330 typename COVMATRIX::PlainObject>::RealVectorType eigVals =
331 eigensolver.eigenvalues();
335 eigVals = eigVals.array().sqrt();
336 for (
typename COVMATRIX::Index i = 0; i < eigVecs.cols(); i++)
337 eigVecs.col(i) *= eigVals[i];
340 ret.resize(desiredSamples);
341 for (
size_t k = 0; k < desiredSamples; k++)
344 for (
size_t i = 0; i < N; i++)
347 for (
size_t d = 0; d < N; d++)
348 ret[k][d] += eigVecs.coeff(d, i) * rnd;
351 for (
size_t d = 0; d < N; d++) ret[k][d] += (*
mean)[d];
366 out_result = in_vector;
367 const size_t N = out_result.size();
395 MAT&
matrix,
const double unif_min = 0,
const double unif_max = 1)
397 for (
size_t r = 0;
r <
matrix.rows();
r++)
398 for (
size_t c = 0;
c <
matrix.cols();
c++)
408 std::vector<T>& v_out,
const T& unif_min = 0,
const T& unif_max = 1)
410 size_t n = v_out.size();
411 for (
size_t r = 0;
r <
n;
r++)
424 for (
size_t r = 0;
r <
matrix.rows();
r++)
425 for (
size_t c = 0;
c <
matrix.cols();
c++)
435 std::vector<T>& v_out,
const T&
mean = 0,
const T&
std = 1)
437 size_t n = v_out.size();
438 for (
size_t r = 0;
r <
n;
r++)
456 const std::vector<T>& in_vector, std::vector<T>& out_result)
474 template <
typename T,
typename MATRIX>
476 const MATRIX&
cov,
size_t desiredSamples, std::vector<std::vector<T>>& ret,
477 std::vector<T>* samplesLikelihoods =
nullptr)
480 ret, desiredSamples,
cov,
static_cast<const std::vector<T>*
>(
nullptr),
489 template <
typename T,
typename MATRIXLIKE>
491 const MATRIXLIKE&
cov,
size_t desiredSamples,
492 std::vector<std::vector<T>>& ret)
502 template <
typename T,
typename MATRIX>
uint32_t drawUniform32bit()
Generate a uniformly distributed pseudo-random number using the MT19937 algorithm, in the whole range of 32-bit integers.
double drawUniform(const double Min, const double Max)
Generate a uniformly distributed pseudo-random number using the MT19937 algorithm, scaled to the selected range.
void permuteVector(const VEC &in_vector, VEC &out_result)
Returns a random permutation of a vector: all the elements of the input vector are in the output but ...
MATRIX drawDefinitePositiveMatrix(const size_t dim, const double std_scale=1.0, const double diagonal_epsilon=1e-8)
Generates a random definite-positive matrix of the given size, using the formula C = v*v^t + epsilon*...
void drawUniformUnsignedInt(uint32_t &ret_number)
You can call this overloaded method with either 32 or 64bit unsigned ints for the sake of general cod...
void shuffle(RandomIt first, RandomIt last, URBG &&g)
Uniform shuffle a sequence.
void drawGaussian1DVector(VEC &v, const double mean=0, const double std=1)
Fills the given vector with independent, 1D-normally distributed samples.
void randomize(const uint32_t seed)
Initialize the PRNG from the given random seed.
void drawGaussian1DMatrix(MAT &matrix, const double mean=0, const double std=1)
Fills the given matrix with independent, 1D-normally distributed samples.
std::uniform_int_distribution< uint64_t > m_uint64
std::normal_distribution< double > m_normdistribution
A thred-safe pseudo random number generator, based on an internal MT19937 randomness generator...
void randomPermutation(const std::vector< T > &in_vector, std::vector< T > &out_result)
Returns a random permutation of a vector: all the elements of the input vector are in the output but ...
void drawUniformUnsignedIntRange(T &ret_number, const U min_val, const V max_val)
Return a uniform unsigned integer in the range [min_val,max_val] (both inclusive) ...
double drawGaussian1D(const double mean, const double std)
Generate a normally distributed pseudo-random number.
void drawUniformMatrix(MAT &matrix, const double unif_min=0, const double unif_max=1)
Fills the given matrix with independent, uniformly distributed samples.
void drawGaussianMultivariate(VECTORLIKE &out_result, const COVMATRIX &cov, const VECTORLIKE *mean=nullptr)
Generate multidimensional random samples according to a given covariance matrix.
void drawGaussianMultivariateMany(VECTOR_OF_VECTORS &ret, size_t desiredSamples, const COVMATRIX &cov, const typename VECTOR_OF_VECTORS::value_type *mean=nullptr)
Generate a given number of multidimensional random samples according to a given covariance matrix...
void Randomize(const uint32_t seed)
Randomize the generators.
std::mt19937_64 m_MT19937
Data used internally by the MT19937 PRNG algorithm.
void vectorRandomNormal(std::vector< T > &v_out, const T &mean=0, const T &std=1)
Generates a random vector with independent, normally distributed samples.
CRandomGenerator(const uint32_t seed)
Constructor for providing a custom random seed to initialize the PRNG.
ptrdiff_t random_generator_for_STL(ptrdiff_t i)
A random number generator for usage in STL algorithms expecting a function like this (eg...
uint64_t drawUniform64bit()
Returns a uniformly distributed pseudo-random number by joining two 32bit numbers from drawUniform32b...
unsigned __int64 uint64_t
This is the global namespace for all Mobile Robot Programming Toolkit (MRPT) libraries.
GLdouble GLdouble GLdouble r
void drawGaussianMultivariate(std::vector< T > &out_result, const MATRIX &cov, const std::vector< T > *mean=nullptr)
Generate multidimensional random samples according to a given covariance matrix.
void matrixRandomNormal(MAT &matrix, const double mean=0, const double std=1)
Fills the given matrix with independent, normally distributed samples.
void matrixRandomUni(MAT &matrix, const double unif_min=0, const double unif_max=1)
Fills the given matrix with independent, uniformly distributed samples.
void randomize()
Randomize the generators, based on std::random_device.
void MT19937_initializeGenerator(const uint32_t &seed)
std::uniform_int_distribution< uint32_t > m_uint32
CRandomGenerator()
Default constructor: initialize random seed based on current time.
void vectorRandomUni(std::vector< T > &v_out, const T &unif_min=0, const T &unif_max=1)
Fills the given matrix with independent, uniformly distributed samples.
Eigen::Matrix< typename MATRIX::Scalar, MATRIX::ColsAtCompileTime, MATRIX::ColsAtCompileTime > cov(const MATRIX &v)
Computes the covariance matrix from a list of samples in an NxM matrix, where each row is a sample...
void randomNormalMultiDimensionalMany(const MATRIX &cov, size_t desiredSamples, std::vector< std::vector< T >> &ret, std::vector< T > *samplesLikelihoods=nullptr)
Generate a given number of multidimensional random samples according to a given covariance matrix...
void randomNormalMultiDimensional(const MATRIX &cov, std::vector< T > &out_result)
Generate multidimensional random samples according to a given covariance matrix.
unsigned __int32 uint32_t
CRandomGenerator & getRandomGenerator()
A static instance of a CRandomGenerator class, for use in single-thread applications.
void drawUniformVector(VEC &v, const double unif_min=0, const double unif_max=1)
Fills the given vector with independent, uniformly distributed samples.
GLuint GLuint GLsizei GLenum type
EIGEN_STRONG_INLINE double mean() const
Computes the mean of the entire matrix.
void drawUniformUnsignedInt(uint64_t &ret_number)
double drawGaussian1D_normalized()
Generate a normalized (mean=0, std=1) normally distributed sample.