CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
cbl::statistics::Posterior Class Reference

The class Posterior. More...

#include "Headers/Posterior.h"

Inheritance diagram for cbl::statistics::Posterior:
Collaboration diagram for cbl::statistics::Posterior:

Public Member Functions

std::shared_ptr< ModelParametersparameters () const
 return the posterior parameters More...
 
double operator() (std::vector< double > &pp) const
 the un-normalized posterior More...
 
std::vector< double > weight (const int start=0, const int thin=1) const
 return the internal member m_weights More...
 
double log (std::vector< double > &pp) const
 the logarithm of the un-normalized posterior More...
 
double chi2 (const std::vector< double > parameter={}) const
 the \(\chi^2\) More...
 
void set_model (std::shared_ptr< Model > model=NULL, const std::shared_ptr< ModelParameters > model_parameters=NULL)
 set the model for the posterior analysis More...
 
void set (const std::vector< std::shared_ptr< PriorDistribution >> prior_distributions, const std::shared_ptr< data::Data > data, const std::shared_ptr< Model > model, const LikelihoodType likelihood_type, const std::vector< size_t > x_index, const int w_index, const int seed)
 set the posterior type using the LikelihoodType object More...
 
int get_Nparameters ()
 get the number of parameters used More...
 
int m_get_seed ()
 get the value of the internal variable m_seed More...
 
bool get_m_use_grid ()
 get the value of the internal variable m_use_grid More...
 
std::vector< std::vector< double > > get_parameters ()
 get the values of the internal variable m_parameters More...
 
std::vector< double > get_log_posterior ()
 get the values of the internal variable m_log_posterior More...
 
std::shared_ptr< void > get_m_likelihood_inputs ()
 get the values of the internal variable m_likelihood_inputs More...
 
std::shared_ptr< cbl::statistics::Priorget_m_prior ()
 get the values of the internal variable m_prior More...
 
LogLikelihood_function get_m_log_likelihood_function ()
 get the values of the internal variable m_log_likelihood_function More...
 
Likelihood_function get_m_likelihood_function ()
 get the values of the internal variable m_likelihood_function More...
 
Likelihood_function get_m_likelihood_function_grid ()
 get the values of the internal variable m_likelihood_function_grid More...
 
LogLikelihood_function get_m_log_likelihood_function_grid ()
 get the values of the internal variable m_log_likelihood_function_grid More...
 
void set_response_function (std::shared_ptr< statistics::Model > response)
 set the response function used to compute the super-sample covariance More...
 
std::shared_ptr< statistics::Modelget_response_function ()
 return the response function used to compute the super-sample covariance More...
 
void maximize (const std::vector< double > start, const std::vector< std::vector< double >> parameter_limits, const unsigned int max_iter=10000, const double tol=1.e-6, const double epsilon=1.e-3)
 function that maximizes the posterior, finds the best-fit parameters and store them in the model More...
 
void maximize (const std::vector< double > start, const unsigned int max_iter=10000, const double tol=1.e-6, const double epsilon=1.e-4)
 function that maximize the posterior, find the best-fit parameters and store them in model More...
 
void initialize_chains (const int chain_size, const int n_walkers)
 initialize the chains by drawing from the prior distributions More...
 
void initialize_chains (const int chain_size, const int n_walkers, const double radius, const std::vector< double > start, const unsigned int max_iter=10000, const double tol=1.e-6, const double epsilon=1.e-3)
 initialize the chains in a ball around the posterior best-fit parameter values More...
 
void initialize_chains (const int chain_size, const int n_walkers, std::vector< double > &value, const double radius)
 initialize the chains in a ball around the input parameter values More...
 
void initialize_chains (const int chain_size, const std::vector< std::vector< double >> chain_value)
 initialize the chains with input values More...
 
void initialize_chains (const int chain_size, const int n_walkers, const std::string input_dir, const std::string input_file)
 initialize the chains reading from an input file More...
 
void sample_stretch_move (const double aa=2, const bool parallel=true, const std::string outputFile=par::defaultString, const int start=0, const int thin=1, const int nbins=50)
 sample the posterior using the stretch-move sampler (Foreman-Mackey et al. 2012) More...
 
void importance_sampling (const std::string input_dir, const std::string input_file, const std::vector< size_t > column={}, const int header_lines_to_skip=1, const bool is_FITS_format=false, const bool apply_to_likelihood=false, const int n_walkers=100)
 perform importance sampling More...
 
void read_chain (const std::string input_dir, const std::string input_file, const int n_walkers, const std::vector< size_t > column={}, const int header_lines_to_skip=1, const bool is_FITS_format=false)
 read the chains More...
 
void read_chain_ascii (const std::string input_dir, const std::string input_file, const int n_walkers, const std::vector< size_t > column={}, const int header_lines_to_skip=1)
 read the chains from an ascii file More...
 
void read_chain_fits (const std::string input_dir, const std::string input_file, const int n_walkers, const std::vector< size_t > column)
 read the chains from an ascii file More...
 
void write_chain (const std::string output_dir, const std::string output_file, const int start=0, const int thin=1, const bool is_FITS_format=false, const int prec=5, const int ww=14)
 write the chains obtained after the MCMC sampling More...
 
void write_chain_ascii (const std::string output_dir, const std::string output_file, const int start=0, const int thin=1, const int prec=5, const int ww=14)
 write the chains obtained after the MCMC sampling on an ascii file More...
 
void write_chain_fits (const std::string output_dir, const std::string output_file, const int start=0, const int thin=1)
 write the chains obtained after the MCMC sampling on a FITS file More...
 
void write_maximization_results (const std::string output_dir, const std::string root_file)
 write maximization results on a file More...
 
void show_results (const int start, const int thin, const int nbins=50, const bool show_mode=false, const int ns=-1, const int nb=-1)
 show the results of the MCMC sampling on screen More...
 
void write_results (const std::string output_dir, const std::string root_file, const int start=0, const int thin=1, const int nbins=50, const bool fits=false, const bool compute_mode=false, const int ns=-1, const int nb=-1)
 store the results of the MCMC sampling to file More...
 
void write_model_from_chain (const std::string output_dir, const std::string output_file, const std::vector< double > xx={}, const std::vector< double > yy={}, const int start=0, const int thin=1)
 write the model at xx, yy computing 16th, 50th and 84th percentiles from the MCMC chains More...
 
Constructors/destructors
 Posterior ()=default
 default constructor
 
 Posterior (const std::vector< std::shared_ptr< PriorDistribution >> prior_distributions, const Likelihood &likelihood, const int seed=5341)
 constructor More...
 
 Posterior (const std::vector< std::shared_ptr< PriorDistribution >> prior_distributions, const std::shared_ptr< data::Data > data, const std::shared_ptr< Model > model, const LikelihoodType likelihood_type, const std::vector< size_t > x_index, const int w_index, const int seed=5341)
 constructor More...
 
 Posterior (const std::string file_name, const std::string path_file, const std::vector< int > usecols, const std::vector< std::string > parameter_names, const int skip_header)
 constructor More...
 
 ~Posterior ()=default
 default destructor
 
- Public Member Functions inherited from cbl::statistics::Likelihood
std::shared_ptr< ModelParametersparameters () const
 return the likelihood parameters More...
 
double operator() (std::vector< double > &pp) const
 evaluate the likelihood More...
 
double log (std::vector< double > &parameter) const
 evaluate the logarithm of the likelihood for the input parameters More...
 
void set_data (std::shared_ptr< data::Data > data)
 set the data for the likelihood analysis More...
 
void set_model (std::shared_ptr< Model > model=NULL, const std::shared_ptr< ModelParameters > model_parameters=NULL)
 set the model for the likelihood analysis More...
 
void unset_grid ()
 set the likelihood function with internal values of LikelihoodType
 
std::shared_ptr< Modelget_m_model ()
 get the values of the internal variable m_model More...
 
std::shared_ptr< ModelParametersget_model_parameters ()
 get the values of the internal variable m_model More...
 
std::shared_ptr< data::Dataget_m_data ()
 get the values of the internal variable m_data More...
 
void set_grid (const int npoints, const std::vector< std::vector< double >> parameter_limits, const std::string file, const bool read=false)
 set the likelihood function as a grid, to speed up computation: this only works for one or two free parameters More...
 
void set_function (const LikelihoodType likelihood_type, const std::vector< size_t > x_index={0, 2}, const int w_index=-1, const double prec=1.e-10, const int Nres=-1)
 set the likelihood type using the LikelihoodType object More...
 
void set_function (const LogLikelihood_function likelihood_function)
 set the likelihood function More...
 
void set_log_function (const LogLikelihood_function loglikelihood_function)
 set the natural logarithm of the likelihood function More...
 
void write_results (const std::string dir_output, const std::string file)
 write best-fit results on a file More...
 
void maximize (const std::vector< double > start, const std::vector< std::vector< double >> parameter_limits, const unsigned int max_iter=10000, const double tol=1.e-6, const double epsilon=1.e-3)
 function that maximizes the likelihood, finds the best-fit parameters and stores them in the model More...
 
void write_model (const std::string output_dir, const std::string output_file, const std::vector< double > parameters, const std::vector< double > xx={}, const std::vector< double > yy={})
 write the model at xx at bestfit More...
 
void write_model_at_bestfit (const std::string output_dir, const std::string output_file, const std::vector< double > xx={}, const std::vector< double > yy={})
 write the model at xx at bestfit More...
 
 Likelihood ()
 default constructor
 
 Likelihood (const std::shared_ptr< data::Data > data, const std::shared_ptr< Model > model, const LikelihoodType likelihood_type, const std::vector< size_t > x_index={0, 2}, const int w_index=-1, const std::shared_ptr< ModelParameters > model_parameters=NULL, const double prec=1.e-10, const int Nres=-1)
 constructor More...
 
 Likelihood (const std::shared_ptr< data::Data > data, const std::shared_ptr< Model > model, const Likelihood_function log_likelihood_function, const std::shared_ptr< ModelParameters > model_parameters=NULL)
 constructor More...
 
virtual ~Likelihood ()=default
 default destructor
 

Protected Member Functions

void m_set_seed (const int seed)
 set the internal attribute m_seed and related quantities More...
 
int m_generate_seed ()
 generate a seed More...
 
void m_isSet_response ()
 check if the response function used to compute the super-sample covariance is set
 
- Protected Member Functions inherited from cbl::statistics::Likelihood
void m_set_grid_likelihood_1D (const int npoints, const std::vector< std::vector< double >> parameter_limits, const std::string output_file)
 set the likelihood grid and write the grid on a file More...
 
void m_set_grid_likelihood_1D (const std::string input_file)
 set the likelihood grid from file More...
 
void m_set_grid_likelihood_2D (const int npoints, const std::vector< std::vector< double >> parameter_limits, const std::string output_file)
 set the likelihood grid and write the grid on a file More...
 
void m_set_grid_likelihood_2D (const std::string input_file)
 set the likelihood grid from file More...
 

Protected Attributes

std::shared_ptr< cbl::statistics::Priorm_prior
 the prior distribution
 
std::shared_ptr< statistics::Modelm_response_func = NULL
 response function for the computation of the super-sample covariance
 
std::vector< double > m_log_likelihood
 the logarithm of the likelihood
 
std::vector< double > m_acceptance
 the MCMC acceptance rate
 
int m_seed
 general seed for prior/posterior distribution and sampler
 
std::shared_ptr< cbl::random::UniformRandomNumbers_Intm_seed_generator
 seed generator
 
std::vector< double > m_log_posterior
 the logarithm of the posterior
 
int m_Nparameters
 the number of parameters
 
std::vector< std::vector< double > > m_parameters
 the parameters
 
std::vector< double > m_weight
 the chain weight
 
- Protected Attributes inherited from cbl::statistics::Likelihood
std::shared_ptr< data::Datam_data
 data containers
 
std::shared_ptr< Modelm_model
 model to test
 
std::shared_ptr< void > m_likelihood_inputs
 likelihood inputs
 
std::shared_ptr< ModelParametersm_model_parameters
 likelihood parameters
 
LikelihoodType m_likelihood_type = LikelihoodType::_NotSet_
 type of the likelihood
 
LogLikelihood_function m_log_likelihood_function
 log-likelihood function
 
Likelihood_function m_likelihood_function
 likelihood function
 
LogLikelihood_function m_likelihood_function_grid
 likelihood function on a grid
 
LogLikelihood_function m_log_likelihood_function_grid
 log-likelihood function on a grid
 
std::vector< size_t > m_x_index
 
int m_w_index
 the index in extra info where the bin weight is stored
 
bool m_use_grid = false
 use grid
 

Detailed Description

The class Posterior.

This class is used to define the distribution

Examples
fit.cpp.

Definition at line 55 of file Posterior.h.

Constructor & Destructor Documentation

◆ Posterior() [1/3]

cbl::statistics::Posterior::Posterior ( const std::vector< std::shared_ptr< PriorDistribution >>  prior_distributions,
const Likelihood likelihood,
const int  seed = 5341 
)

constructor

Parameters
prior_distributionsvector containing the priors for the likelihood parameters
likelihoodpointer to an object of type Likelihood
seedgeneral seed for prior/posterior distribution and sampler

Definition at line 67 of file Posterior.cpp.

◆ Posterior() [2/3]

cbl::statistics::Posterior::Posterior ( const std::vector< std::shared_ptr< PriorDistribution >>  prior_distributions,
const std::shared_ptr< data::Data data,
const std::shared_ptr< Model model,
const LikelihoodType  likelihood_type,
const std::vector< size_t >  x_index,
const int  w_index,
const int  seed = 5341 
)

constructor

Parameters
prior_distributionsvector containing the priors for the likelihood parameters
dataobject of type data
modelobject of type model
likelihood_typetype of the likelihood
x_indexindex(s) of the extra info std::vector containing the point(s) where the model is evaluated
w_indexstd::vector containing the data point weight
seedgeneral seed for prior/posterior distribution and sampler

Definition at line 58 of file Posterior.cpp.

◆ Posterior() [3/3]

cbl::statistics::Posterior::Posterior ( const std::string  file_name,
const std::string  path_file,
const std::vector< int >  usecols,
const std::vector< std::string >  parameter_names,
const int  skip_header 
)

constructor

Parameters
file_namethe name of the input file
path_filethe path of the input file
parameter_namesthe vector containing the names of parameters
usecolsthe columns to read
skip_headerthe number of lines to skip

Definition at line 84 of file Posterior.cpp.

Member Function Documentation

◆ chi2()

double cbl::statistics::Posterior::chi2 ( const std::vector< double >  parameter = {}) const

the \(\chi^2\)

this function computes

\[ \chi^2 \equiv -2\log(\frac{P(\vec{\theta} | \vec{d})}{Pr(\vec{\theta})}) \]

where \(P(\vec{\theta} | \vec{d})\) is the posterior and \(Pr(\vec{\theta})\) is the prior

Parameters
parametervector containing the model parameters; if not provided in input, the chi2 is evaluated at best-fit parameter values
Returns
\(\chi^2\)

Definition at line 155 of file Posterior.cpp.

◆ get_log_posterior()

std::vector<double> cbl::statistics::Posterior::get_log_posterior ( )
inline

get the values of the internal variable m_log_posterior

Returns
vector of double containing the logPosterior for each MCMC step

Definition at line 334 of file Posterior.h.

◆ get_m_likelihood_function()

Likelihood_function cbl::statistics::Posterior::get_m_likelihood_function ( )
inline

get the values of the internal variable m_likelihood_function

Returns
m_likelihood_function

Definition at line 362 of file Posterior.h.

◆ get_m_likelihood_function_grid()

Likelihood_function cbl::statistics::Posterior::get_m_likelihood_function_grid ( )
inline

get the values of the internal variable m_likelihood_function_grid

Returns
m_likelihood_function_grid

Definition at line 369 of file Posterior.h.

◆ get_m_likelihood_inputs()

std::shared_ptr<void> cbl::statistics::Posterior::get_m_likelihood_inputs ( )
inline

get the values of the internal variable m_likelihood_inputs

Returns
pointer containing m_likelihood_inputs

Definition at line 341 of file Posterior.h.

◆ get_m_log_likelihood_function()

LogLikelihood_function cbl::statistics::Posterior::get_m_log_likelihood_function ( )
inline

get the values of the internal variable m_log_likelihood_function

Returns
m_log_likelihood_function

Definition at line 355 of file Posterior.h.

◆ get_m_log_likelihood_function_grid()

LogLikelihood_function cbl::statistics::Posterior::get_m_log_likelihood_function_grid ( )
inline

get the values of the internal variable m_log_likelihood_function_grid

Returns
m_log_likelihood_function_grid

Definition at line 376 of file Posterior.h.

◆ get_m_prior()

std::shared_ptr<cbl::statistics::Prior> cbl::statistics::Posterior::get_m_prior ( )
inline

get the values of the internal variable m_prior

Returns
pointer containing the prior

Definition at line 348 of file Posterior.h.

◆ get_m_use_grid()

bool cbl::statistics::Posterior::get_m_use_grid ( )
inline

get the value of the internal variable m_use_grid

Returns
the boolean value of m_use_grid

Definition at line 320 of file Posterior.h.

◆ get_Nparameters()

int cbl::statistics::Posterior::get_Nparameters ( )
inline

get the number of parameters used

Returns
the value of the internal variable m_Nparameters

Definition at line 306 of file Posterior.h.

◆ get_parameters()

std::vector<std::vector<double> > cbl::statistics::Posterior::get_parameters ( )
inline

get the values of the internal variable m_parameters

Returns
vector of vector containing parameters for each MCMC step

Definition at line 327 of file Posterior.h.

◆ get_response_function()

std::shared_ptr<statistics::Model> cbl::statistics::Posterior::get_response_function ( )
inline

return the response function used to compute the super-sample covariance

Returns
pointer to the response function

Definition at line 392 of file Posterior.h.

◆ importance_sampling()

void cbl::statistics::Posterior::importance_sampling ( const std::string  input_dir,
const std::string  input_file,
const std::vector< size_t >  column = {},
const int  header_lines_to_skip = 1,
const bool  is_FITS_format = false,
const bool  apply_to_likelihood = false,
const int  n_walkers = 100 
)

perform importance sampling

Importance sampling is a convenient technique to join independet datasets.

This function takes in input a chain and computes the posterior looping over all entries. It's possible to specify the columns, in case the input chain has different ordering, or larger number of parameters.

It's possible either to sum or overwrite the log likelihood.

Parameters
input_dirthe input directory
input_filethe input file
columnthe column of the input file to be read
header_lines_to_skipthe lines to be skipped in the chain file
is_FITS_formattrue \(\rightarrow\) the format of the input file is FITS; false \(\rightarrow\) the format of the input file is ASCII
apply_to_likelihoodtrue \(\rightarrow\) the likelihood ratio is used as weight; false \(\rightarrow\) the posterior ratio is used as weight
n_walkersthis parameter is used here only to parallelize the computation
Warning
column is used only for ASCII chain files

Definition at line 312 of file Posterior.cpp.

◆ initialize_chains() [1/5]

void cbl::statistics::Posterior::initialize_chains ( const int  chain_size,
const int  n_walkers 
)

initialize the chains by drawing from the prior distributions

the starting values of the chain are extracted from the (possibly different) distributions of the priors

Parameters
chain_sizethe chain lenght
n_walkersthe number of parallel chains

Definition at line 615 of file Posterior.cpp.

◆ initialize_chains() [2/5]

void cbl::statistics::Posterior::initialize_chains ( const int  chain_size,
const int  n_walkers,
const double  radius,
const std::vector< double >  start,
const unsigned int  max_iter = 10000,
const double  tol = 1.e-6,
const double  epsilon = 1.e-3 
)

initialize the chains in a ball around the posterior best-fit parameter values

the starting values of the chain are extracted from uniform distributions in the range [parameter-radius, parameter+radius] (for each likelihood parameter)

this function first maximizes the posterior, starting the computation at the values of the input vector 'start', then it inizializes the chain

Parameters
chain_sizethe chain lenght
n_walkersthe number of parallel chains
radiusradius of the ball in parameter space
startstd::vector containing initial values for the posterior maximization
max_iterthe maximum number of iterations
tolthe tolerance in finding convergence
epsilonthe simplex side

Definition at line 625 of file Posterior.cpp.

◆ initialize_chains() [3/5]

void cbl::statistics::Posterior::initialize_chains ( const int  chain_size,
const int  n_walkers,
const std::string  input_dir,
const std::string  input_file 
)

initialize the chains reading from an input file

the starting values of the chain are get from the last lines of an input chain file; it can be used to continue an MCMC sampling computation

Parameters
chain_sizethe chain lenght
n_walkersthe number of parallel chains
input_dirthe input directory
input_filethe input file

Definition at line 661 of file Posterior.cpp.

◆ initialize_chains() [4/5]

void cbl::statistics::Posterior::initialize_chains ( const int  chain_size,
const int  n_walkers,
std::vector< double > &  value,
const double  radius 
)

initialize the chains in a ball around the input parameter values

the starting values of the chain are extracted from uniform distributions in the range [value[i]-radius, value[i]+radius] (for each i-th likelihood parameter)

Parameters
chain_sizethe chain lenght
n_walkersthe number of parallel chains
valuevector containing the input values, centres of the ball in the parameter space
radiusradius of the ball in the parameter space

Definition at line 637 of file Posterior.cpp.

◆ initialize_chains() [5/5]

void cbl::statistics::Posterior::initialize_chains ( const int  chain_size,
const std::vector< std::vector< double >>  chain_value 
)

initialize the chains with input values

the starting values of the chain are the elements of the input matrix 'chain_values'

Parameters
chain_sizethe chain lenght
chain_valuematrix of size (n_walkers, n_parameters), starting values of the chain

Definition at line 647 of file Posterior.cpp.

◆ log()

double cbl::statistics::Posterior::log ( std::vector< double > &  pp) const

the logarithm of the un-normalized posterior

\[ P(\vec{\theta} | \vec{d}) = \mathcal{L}(\vec{d}|\vec{\theta}) \cdot Pr(\vec{\theta}) \]

where \(P(\vec{\theta} | \vec{d})\) is the posterior, \(\mathcal{L}(\vec{d}|\vec{\theta})\) is the likelihood, \(Pr(\vec{\theta})\) is the prior, \(\vec{\theta}\) is the set of model parameters and \(\vec{d}\) is the dataset

Parameters
ppthe parameters
Returns
the logarithm of the un-normalized posterior

Definition at line 137 of file Posterior.cpp.

◆ m_generate_seed()

int cbl::statistics::Posterior::m_generate_seed ( )
inlineprotected

generate a seed

Returns
a seed generated from m_seed_generator

Definition at line 104 of file Posterior.h.

◆ m_get_seed()

int cbl::statistics::Posterior::m_get_seed ( )
inline

get the value of the internal variable m_seed

Returns
the value of the internal variable m_seed

Definition at line 313 of file Posterior.h.

◆ m_set_seed()

void cbl::statistics::Posterior::m_set_seed ( const int  seed)
protected

set the internal attribute m_seed and related quantities

Parameters
seedthe seed

Definition at line 47 of file Posterior.cpp.

◆ maximize() [1/2]

void cbl::statistics::Posterior::maximize ( const std::vector< double >  start,
const std::vector< std::vector< double >>  parameter_limits,
const unsigned int  max_iter = 10000,
const double  tol = 1.e-6,
const double  epsilon = 1.e-3 
)
inline

function that maximizes the posterior, finds the best-fit parameters and store them in the model

this function exploits the Nelder-Mead method https://en.wikipedia.org/wiki/Nelder%E2%80%93Mead_method

the algorithm defines a simplex (i.e a k-dimensional polytope which is the convex hull of its k+1 vertices) in the parameter space. At each step, it identifies the simplex vertex at which the function to be minimised (i.e. the negative posterior in this case) has the greatest value, and moves it, via reflections and scaling, to a new position in which the function has a lower value. This iteration stops when the simplex area becomes lower than the tolerance. For instance, in 2D, the starting vertices of the simplex (a triangle in 2D) are the following: (start[0], start[1]) ; (start[0]+epsilon, start[1]) ; (start[0], start[1]+epsilon)

Parameters
startstd::vector containing initial values for the posterior maximization
parameter_limitslimits for the parameters
max_iterthe maximum number of iterations
tolthe tolerance in finding convergence
epsilonthe simplex side

Definition at line 428 of file Posterior.h.

◆ maximize() [2/2]

void cbl::statistics::Posterior::maximize ( const std::vector< double >  start,
const unsigned int  max_iter = 10000,
const double  tol = 1.e-6,
const double  epsilon = 1.e-4 
)

function that maximize the posterior, find the best-fit parameters and store them in model

Parameters
startstd::vector containing initial values for the posterior maximization
max_iterthe maximum number of iterations
tolthe tolerance in finding convergence
epsilonthe simplex side

Definition at line 212 of file Posterior.cpp.

◆ operator()()

double cbl::statistics::Posterior::operator() ( std::vector< double > &  pp) const

the un-normalized posterior

\[ P((\vec{\theta} | \vec{d}) = \mathcal{L}(\vec{d}|\vec{\theta}) \cdot Pr(\vec{\theta}) \]

where \(P\) is the posterior, \(\mathcal{L}(\vec{d}|\vec{\theta})\) is the likelihood and \(Pr(\vec{\theta})\) is the prior

Parameters
ppthe parameters
Returns
the value of the un-normalized posterior

Definition at line 105 of file Posterior.cpp.

◆ parameters()

std::shared_ptr<ModelParameters> cbl::statistics::Posterior::parameters ( ) const
inline

return the posterior parameters

Returns
pointer containing the posterior parameters

Definition at line 197 of file Posterior.h.

◆ read_chain()

void cbl::statistics::Posterior::read_chain ( const std::string  input_dir,
const std::string  input_file,
const int  n_walkers,
const std::vector< size_t >  column = {},
const int  header_lines_to_skip = 1,
const bool  is_FITS_format = false 
)

read the chains

Parameters
input_dirthe input directory
input_filethe input file
n_walkersthe number of parallel chains
columnthe columns of the input file to be read
header_lines_to_skipthe lines to be skipped in the chain file
is_FITS_formattrue \(\rightarrow\) the format of the input file is FITS; false \(\rightarrow\) the format of the input file is ASCII

Definition at line 361 of file Posterior.cpp.

◆ read_chain_ascii()

void cbl::statistics::Posterior::read_chain_ascii ( const std::string  input_dir,
const std::string  input_file,
const int  n_walkers,
const std::vector< size_t >  column = {},
const int  header_lines_to_skip = 1 
)

read the chains from an ascii file

Parameters
input_dirthe intput directory
input_filethe intput file
n_walkersthe number of parallel chains
columnthe columns of the input file to be read.
header_lines_to_skipthe lines to be skipped in the chain file

Definition at line 373 of file Posterior.cpp.

◆ read_chain_fits()

void cbl::statistics::Posterior::read_chain_fits ( const std::string  input_dir,
const std::string  input_file,
const int  n_walkers,
const std::vector< size_t >  column 
)

read the chains from an ascii file

Parameters
input_dirthe input directory
input_filethe input file
n_walkersthe number of parallel chains
columnthe columns of the input file to be read

Definition at line 446 of file Posterior.cpp.

◆ sample_stretch_move()

void cbl::statistics::Posterior::sample_stretch_move ( const double  aa = 2,
const bool  parallel = true,
const std::string  outputFile = par::defaultString,
const int  start = 0,
const int  thin = 1,
const int  nbins = 50 
)

sample the posterior using the stretch-move sampler (Foreman-Mackey et al. 2012)

Parameters
aathe parameter of the \(g(z)\) distribution
parallelfalse \(\rightarrow\) non-parallel sampler; true \(\rightarrow\) parallel sampler
outputFileoutput file where the chains are written during run-time. Leave it to default value to have no output. WARNING: this option is intended for debug. It only works for the non-parallelized stretch-move algorithm. The chain in output will be written in a different format with respect to the method cbl::statistics::Posterior::write_chain::ascii col1) chain step col2) walker index col3-npar) parameter values col npar+3) value of the posterior
startthe minimum chain position used to compute the median
thinthe step used for chain dilution
nbinsthe number of bins to estimate the posterior distribution, used to assess its properties
Warning
if parallel is set true, than pointers cannot be used inside the posterior function

Definition at line 268 of file Posterior.cpp.

◆ set()

void cbl::statistics::Posterior::set ( const std::vector< std::shared_ptr< PriorDistribution >>  prior_distributions,
const std::shared_ptr< data::Data data,
const std::shared_ptr< Model model,
const LikelihoodType  likelihood_type,
const std::vector< size_t >  x_index,
const int  w_index,
const int  seed 
)

set the posterior type using the LikelihoodType object

Parameters
prior_distributionsvector containing the priors for the likelihood parameters
datapointer to an object of type Data
modelpointer to an object of type model
likelihood_typethe likelihood type, specified with the LikelihoodType object
x_indexindex(s) of the extra info std::vector containing the point(s) where to evaluate the model
w_indexstd::vector containing the data point weight
seedthe seed

Definition at line 190 of file Posterior.cpp.

◆ set_model()

void cbl::statistics::Posterior::set_model ( std::shared_ptr< Model model = NULL,
const std::shared_ptr< ModelParameters model_parameters = NULL 
)

set the model for the posterior analysis

Parameters
modelpointer to the model
model_parametersparameters of the model

Definition at line 166 of file Posterior.cpp.

◆ set_response_function()

void cbl::statistics::Posterior::set_response_function ( std::shared_ptr< statistics::Model response)
inline

set the response function used to compute the super-sample covariance

Parameters
responsepointer to the response statistics::Model object

Definition at line 384 of file Posterior.h.

◆ show_results()

void cbl::statistics::Posterior::show_results ( const int  start,
const int  thin,
const int  nbins = 50,
const bool  show_mode = false,
const int  ns = -1,
const int  nb = -1 
)

show the results of the MCMC sampling on screen

if the covariance matrix has been estimated from a set of mock catalogues, and the input parameters ns (number of samples used to estimate the covariance matrix) and nb (number of data measurements, e.g. the bins of the dataset) are provided (>0), then the parameter errors ( \(\sigma_p\)) will be corrected to take into account the uncertainities in the covariance estimate (Percival et al. 2014):

\[ \sigma_p = \sqrt{\frac{1+B(n_b-n_p)}{1+A+B(n_p+1)}} \]

where

\[ A = \frac{2}{(n_s-n_b-1)(n_s-n_b-4)} \,, \]

\[ B = \frac{(n_s-n_b-2)}{(n_s-n_b-1)(n_s-n_b-4)} \,. \]

this correction can be applied only if the likelihood is Gaussian. Morever, the inverce covariance matrix estimator has to be corrected to take into account the inverse Wishart distribution (Hartlap, Simon and Schneider 2006).

Parameters
startthe minimum chain position to be written
thinthe step used for dilution on screen
nbinsthe number of bins to estimate the posterior distribution, used to assess its properties
show_modetrue \(\rightarrow\) show the posterior mode; false \(\rightarrow\) do not show the posterior mode
nsnumber of samples used to estimate the covariance matrix
nbnumber of data measurements, e.g. the bins of the dataset

Definition at line 732 of file Posterior.cpp.

◆ weight()

vector< double > cbl::statistics::Posterior::weight ( const int  start = 0,
const int  thin = 1 
) const

return the internal member m_weights

Parameters
startthe minimum chain position to be written
thinthe step used for dilution
Returns
vector containing the posterior weights

Definition at line 120 of file Posterior.cpp.

◆ write_chain()

void cbl::statistics::Posterior::write_chain ( const std::string  output_dir,
const std::string  output_file,
const int  start = 0,
const int  thin = 1,
const bool  is_FITS_format = false,
const int  prec = 5,
const int  ww = 14 
)

write the chains obtained after the MCMC sampling

Parameters
output_dirthe output directory
output_filethe output file
startthe minimum chain position to be written
thinthe step used for dilution
is_FITS_formattrue \(\rightarrow\) the format of the input file is FITS; false \(\rightarrow\) the format of the input file is ASCII
precdecimal precision
wwnumber of characters to be used as field width
Warning
column only work for ascii chain file

Definition at line 490 of file Posterior.cpp.

◆ write_chain_ascii()

void cbl::statistics::Posterior::write_chain_ascii ( const std::string  output_dir,
const std::string  output_file,
const int  start = 0,
const int  thin = 1,
const int  prec = 5,
const int  ww = 14 
)

write the chains obtained after the MCMC sampling on an ascii file

Parameters
output_dirthe output directory
output_filethe output file
startthe minimum chain position to be written
thinthe step used for dilution
precdecimal precision
wwnumber of characters to be used as field width

Definition at line 502 of file Posterior.cpp.

◆ write_chain_fits()

void cbl::statistics::Posterior::write_chain_fits ( const std::string  output_dir,
const std::string  output_file,
const int  start = 0,
const int  thin = 1 
)

write the chains obtained after the MCMC sampling on a FITS file

Parameters
output_dirthe output directory
output_filethe output file
startthe minimum chain position to be written
thinthe step used for dilution

Definition at line 560 of file Posterior.cpp.

◆ write_maximization_results()

void cbl::statistics::Posterior::write_maximization_results ( const std::string  output_dir,
const std::string  root_file 
)

write maximization results on a file

Parameters
output_dirthe output directory
root_filethe root of the output file to be written

Definition at line 699 of file Posterior.cpp.

◆ write_model_from_chain()

void cbl::statistics::Posterior::write_model_from_chain ( const std::string  output_dir,
const std::string  output_file,
const std::vector< double >  xx = {},
const std::vector< double >  yy = {},
const int  start = 0,
const int  thin = 1 
)

write the model at xx, yy computing 16th, 50th and 84th percentiles from the MCMC chains

Parameters
output_dirthe output directory
output_filethe output file
xxvector of points at which the model is computed
yyvector of points at which the model is computed
startthe minimum chain position to be written
thinthe step used for dilution on screen

Definition at line 753 of file Posterior.cpp.

◆ write_results()

void cbl::statistics::Posterior::write_results ( const std::string  output_dir,
const std::string  root_file,
const int  start = 0,
const int  thin = 1,
const int  nbins = 50,
const bool  fits = false,
const bool  compute_mode = false,
const int  ns = -1,
const int  nb = -1 
)

store the results of the MCMC sampling to file

this function stores to file the posterior mean, the posterior standard deviation, the posterior median, 18th and 82th posterior percentiles, and, optionally, the posterior mode.

If the covariance matrix has been estimated from a set of mock catalogues, and the input parameters ns (number of samples used to estimate the covariance matrix) and nb (number of data measurements, e.g. the bins of the dataset) are provided (>0), then the parameter errors ( \(\sigma_p\)) will be corrected to take into account the uncertainities in the covariance estimate (Percival et al. 2014):

\[ \sigma_p = \sqrt{\frac{1+B(n_b-n_p)}{1+A+B(n_p+1)}} \]

where

\[ A = \frac{2}{(n_s-n_b-1)(n_s-n_b-4)} \,, \]

\[ B = \frac{(n_s-n_b-2)}{(n_s-n_b-1)(n_s-n_b-4)} \,. \]

this correction can be applied only if the likelihood is Gaussian. Morever, the inverce covariance matrix estimator has to be corrected to take into account the inverse Wishart distribution (Hartlap, Simon and Schneider 2006).

Parameters
output_dirthe output directory
root_filethe root of the output file to be written
startthe minimum chain position to be written
thinthe step used for dilution on screen
nbinsthe number of bins to estimate the posterior distribution, used to assess its properties
fitsfalse \(\rightarrow\) ascii file; true \(\rightarrow\) fits file
compute_modetrue \(\rightarrow\) compute the posterior mode; false \(\rightarrow\) do not compute the posterior mode
nsnumber of samples used to estimate the covariance matrix
nbnumber of data measurements, e.g. the bins of the dataset

Definition at line 741 of file Posterior.cpp.


The documentation for this class was generated from the following files: