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

The class PosteriorParameters. More...

#include <PosteriorParameters.h>

Inheritance diagram for cbl::statistics::PosteriorParameters:
Collaboration diagram for cbl::statistics::PosteriorParameters:

Public Member Functions

size_t nparameters_free () const override
 return the number of free parameters More...
 
std::string status (const int p) const
 return the model parameter status More...
 
std::vector< std::string > status () const
 return all the model parameter status More...
 
std::vector< unsigned int > free_parameter () const override
 return the private member m_free_parameter More...
 
size_t nparameters_fixed () const override
 return the number of fixed parameters More...
 
std::vector< unsigned int > fixed_parameter () const
 return the private member m_fixed_parameter More...
 
std::vector< double > full_parameter (const std::vector< double > parameter_value) const override
 return all the model parameters More...
 
void set_parameters (const size_t nparameters, const std::vector< std::shared_ptr< PriorDistribution >> priorDistributions, std::vector< ParameterType > parameterTypes, std::vector< std::string > parameterNames) override
 set the parameter More...
 
void show_results (const int start, const int thin, const int nbins, const int seed=34121, const bool show_mode=false, const int ns=-1, const int nb=-1, const std::vector< double > weight={})
 show the results on the standard output More...
 
void write_results (const std::string dir, const std::string file, const int start, const int thin, const int nbins, const int seed=34121, const bool compute_mode=false, const int ns=-1, const int nb=-1, const std::vector< double > weight={})
 store the results to file More...
 
size_t chain_size () const
 return the private member m_chain_size More...
 
size_t chain_nwalkers () const
 return the private member m_chain_nwalkers More...
 
void set_chain (const size_t size, const size_t nwalkers)
 set the chain More...
 
void reset_chain ()
 reset the chain using m_size and m_nwalkers
 
void expand_chain (const int append)
 expand the already existing chain More...
 
double chain_value (const int param, const int pos, const int ww) const
 return the private member m_chain_value at the pos step for the ww-th walker, for the chosen parameter More...
 
std::vector< double > chain_value_parameter (const int pos, const int ww) const
 return the private member m_values at the pp-th step for the ww-th step for all the parameters More...
 
std::vector< double > parameter_chain_values (const int param, const int start=0, const int thin=1) const
 return all the chain values for a parameter More...
 
void set_chain_value (const int param, const int pos, const int ww, const double value)
 set the private member m_chain_value at the pp-th step for the ww-th step More...
 
void set_chain_values (const std::vector< std::vector< double >> values, const int nwalkers)
 set the chain values More...
 
void set_chain_values (const std::vector< std::vector< std::vector< double >>> values)
 set the chain values More...
 
void initialize_chain (const std::vector< std::vector< double >> values)
 initialize the chain values More...
 
void initialize_chain_from_prior ()
 initialize the chain values random sampling the parameter priors
 
void initialize_chain_ball (const std::vector< double > center, const double radius, const double seed)
 initialize the chain values More...
 
void initialize_chain_ball_bestfit (const double radius, const double seed)
 initialize the chain values around bestfit values More...
 
Constructors/destructors
 PosteriorParameters ()=default
 default constructor
 
 PosteriorParameters (const size_t nparameters, const std::vector< std::shared_ptr< PriorDistribution >> priorDistributions, std::vector< ParameterType > parameterTypes, std::vector< std::string > parameterNames)
 constructor for PosteriorParameters More...
 
 ~PosteriorParameters ()=default
 default destructor
 
Member functions used to set private/protected of the PosteriorParameters
void set_parameter_covariance (const int start=0, const int thin=1)
 set the internal method m_parameter_covariance More...
 
double parameter_covariance (const int i, const int j) const
 return the protected member m_parameter_covariance More...
 
std::vector< std::vector< double > > parameter_covariance () const
 return the protected member m_parameter_covariance More...
 
Member functions used to interact with prior distribution
void set_prior_distribution (const int p, const std::shared_ptr< PriorDistribution > priorDistribution)
 set the prior distribution for the p-th parameter More...
 
void set_prior_distribution (const std::vector< std::shared_ptr< PriorDistribution >> priorDistribution)
 set the prior distributions for the parameters More...
 
void set_prior_distribution_seed (const std::shared_ptr< random::UniformRandomNumbers_Int > ran_generator)
 set the prior distributions for the parameters More...
 
std::shared_ptr< PriorDistributionprior_distribution (const int p) const
 get the prior distribution for the p-th parameter More...
 
std::vector< std::shared_ptr< PriorDistribution > > prior_distribution () const
 get the prior distribution for the p-th parameter More...
 
std::shared_ptr< Priorprior () const
 get the prior function More...
 
Member functions used to get/set best fit values
double bestfit_value (const int p) const
 get the protected member m_value More...
 
std::vector< double > bestfit_value () const
 get the protected member m_parameter_bestfit_value More...
 
void set_bestfit_values (const std::vector< double > bestfit_value)
 set the protected member m_bestfit_value More...
 
void set_bestfit_values (const int start, const int thin, const int nbins, const int seed) override
 set the protected member m_bestfit_value More...
 
void write_bestfit_info ()
 write bestfit info
 
Member functions used to interact with posterior distribution
void set_posterior_distribution (const int start, const int thin, const int nbins, const int seed=34121, const std::vector< double > weight={})
 set the posterior distribution from the chains More...
 
std::shared_ptr< PosteriorDistributionposterior_distribution (const int par) const
 get the posterior distribution for the chosen parameter More...
 
- Public Member Functions inherited from cbl::statistics::ModelParameters
virtual void reset ()
 reset the parameter vectors
 
virtual std::vector< double > chain_value_parameters (const int pos, const int ww) const
 return the private member m_values at the pp-th step for the ww-th step for all the parameters More...
 
 ModelParameters ()=default
 default constructor
 
 ModelParameters (const size_t nparameters, std::vector< ParameterType > parameterTypes, std::vector< std::string > parameterNames)
 constructor for ModelParameters More...
 
 ModelParameters (const size_t nparameters, std::vector< std::string > parameterNames)
 constructor for ModelParameters More...
 
 ~ModelParameters ()=default
 default destructor
 
size_t nparameters () const
 return the total number of parameters More...
 
size_t nparameters_base () const
 return the number of base parameters More...
 
std::vector< unsigned int > base_parameter () const
 return the private member m_base_parameters More...
 
size_t nparameters_derived () const
 return the number of derived parameters More...
 
size_t nparameters_correlated () const
 return the number of correlated parameters More...
 
std::vector< unsigned int > derived_parameter () const
 return the private member m_derived_parameter More...
 
ParameterType type (const int p) const
 return the model parameter type More...
 
std::vector< ParameterTypetype () const
 return all the model parameter names More...
 
std::string name (const int p) const
 return the model parameter name More...
 
std::vector< std::string > name () const
 return all the model parameter names More...
 
virtual void set_parameters (const size_t nparameters, std::vector< ParameterType > parameterTypes, std::vector< std::string > parameterNames)
 set the parameter More...
 
virtual void set_parameters (const size_t nparameters, std::vector< std::string > parameterNames)
 set the parameter More...
 
virtual void free (const int p)
 set m_fixed to false More...
 
virtual void fix (const int p, const double value)
 fix the parameter at the input value More...
 
virtual void fix_at_bestfit (const int p)
 fix the parameter at the bestfit value, contained in m_bestfit_value More...
 

Protected Member Functions

void m_set_parameter_type () override
 private member to set the parameter types
 
int m_inds_to_index (const int pp, const int ww) const
 get the position in the vector m_values from position index and walker index More...
 

Protected Attributes

size_t m_nparameters_free = 0
 the numbers of free parameters
 
size_t m_nparameters_fixed = 0
 the number of fixed parameters
 
std::vector< unsigned int > m_fixed_parameter
 the indexes of fixed parameters
 
std::vector< unsigned int > m_free_parameter
 the indexes of the free parameters
 
std::vector< std::shared_ptr< PriorDistribution > > m_parameter_prior
 the parameter prior distributions
 
std::vector< std::shared_ptr< PriorDistribution > > m_correlated_parameter_prior
 the correlated parameter prior distributions
 
std::vector< std::shared_ptr< PosteriorDistribution > > m_posterior_distribution
 the parameter posterior distributions
 
std::vector< double > m_parameter_bestfit_value
 the best-fit parameter values, either the medians of the MCMC chain or the maximum of the posterior (depending of what has been calculated)
 
std::vector< std::vector< double > > m_parameter_covariance
 the function parameter covariance matrix
 
size_t m_chain_size
 the lenght of the chain
 
size_t m_chain_nwalkers
 the number of parallel walkers
 
std::vector< size_t > m_multipriors
 the lenght of each set of multidimensional prior
 
std::vector< std::vector< double > > m_chain_value
 content of the chain
 
- Protected Attributes inherited from cbl::statistics::ModelParameters
std::vector< ParameterTypem_parameter_type
 model parameters
 
std::vector< std::string > m_parameter_name
 model parameter types
 
size_t m_nparameters = 0
 number of parameters
 
size_t m_nparameters_base = 0
 number of base parameters
 
size_t m_nparameters_derived = 0
 number of derived parameters
 
size_t m_nparameters_correlated = 0
 number of correlated parameters
 
std::vector< unsigned int > m_base_parameter
 indexes of base parameters
 
std::vector< unsigned int > m_derived_parameter
 indexes of the derived parameters
 

Detailed Description

The class PosteriorParameters.

"Headers/PosteriorParameters.h"

This class is used to define the model parameters in a bayesian framework

Definition at line 53 of file PosteriorParameters.h.

Constructor & Destructor Documentation

◆ PosteriorParameters()

cbl::statistics::PosteriorParameters::PosteriorParameters ( const size_t  nparameters,
const std::vector< std::shared_ptr< PriorDistribution >>  priorDistributions,
std::vector< ParameterType parameterTypes,
std::vector< std::string >  parameterNames 
)

constructor for PosteriorParameters

Parameters
nparametersthe number of parameters
priorDistributionsthe priorDistributions
parameterTypesthe parameter types
parameterNamesthe parameter names

Definition at line 193 of file PosteriorParameters.cpp.

Member Function Documentation

◆ bestfit_value() [1/2]

std::vector< double > cbl::statistics::PosteriorParameters::bestfit_value ( ) const
virtual

get the protected member m_parameter_bestfit_value

Returns
the parameter bestfit values

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 437 of file PosteriorParameters.cpp.

◆ bestfit_value() [2/2]

double cbl::statistics::PosteriorParameters::bestfit_value ( const int  p) const
virtual

get the protected member m_value

Parameters
pthe p-th parameter
Returns
the bestfit value of the parameter

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 426 of file PosteriorParameters.cpp.

◆ chain_nwalkers()

size_t cbl::statistics::PosteriorParameters::chain_nwalkers ( ) const
inlinevirtual

return the private member m_chain_nwalkers

Returns
the chain size

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 536 of file PosteriorParameters.h.

◆ chain_size()

size_t cbl::statistics::PosteriorParameters::chain_size ( ) const
inlinevirtual

return the private member m_chain_size

Returns
the chain size

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 529 of file PosteriorParameters.h.

◆ chain_value()

double cbl::statistics::PosteriorParameters::chain_value ( const int  param,
const int  pos,
const int  ww 
) const
inlinevirtual

return the private member m_chain_value at the pos step for the ww-th walker, for the chosen parameter

Parameters
paramthe parameter index
posthe position in the chain
wwthe walker index
Returns
the chain value

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 577 of file PosteriorParameters.h.

◆ chain_value_parameter()

std::vector< double > cbl::statistics::PosteriorParameters::chain_value_parameter ( const int  pos,
const int  ww 
) const

return the private member m_values at the pp-th step for the ww-th step for all the parameters

Parameters
posthe position in the chain
wwthe walker index
Returns
the chain value

Definition at line 581 of file PosteriorParameters.cpp.

◆ expand_chain()

void cbl::statistics::PosteriorParameters::expand_chain ( const int  append)
virtual

expand the already existing chain

Parameters
appendthe lenght of the empty chunk of the chain

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 562 of file PosteriorParameters.cpp.

◆ fixed_parameter()

std::vector<unsigned int> cbl::statistics::PosteriorParameters::fixed_parameter ( ) const
inlinevirtual

return the private member m_fixed_parameter

Returns
the private member m_fixed_parameter

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 199 of file PosteriorParameters.h.

◆ free_parameter()

std::vector<unsigned int> cbl::statistics::PosteriorParameters::free_parameter ( ) const
inlineoverridevirtual

return the private member m_free_parameter

Returns
the private member m_free_parameter

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 184 of file PosteriorParameters.h.

◆ full_parameter()

std::vector< double > cbl::statistics::PosteriorParameters::full_parameter ( const std::vector< double >  parameter_value) const
overridevirtual

return all the model parameters

Parameters
parameter_valuevector of free parameters
Returns
all the parameter values

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 105 of file PosteriorParameters.cpp.

◆ initialize_chain()

void cbl::statistics::PosteriorParameters::initialize_chain ( const std::vector< std::vector< double >>  values)
virtual

initialize the chain values

Parameters
valuesthe starting values

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 649 of file PosteriorParameters.cpp.

◆ initialize_chain_ball()

void cbl::statistics::PosteriorParameters::initialize_chain_ball ( const std::vector< double >  center,
const double  radius,
const double  seed 
)
virtual

initialize the chain values

Parameters
centerthe ball center
radiusthe ball radius
seedthe random number generator seed

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 683 of file PosteriorParameters.cpp.

◆ initialize_chain_ball_bestfit()

void cbl::statistics::PosteriorParameters::initialize_chain_ball_bestfit ( const double  radius,
const double  seed 
)
virtual

initialize the chain values around bestfit values

Parameters
radiusthe ball radius
seedthe random number generator seed

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 714 of file PosteriorParameters.cpp.

◆ m_inds_to_index()

int cbl::statistics::PosteriorParameters::m_inds_to_index ( const int  pp,
const int  ww 
) const
inlineprotected

get the position in the vector m_values from position index and walker index

Parameters
ppthe position in the chain
wwthe walker
Returns
object of class Chain.

Definition at line 113 of file PosteriorParameters.h.

◆ nparameters_fixed()

size_t cbl::statistics::PosteriorParameters::nparameters_fixed ( ) const
inlineoverridevirtual

return the number of fixed parameters

Returns
the number of fixed parameters

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 192 of file PosteriorParameters.h.

◆ nparameters_free()

size_t cbl::statistics::PosteriorParameters::nparameters_free ( ) const
inlineoverridevirtual

return the number of free parameters

Returns
the number of free parameters

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 161 of file PosteriorParameters.h.

◆ parameter_chain_values()

std::vector< double > cbl::statistics::PosteriorParameters::parameter_chain_values ( const int  param,
const int  start = 0,
const int  thin = 1 
) const
virtual

return all the chain values for a parameter

Parameters
paramthe parameter index
startthe starting position
thinnumber of jumped indexes in the chain
Returns
the chain value

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 592 of file PosteriorParameters.cpp.

◆ parameter_covariance() [1/2]

std::vector< std::vector< double > > cbl::statistics::PosteriorParameters::parameter_covariance ( ) const
virtual

return the protected member m_parameter_covariance

Returns
vector containing the parameter covariance matrix

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 284 of file PosteriorParameters.cpp.

◆ parameter_covariance() [2/2]

double cbl::statistics::PosteriorParameters::parameter_covariance ( const int  i,
const int  j 
) const
virtual

return the protected member m_parameter_covariance

Parameters
iindex i
jindex j
Returns
the value of the parameter covariance at i,j

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 275 of file PosteriorParameters.cpp.

◆ posterior_distribution()

std::shared_ptr<PosteriorDistribution> cbl::statistics::PosteriorParameters::posterior_distribution ( const int  par) const
inlinevirtual

get the posterior distribution for the chosen parameter

Parameters
parthe index of the parameter
Returns
the protected member m_posterior_distribution[par]

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 416 of file PosteriorParameters.h.

◆ prior()

std::shared_ptr< cbl::statistics::Prior > cbl::statistics::PosteriorParameters::prior ( ) const
virtual

get the prior function

Returns
pointer to a class Prior

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 385 of file PosteriorParameters.cpp.

◆ prior_distribution() [1/2]

std::vector<std::shared_ptr<PriorDistribution> > cbl::statistics::PosteriorParameters::prior_distribution ( ) const
inlinevirtual

get the prior distribution for the p-th parameter

Returns
vector containing pointers to the prior distributions

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 318 of file PosteriorParameters.h.

◆ prior_distribution() [2/2]

std::shared_ptr<PriorDistribution> cbl::statistics::PosteriorParameters::prior_distribution ( const int  p) const
inlinevirtual

get the prior distribution for the p-th parameter

Parameters
pthe p-th parameter
Returns
pointer to the prior distribution of the p-th parameter

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 310 of file PosteriorParameters.h.

◆ set_bestfit_values() [1/2]

void cbl::statistics::PosteriorParameters::set_bestfit_values ( const int  start,
const int  thin,
const int  nbins,
const int  seed 
)
overridevirtual

set the protected member m_bestfit_value

Parameters
startthe starting position
thinnumber of jumped indexes in the chain
nbinsthe number of bins
seedseed for random extraction

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 464 of file PosteriorParameters.cpp.

◆ set_bestfit_values() [2/2]

void cbl::statistics::PosteriorParameters::set_bestfit_values ( const std::vector< double >  bestfit_value)
virtual

set the protected member m_bestfit_value

Parameters
bestfit_valueparameter bestfit values

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 449 of file PosteriorParameters.cpp.

◆ set_chain()

void cbl::statistics::PosteriorParameters::set_chain ( const size_t  size,
const size_t  nwalkers 
)
virtual

set the chain

Parameters
sizethe chain lenght
nwalkersthe number of parallel walkers

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 541 of file PosteriorParameters.cpp.

◆ set_chain_value()

void cbl::statistics::PosteriorParameters::set_chain_value ( const int  param,
const int  pos,
const int  ww,
const double  value 
)
inlinevirtual

set the private member m_chain_value at the pp-th step for the ww-th step

Parameters
paramthe parameter index
posthe position in the chain
wwthe walker index
valuethe chain value

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 619 of file PosteriorParameters.h.

◆ set_chain_values() [1/2]

void cbl::statistics::PosteriorParameters::set_chain_values ( const std::vector< std::vector< double >>  values,
const int  nwalkers 
)
virtual

set the chain values

Parameters
valuesthe input chain values
nwalkersthe number of parallel walkers

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 607 of file PosteriorParameters.cpp.

◆ set_chain_values() [2/2]

void cbl::statistics::PosteriorParameters::set_chain_values ( const std::vector< std::vector< std::vector< double >>>  values)
virtual

set the chain values

Parameters
valuesthe input chain values

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 626 of file PosteriorParameters.cpp.

◆ set_parameter_covariance()

void cbl::statistics::PosteriorParameters::set_parameter_covariance ( const int  start = 0,
const int  thin = 1 
)
virtual

set the internal method m_parameter_covariance

Parameters
startthe starting position
thinnumber of jumped indexes in the chain

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 261 of file PosteriorParameters.cpp.

◆ set_parameters()

void cbl::statistics::PosteriorParameters::set_parameters ( const size_t  nparameters,
const std::vector< std::shared_ptr< PriorDistribution >>  priorDistributions,
std::vector< ParameterType parameterTypes,
std::vector< std::string >  parameterNames 
)
overridevirtual

set the parameter

Parameters
nparametersthe number of parameters
priorDistributionsvector containing the parameter priors
parameterTypesthe parameter types
parameterNamesthe parameter names

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 202 of file PosteriorParameters.cpp.

◆ set_posterior_distribution()

void cbl::statistics::PosteriorParameters::set_posterior_distribution ( const int  start,
const int  thin,
const int  nbins,
const int  seed = 34121,
const std::vector< double >  weight = {} 
)
virtual

set the posterior distribution from the chains

Parameters
startthe starting position
thinnumber of jumped indexes in the chain
nbinsthe number of bins
seedseed for random extraction
weightchain weight

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 723 of file PosteriorParameters.cpp.

◆ set_prior_distribution() [1/2]

void cbl::statistics::PosteriorParameters::set_prior_distribution ( const int  p,
const std::shared_ptr< PriorDistribution priorDistribution 
)
virtual

set the prior distribution for the p-th parameter

Parameters
pthe p-th parameter
priorDistributionthe prior distribution

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 293 of file PosteriorParameters.cpp.

◆ set_prior_distribution() [2/2]

void cbl::statistics::PosteriorParameters::set_prior_distribution ( const std::vector< std::shared_ptr< PriorDistribution >>  priorDistribution)
virtual

set the prior distributions for the parameters

Parameters
priorDistributionthe prior distributions

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 319 of file PosteriorParameters.cpp.

◆ set_prior_distribution_seed()

void cbl::statistics::PosteriorParameters::set_prior_distribution_seed ( const std::shared_ptr< random::UniformRandomNumbers_Int ran_generator)
virtual

set the prior distributions for the parameters

Parameters
ran_generatorthe random generator

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 375 of file PosteriorParameters.cpp.

◆ show_results()

void cbl::statistics::PosteriorParameters::show_results ( const int  start,
const int  thin,
const int  nbins,
const int  seed = 34121,
const bool  show_mode = false,
const int  ns = -1,
const int  nb = -1,
const std::vector< double >  weight = {} 
)
virtual

show the results on the standard output

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 starting position
thinnumber of jumped indexes in the chain
nbinsthe number of bins
seedseed for random extraction
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
weightchain weight

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 746 of file PosteriorParameters.cpp.

◆ status() [1/2]

vector< string > cbl::statistics::PosteriorParameters::status ( ) const
virtual

return all the model parameter status

Returns
vector containing all the parameter statuss

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 76 of file PosteriorParameters.cpp.

◆ status() [2/2]

string cbl::statistics::PosteriorParameters::status ( const int  p) const
virtual

return the model parameter status

Parameters
pthe index of the parameter
Returns
the parameter status

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 45 of file PosteriorParameters.cpp.

◆ write_results()

void cbl::statistics::PosteriorParameters::write_results ( const std::string  dir,
const std::string  file,
const int  start,
const int  thin,
const int  nbins,
const int  seed = 34121,
const bool  compute_mode = false,
const int  ns = -1,
const int  nb = -1,
const std::vector< double >  weight = {} 
)
virtual

store the results to file

this function stores to file the posterior mean, the posterior standard deviation, the posterior median, 16th and 84th 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
dirname of the output folder
filename of the output file
startthe starting position
thinnumber of jumped indexes in the chain
nbinsthe number of bins
seedseed for random extraction
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
weightchain weight

Reimplemented from cbl::statistics::ModelParameters.

Definition at line 834 of file PosteriorParameters.cpp.


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