26 #ifndef O2SCL_PROB_DENS_FUNC_H 27 #define O2SCL_PROB_DENS_FUNC_H 29 #include <gsl/gsl_rng.h> 30 #include <gsl/gsl_randist.h> 31 #include <gsl/gsl_cdf.h> 35 #include <boost/numeric/ublas/vector.hpp> 37 #include <o2scl/hist.h> 38 #include <o2scl/rng_gsl.h> 39 #include <o2scl/search_vec.h> 40 #include <o2scl/cholesky.h> 42 #ifndef DOXYGEN_NO_O2NS 69 virtual double pdf(
double x)
const {
79 virtual double cdf(
double x)
const {
137 O2SCL_ERR2(
"Tried to create a Gaussian dist. with sigma",
138 "<0 in prob_dens_gaussian::prob_dens_gaussian().",
182 "in prob_dens_gaussian::prob_dens_gaussian().",
197 O2SCL_ERR2(
"Width not set in prob_dens_gaussian::",
206 O2SCL_ERR2(
"Width not set in prob_dens_gaussian::",
209 return cent_+gsl_ran_gaussian(&r,sigma_);
213 virtual double pdf(
double x)
const {
215 O2SCL_ERR2(
"Width not set in prob_dens_gaussian::",
218 return gsl_ran_gaussian_pdf(x-cent_,sigma_);
224 O2SCL_ERR2(
"Width not set in prob_dens_gaussian::",
227 return log(gsl_ran_gaussian_pdf(x-cent_,sigma_));
231 virtual double cdf(
double x)
const {
233 O2SCL_ERR2(
"Width not set in prob_dens_gaussian::",
236 return gsl_cdf_gaussian_P(x-cent_,sigma_);
242 O2SCL_ERR2(
"Width not set in prob_dens_gaussian::",
245 if (in_cdf<0.0 || in_cdf>1.0) {
246 O2SCL_ERR2(
"Requested cdf inverse outside of [0,1] in ",
247 "prob_dens_gaussian::invert_cdf().",
exc_einval);
249 return gsl_cdf_gaussian_Pinv(in_cdf,sigma_)+cent_;
255 O2SCL_ERR2(
"Width not set in prob_dens_gaussian::",
273 virtual double lower_limit()
const=0;
276 virtual double upper_limit()
const=0;
336 if (
this==&pdg)
return *
this;
365 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
374 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
383 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
386 return gsl_ran_flat(&r,ll,ul);
390 virtual double pdf(
double x)
const {
392 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
395 if (x<ll || x>ul)
return 0.0;
396 return gsl_ran_flat_pdf(x,ll,ul);
402 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
405 if (x<ll || x>ul)
return 0.0;
406 return log(gsl_ran_flat_pdf(x,ll,ul));
410 virtual double cdf(
double x)
const {
412 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
415 if (x<ll)
return 0.0;
416 if (x>ul)
return 1.0;
417 return gsl_cdf_flat_P(x,ll,ul);
423 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
426 if (in_cdf<0.0 || in_cdf>1.0) {
427 O2SCL_ERR2(
"Requested cdf inverse outside of [0,1] in ",
428 "prob_dens_uniform::invert_cdf().",
exc_einval);
430 return gsl_cdf_flat_Pinv(in_cdf,ll,ul);
436 O2SCL_ERR2(
"Limits not set in prob_dens_uniform::",
498 O2SCL_ERR2(
"Tried to create log normal dist. with mu or sigma",
499 "<0 in prob_dens_lognormal::prob_dens_lognormal().",
518 if (
this==&pdg)
return *
this;
528 O2SCL_ERR2(
"Tried to set mu or sigma negative",
529 "in prob_dens_lognormal::prob_dens_lognormal().",
543 return gsl_ran_lognormal(&r,mu_,sigma_);
547 virtual double pdf(
double x)
const {
551 return gsl_ran_lognormal_pdf(x,mu_,sigma_);
559 return log(gsl_ran_lognormal_pdf(x,mu_,sigma_));
563 virtual double cdf(
double x)
const {
567 return gsl_cdf_lognormal_P(x,mu_,sigma_);
572 if (in_cdf<0.0 || in_cdf>1.0) {
573 O2SCL_ERR2(
"Requested cdf inverse outside of [0,1] in ",
574 "prob_dens_lognormal::invert_cdf().",
exc_einval);
576 return gsl_cdf_lognormal_Pinv(in_cdf,mu_,sigma_);
582 O2SCL_ERR2(
"Parameters not set in prob_dens_lognormal::",
636 virtual double lower_limit()
const;
639 virtual double upper_limit()
const;
642 virtual double pdf(
double x)
const;
645 virtual double log_pdf(
double x)
const;
648 virtual double cdf(
double x)
const;
664 template<
class vec_t=boost::numeric::ublas::vector<
double> >
670 virtual size_t dim()
const {
675 virtual double pdf(
const vec_t &x)
const {
680 virtual double log_pdf(
const vec_t &x)
const {
694 template<
class vec_t=boost::numeric::ublas::vector<
double> >
700 std::vector<prob_dens_func>
list;
709 virtual size_t dim()
const {
714 virtual double pdf(
const vec_t &x)
const {
716 for(
size_t i=0;i<list.size();i++) ret*=list[i].
pdf(x[i]);
721 virtual double log_pdf(
const vec_t &x)
const {
723 for(
size_t i=0;i<list.size();i++) ret*=list[i].
pdf(x[i]);
729 for(
size_t i=0;i<list.size();i++) x[i]=list[i]();
739 template<
class vec_t=boost::numeric::ublas::vector<
double>,
740 class mat_t=boost::numeric::ublas::matrix<
double> >
772 virtual size_t dim()
const {
783 set(p_ndim,p_peak,covar);
788 void set(
size_t p_ndim, vec_t &p_peak, mat_t &covar) {
790 O2SCL_ERR(
"Zero dimension in prob_dens_mdim_gaussian::set().",
796 for(
size_t i=0;i<ndim;i++) peak[i]=p_peak[i];
806 o2scl_linalg::cholesky_invert<mat_t>(ndim,covar_inv);
810 for(
size_t i=0;i<ndim;i++) {
812 for(
size_t j=0;j<ndim;j++) {
813 if (i<j) chol(i,j)=0.0;
823 virtual double pdf(
const vec_t &x)
const {
825 O2SCL_ERR2(
"Distribution not set in prob_dens_mdim_gaussian::",
829 for(
size_t i=0;i<ndim;i++) q[i]=x[i]-peak[i];
830 vtmp=prod(covar_inv,q);
831 ret*=exp(-0.5*inner_prod(q,vtmp));
836 virtual double log_pdf(
const vec_t &x)
const {
838 O2SCL_ERR2(
"Distribution not set in prob_dens_mdim_gaussian::",
841 double ret=log(norm);
842 for(
size_t i=0;i<ndim;i++) q[i]=x[i]-peak[i];
843 vtmp=prod(covar_inv,q);
844 ret+=-0.5*inner_prod(q,vtmp);
851 O2SCL_ERR2(
"Distribution not set in prob_dens_mdim_gaussian::",
854 for(
size_t i=0;i<ndim;i++) q[i]=pdg();
856 for(
size_t i=0;i<ndim;i++) x[i]=peak[i]+vtmp[i];
866 template<
class vec_t=boost::numeric::ublas::vector<
double> >
872 virtual size_t dim()
const {
877 virtual double pdf(
const vec_t &x,
const vec_t &
x2)
const=0;
880 virtual double log_pdf(
const vec_t &x,
const vec_t &x2)
const=0;
883 virtual void operator()(
const vec_t &x, vec_t &x2)
const=0;
904 template<
class vec_t=boost::numeric::ublas::vector<
double> >
911 std::random_device
rd;
937 (vec_t &step, vec_t &low, vec_t &high) {
939 for(
size_t i=0;i<step.size();i++) {
940 u_step.push_back(step[i]);
941 if (low[i]>high[i]) {
946 u_low.push_back(low[i]);
947 u_high.push_back(high[i]);
948 d_pdf/=high[i]-low[i];
953 virtual size_t dim()
const {
954 return u_step.size();
958 virtual double pdf(
const vec_t &x,
const vec_t &
x2)
const {
963 virtual double log_pdf(
const vec_t &x,
const vec_t &
x2)
const {
969 size_t nv=u_step.size();
970 for(
size_t i=0;i<nv;i++) {
971 while (x2[i]<u_low[i] || x2[i]>u_high[i]) {
972 x2[i]=x[i]+u_step[i]*(rg.
random()*2.0-1.0);
985 template<
class vec_t=boost::numeric::ublas::vector<
double> >
998 virtual size_t dim()
const {
1003 virtual double pdf(
const vec_t &x,
const vec_t &
x2)
const {
1004 return base.
pdf(x2);
1008 virtual double log_pdf(
const vec_t &x,
const vec_t &
x2)
const {
1024 template<
class vec_t=boost::numeric::ublas::vector<
double>,
1025 class mat_t=boost::numeric::ublas::matrix<
double> >
1072 void set(
size_t p_ndim, mat_t &covar) {
1074 O2SCL_ERR(
"Zero dimension in prob_cond_mdim_gaussian::set().",
1088 o2scl_linalg::cholesky_invert<mat_t>(ndim,covar_inv);
1092 for(
size_t i=0;i<ndim;i++) {
1094 for(
size_t j=0;j<ndim;j++) {
1095 if (i<j) chol(i,j)=0.0;
1105 virtual double pdf(
const vec_t &x,
const vec_t &
x2)
const {
1107 O2SCL_ERR2(
"Distribution not set in prob_cond_mdim_gaussian::",
1111 for(
size_t i=0;i<ndim;i++) q[i]=x2[i]-x[i];
1112 vtmp=prod(covar_inv,q);
1113 ret*=exp(-0.5*inner_prod(q,vtmp));
1118 virtual double log_pdf(
const vec_t &x,
const vec_t &
x2)
const {
1120 O2SCL_ERR2(
"Distribution not set in prob_cond_mdim_gaussian::",
1123 double ret=log(norm);
1124 for(
size_t i=0;i<ndim;i++) q[i]=x2[i]-x[i];
1125 vtmp=prod(covar_inv,q);
1126 ret+=-0.5*inner_prod(q,vtmp);
1133 O2SCL_ERR2(
"Distribution not set in prob_cond_mdim_gaussian::",
1136 for(
size_t i=0;i<ndim;i++) q[i]=pdg();
1138 for(
size_t i=0;i<ndim;i++) x2[i]=x[i]+vtmp[i];
1145 #ifndef DOXYGEN_NO_O2NS vec_t peak
Location of the peak.
double sigma_
Width parameter.
virtual double entropy() const
The inverse cumulative distribution function.
virtual double log_pdf(const vec_t &x, const vec_t &x2) const
The log of the normalized density.
The main O<span style='position: relative; top: 0.3em; font-size: 0.8em'>2</span>scl O$_2$scl names...
virtual double cdf(double x) const
The cumulative distribution function (from the lower tail)
prob_cond_mdim_gaussian(size_t p_ndim, mat_t &covar)
Create a distribution from the covariance matrix.
double mean()
Get the center.
virtual double log_pdf(const vec_t &x) const
The log of the normalized density.
virtual double invert_cdf(double in_cdf) const
The inverse cumulative distribution function.
virtual size_t dim() const
The dimensionality.
std::vector< double > u_step
Desc.
virtual double pdf(const vec_t &x) const
The normalized density.
virtual size_t dim() const
The dimensionality.
void set_seed(unsigned long int s)
Set the seed.
prob_dens_gaussian & operator=(const prob_dens_gaussian &pdg)
Copy constructor with operator=.
o2scl::prob_dens_gaussian pdg
Standard normal.
virtual double pdf(const vec_t &x, const vec_t &x2) const
The normalized density.
virtual double log_pdf(const vec_t &x, const vec_t &x2) const
The log of the normalized density.
void set_sigma(double sigma)
Set the Gaussian width (must be positive)
invalid argument supplied by user
prob_dens_mdim_gaussian(size_t p_ndim, vec_t &p_peak, mat_t &covar)
Create a distribution from the covariance matrix.
ubvector sum
Normalized partial sum of histogram bins.
virtual size_t dim() const
The dimensionality.
size_t n
Number of original histogram bins.
vec_t vtmp
Temporary storage 2.
virtual double entropy() const
Inverse cumulative distribution function (from the lower tail)
prob_dens_lognormal & operator=(const prob_dens_lognormal &pdg)
Copy constructor with operator=.
virtual void operator()(vec_t &x) const
Sample the distribution.
prob_dens_lognormal()
Create a blank lognormal distribution.
o2scl::rng_gsl r
Base GSL random number generator.
o2scl::prob_dens_gaussian pdg
Standard normal.
double random()
Return a random number in .
std::random_device rd
Desc.
A one-dimensional Gaussian probability density.
A multi-dimensional conditional probability density function.
virtual double cdf(double x) const
The cumulative distribution function (from the lower tail)
rng_gsl r
The GSL random number generator.
prob_dens_gaussian()
Create a standard normal distribution.
mat_t chol
Cholesky decomposition.
virtual double log_pdf(const vec_t &x) const
The log of the normalized density.
mat_t covar_inv
Inverse of the covariance matrix.
virtual double invert_cdf(double cdf) const
The inverse cumulative distribution function.
void set_mu_sigma(double mu, double sigma)
Set the maximum and width of the lognormal distribution.
A one-dimensional histogram class.
prob_dens_gaussian(const prob_dens_gaussian &pdg)
Copy constructor.
virtual double pdf(const vec_t &x) const
The normalized density.
prob_dens_lognormal(const prob_dens_lognormal &pdg)
Copy constructor.
virtual void operator()(const vec_t &x, vec_t &x2) const
Sample the distribution.
vec_t vtmp
Temporary storage 2.
virtual double invert_cdf(double in_cdf) const
The inverse cumulative distribution function.
prob_cond_mdim_gaussian()
Create an empty distribution.
prob_dens_gaussian(double cent, double sigma)
Create a Gaussian distribution with width sigma.
A multi-dimensional Gaussian conditional probability density function.
A multidimensional distribution formed by the product of several one-dimensional distributions.
double norm
Normalization factor.
A one-dimensional probability density function.
virtual double pdf(double x) const
The normalized density.
virtual double entropy() const
The inverse cumulative distribution function.
A one-dimensional probability density over a finite range.
virtual double entropy() const
Entropy of the distribution ( )
virtual void operator()(vec_t &x) const
Sample the distribution.
void cholesky_decomp(const size_t M, mat_t &A)
Compute the in-place Cholesky decomposition of a symmetric positive-definite square matrix...
virtual double metrop_hast(const vec_t &x, vec_t &x2) const
Sample the distribution and return the log of the Metropolis-Hastings ratio.
#define O2SCL_ERR2(d, d2, n)
Set an error, two-string version.
A multi-dimensional probability density function.
virtual double cdf(double x) const
The cumulative distribution function (from the lower tail)
ubvector range
Vector specifying original histogram bins.
A multi-dimensional conditional probability density function independent of the input.
double cent_
Central value.
virtual double pdf(const vec_t &x) const
The normalized density.
virtual size_t dim() const
Return the dimensionality.
search_vec< ubvector > sv
Search through the partial sums.
virtual double log_pdf(double x) const
The log of the normalized density.
virtual double operator()() const
Sample from the specified density.
#define O2SCL_ERR(d, n)
Set an error with message d and code n.
virtual double operator()() const
Sample from the specified density.
double sigma_
Width parameter.
virtual double log_pdf(const vec_t &x, const vec_t &x2) const
The log of the normalized density.
vec_t q
Temporary storage 1.
virtual double log_pdf(double x) const
The log of the normalized density.
A one-dimensional probability density over the positive real numbers.
vec_t q
Temporary storage 1.
virtual size_t dim() const
The dimensionality.
virtual void operator()(const vec_t &x, vec_t &x2) const
Sample the distribution.
virtual double log_pdf(const vec_t &x) const
The log of the normalized density.
A multi-dimensional Gaussian probability density function.
virtual double pdf(const vec_t &x, const vec_t &x2) const
The normalized density.
virtual size_t dim() const
The dimensionality.
Random number generator (GSL)
virtual void operator()(vec_t &x) const
Sample the distribution.
virtual double pdf(const vec_t &x, const vec_t &x2) const
The normalized density.
Searching class for monotonic data with caching.
prob_dens_lognormal(double mu, double sigma)
Create lognormal distribution with mean parameter mu and width parameter sigma.
std::vector< double > u_high
Desc.
double norm
Normalization factor.
A constrained random walk in the shape of a hypercube.
rng_gsl rng
Random number generator.
Probability density function based on a histogram.
size_t ndim
Number of dimensions.
std::vector< prob_dens_func > list
Vector of one-dimensional distributions.
virtual double pdf(double x) const
The normalized density.
virtual size_t dim() const
Return the dimensionality.
virtual double operator()() const
Sample from the specified density.
mat_t covar_inv
Inverse of the covariance matrix.
void set_seed(unsigned long int s)
Set the seed.
void set_seed(unsigned long int s)
Set the seed.
static const double x2[5]
size_t ndim
Number of dimensions.
virtual double pdf(double x) const
The normalized density.
virtual double log_pdf(double x) const
The log of the normalized density.
Lognormal density function.
std::vector< double > u_low
Desc.
double stddev()
Get the Gaussian width.
mat_t chol
Cholesky decomposition.
void set_center(double cent)
Set the center.
virtual void operator()(const vec_t &x, vec_t &x2) const
Sample the distribution.