![]() |
CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
|
The class CombinedPosterior. More...
#include "Headers/CombinedPosterior.h"
Public Member Functions | |
void | set_log_posterior (const std::vector< double > logpostA, const std::vector< double > logpostB) |
set the internal values of m_log_posterior as the concatenation of the logposterior vectors of two cosmological chains More... | |
void | set_parameters (const std::vector< std::vector< double >> parametersA, const std::vector< std::vector< double >> parametersB) |
set the internal values of m_parameters as the concatenation of the parameters vectors of two cosmological chains More... | |
void | set_weight (const std::vector< double > weightsA, const std::vector< double > weightsB) |
set the internal values of m_weight as the concatenation of the weights vectors of two MCMC chains More... | |
void | m_check_repeated_par (int dummy_Nposteriors, std::vector< std::shared_ptr< Posterior >> posteriors, const std::vector< std::string > repeated_par) |
check the repeated parameters More... | |
void | m_check_common_repeated_par (int dummy_Nposteriors, std::vector< std::shared_ptr< Posterior >> posteriors, std::vector< std::string > repeated_par, std::vector< std::vector< std::vector< int >>> common_repeated_par) |
check the common repeated parameters More... | |
void | m_add_prior (bool par_is_repeated, std::vector< std::shared_ptr< Posterior >> posteriors, const int N, const int k, std::vector< std::shared_ptr< cbl::statistics::PriorDistribution >> &prior_distributions) |
add a prior More... | |
void | m_set_common_repeated_par (std::vector< std::shared_ptr< Posterior >> posteriors, const bool is_in_parnames, const int N, const int k, const std::vector< std::string > repeated_par, const std::vector< std::vector< std::vector< int >>> common_repeated_par, std::vector< std::shared_ptr< cbl::statistics::PriorDistribution >> &prior_distributions, std::vector< std::string > ¶meter_names, std::vector< std::string > &original_names, std::vector< ParameterType > ¶meter_types) |
set a common repeated parameter More... | |
void | m_set_repeated_par (std::vector< std::shared_ptr< Posterior >> posteriors, const bool is_in_parnames, const int N, const int k, const std::vector< std::string > repeated_par, const std::vector< std::vector< std::vector< int >>> common_repeated_par, std::vector< std::shared_ptr< cbl::statistics::PriorDistribution >> &prior_distributions, std::vector< std::string > ¶meter_names, std::vector< std::string > &original_names, std::vector< ParameterType > ¶meter_types) |
set a repeated parameter More... | |
void | m_set_parameters_priors (std::vector< std::shared_ptr< Posterior >> posteriors, std::vector< std::string > repeated_par={}, const std::vector< std::vector< std::vector< int >>> common_repeated_par={}) |
set the parameters and the priors More... | |
void | m_set_independent_probes () |
set all the internal variables needed for modelling independent probes | |
void | importance_sampling (const int distNum, const double cell_size, const double rMAX, const double cut_sigma=-1) |
do the importance sampling for two Posterior objects, which has been read externally More... | |
void | importance_sampling (const std::string output_path, const std::string model_nameA, const std::string model_nameB, const std::vector< double > start, const int chain_size, const int nwalkers, const int burn_in, const int thin) |
do the importance sampling for two Posterior objects, that are constructed from dataset and models 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) |
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, const std::string input_dir, const std::string input_file, const int seed) |
initialize the chains reading the input values 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... | |
double | operator() (std::vector< double > &pp) const |
evaluate the un-normalized posterior as the product of the un-normalized entry posteriors. For N Posterior objects, at a given point \( \vec{\theta_i} \) in the parameter space: More... | |
double | log (std::vector< double > &pp) const |
evaluate the logarithm of the un-normalized posterior as the sum of all the logarithm of the un-normalized posteriors of the N entry objects. Considering a given point \( \vec{\theta_i} \) in the parameter space: 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 | 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 | 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 combination on the screen. 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_model_from_chain (const std::string output_dir, const std::string output_file, const int start, const int thin, const std::vector< double > xx={}, const std::vector< double > yy={}) |
write the model computing 16th, 50th and 84th percentiles from the MCMC More... | |
void | write_maximization_results (const std::string dir_output, const std::string file) |
write maximization results on a file More... | |
std::vector< std::vector< double > > | read (const std::string path, const std::string filename) |
read an entire 2d table data (MCMC chain) of unknown dimension More... | |
Constructors/destructors | |
CombinedPosterior ()=default | |
default constructor | |
CombinedPosterior (const std::vector< std::shared_ptr< Posterior >> posteriors, std::vector< std::string > repeated_par={}, const std::vector< std::vector< std::vector< int >>> common_repeated_par={}) | |
Constructor used to set the modelling of statistically independent probes, or used to perform importance sampling. More... | |
CombinedPosterior (const std::vector< std::vector< std::shared_ptr< Posterior >>> posteriors, const std::vector< std::shared_ptr< data::CovarianceMatrix >> covariance, const std::vector< cbl::statistics::LikelihoodType > likelihood_types, const std::vector< std::string > repeated_par={}, const std::vector< std::vector< std::vector< int >>> common_repeated_par={}, const std::vector< std::shared_ptr< cbl::cosmology::SuperSampleCovariance >> SSC={}) | |
Constructor used to set the modelling of statistically dependent probes. For each vector of cbl::statistics::Posterior objects in the posteriors argument, a cbl::data::CovarianceMatrix object must be defined. The probes in each vector within posteriors are described by the same likelihood function. The final log-likelihood is given by the sum of the logarithms of each likelihood describing a set of dependent probes. More... | |
~CombinedPosterior ()=default | |
default destructor | |
![]() | |
std::shared_ptr< ModelParameters > | parameters () 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::Prior > | get_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::Model > | get_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... | |
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 | |
![]() | |
std::shared_ptr< ModelParameters > | parameters () const |
return the likelihood parameters More... | |
double | operator() (std::vector< double > &pp) const |
evaluate the likelihood More... | |
double | log (std::vector< double > ¶meter) 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< Model > | get_m_model () |
get the values of the internal variable m_model More... | |
std::shared_ptr< ModelParameters > | get_model_parameters () |
get the values of the internal variable m_model More... | |
std::shared_ptr< data::Data > | get_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 | |
Private Attributes | |
STR_DependentProbes_data_model | m_data_model |
the container of parameters for statistically dependent probes | |
std::vector< std::shared_ptr< Posterior > > | m_posteriors |
shared pointers vector containing input Posterior objects | |
int | m_Nposteriors |
number of posteriors | |
std::vector< std::vector< int > > | m_parameter_indexes |
indexes of the parameters, for each posterior object, in the total parameter vector | |
std::vector< std::vector< int > > | m_parameter_indexes2 |
indexes of the parameters, for each posterior object, in the total parameter vector | |
std::vector< int > | m_cosmoPar_indexes |
indexes of the cosmological parameters | |
std::vector< std::vector< std::string > > | m_parameter_names |
names of the parameters, for each posterior object | |
bool | impsampling = false |
importance sampling mode | |
std::vector< std::shared_ptr< void > > | m_likelihood_inputs |
likelihood inputs | |
std::vector< LogLikelihood_function > | m_log_likelihood_functions |
log-likelihood functions | |
std::vector< Likelihood_function > | m_likelihood_functions |
likelihood functions | |
std::vector< LogLikelihood_function > | m_likelihood_functions_grid |
likelihood functions on a grid | |
std::vector< LogLikelihood_function > | m_log_likelihood_functions_grid |
log-likelihood functions on a grid | |
std::vector< std::shared_ptr< data::Data > > | m_datasets |
data containers | |
std::vector< std::shared_ptr< Model > > | m_models |
models to test | |
std::vector< bool > | m_use_grid |
use_grid vector | |
Additional Inherited Members | |
![]() | |
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 | |
![]() | |
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... | |
![]() | |
std::shared_ptr< cbl::statistics::Prior > | m_prior |
the prior distribution | |
std::shared_ptr< statistics::Model > | m_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_Int > | m_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 | |
![]() | |
std::shared_ptr< data::Data > | m_data |
data containers | |
std::shared_ptr< Model > | m_model |
model to test | |
std::shared_ptr< void > | m_likelihood_inputs |
likelihood inputs | |
std::shared_ptr< ModelParameters > | m_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 | |
The class CombinedPosterior.
This class is used to define the distribution
Definition at line 105 of file CombinedPosterior.h.
cbl::statistics::CombinedPosterior::CombinedPosterior | ( | const std::vector< std::shared_ptr< Posterior >> | posteriors, |
std::vector< std::string > | repeated_par = {} , |
||
const std::vector< std::vector< std::vector< int >>> | common_repeated_par = {} |
||
) |
Constructor used to set the modelling of statistically independent probes, or used to perform importance sampling.
posteriors | pointers to Posterior objects |
repeated_par | parameters shared by different probes, for which the user wants different posteriors for each probe. For example, if the probes A and B depend on the same parameter \(p\), whose identification string is "par", then "par" must be given in input if the user desires to have different posteriors of \(p\) for the two probes A and B. Every probe depending on \(p\) will provide a different posterior on \(p\). This is useful in the case of astrophysical parameters, e.g. the parameters describing galaxy cluster profiles. |
common_repeated_par | for each argument in repeated_par, a vector of vectors of integers can be defined here. Such vectors define the sets of probes for which the same priors and posteriors are provided for the same parameters. For example, let us consider the probes {A1, A2, A3, A4}, provided in the posteriors parameter. All the probes (A1, A2, A3, A4) depend on the parameter \(p\), whose identification string is "par", and we set repeated_par = {"par"}. If we want the pair of probes {A1, A2}, {A3, A4}, to provide a different posterior on \(p\), we must set common_repeated_par = { { {0,1}, {2,3} } }. This is useful when different probes in the same bin provide constraints on the same parameters. If common_repeated_par is not provided, every probe depending on \(p\) will provide a different posterior on \(p\). If two parameters are provided in repeated_par, e.g. repeated_par = {"par1", "par2"}, and only "par2" must be shared by more than one probe, then leave blank the vector of vectors corresponding to "par1" in common_repeated_par. |
Definition at line 47 of file CombinedPosterior.cpp.
cbl::statistics::CombinedPosterior::CombinedPosterior | ( | const std::vector< std::vector< std::shared_ptr< Posterior >>> | posteriors, |
const std::vector< std::shared_ptr< data::CovarianceMatrix >> | covariance, | ||
const std::vector< cbl::statistics::LikelihoodType > | likelihood_types, | ||
const std::vector< std::string > | repeated_par = {} , |
||
const std::vector< std::vector< std::vector< int >>> | common_repeated_par = {} , |
||
const std::vector< std::shared_ptr< cbl::cosmology::SuperSampleCovariance >> | SSC = {} |
||
) |
Constructor used to set the modelling of statistically dependent probes. For each vector of cbl::statistics::Posterior objects in the posteriors argument, a cbl::data::CovarianceMatrix object must be defined. The probes in each vector within posteriors are described by the same likelihood function. The final log-likelihood is given by the sum of the logarithms of each likelihood describing a set of dependent probes.
If sets of probes are described by user-defined likelihoods, then the cbl::data::CovarianceMatrix objects must not be provided for such sets.
posteriors | vector of vectors of pointers to Posterior objects. In each vector, a set of probes is contained, with a covariance matrix defined by the corresponding cbl::data::CovarianceMatrix object given as input in the second argument of this constructor |
covariance | objects defining the covariance matrices for the Posterior objects given in input |
likelihood_types | likelihood types describing each set of probes |
repeated_par | parameters shared by different probes, for which the user wants different posteriors for each probe. For example, if the probes A and B depend on the same parameter \(p\), whose identification string is "par", then "par" must be given in input if the user desires to have different posteriors of \(p\) for the two probes A and B. This is useful in the case of astrophysical parameters, e.g. the parameters describing galaxy cluster profiles. |
common_repeated_par | for each argument in repeated_par, a vector of vectors of integers can be defined here. Such vectors define the sets of probes for which the same priors and posteriors are provided for the same parameters. For example, let us consider the two sets of probes {A1, A2, A3} and {B1, B2, B3}, provided in the posteriors parameter. All the probes (A1, A2, A3, B1, B2, B3) depend on the parameter \(p\), whose identification string is "par", and we set repeated_par = {"par"}. If we want each pair of probes {A1, B1}, {A2, B2}, {A3, B3}, to provide a different posterior on \(p\), we must set common_repeated_par = { { {0,3}, {1,4}, {2,5} } }. This is useful when different probes in the same bin provide constraints on the same parameters. If common_repeated_par is not provided, every probe depending on \(p\) will provide a different posterior on \(p\). If two parameters are provided in repeated_par, e.g. repeated_par = {"par1", "par2"}, and only "par2" must be shared by more than one probe, then leave blank the vector of vectors corresponding to "par1" in common_repeated_par. |
SSC | vector of pointers to cbl::cosmology::SuperSampleCovariance objects, for the computation of the \(S_{ij}\) matrices. If, for example, the sets of probes A, B, C, are considered, and only for A and C the super-sample covariance must be computed, then set the second element of SSC equal to NULL. |
Definition at line 81 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::importance_sampling | ( | const int | distNum, |
const double | cell_size, | ||
const double | rMAX, | ||
const double | cut_sigma = -1 |
||
) |
do the importance sampling for two Posterior objects, which has been read externally
distNum | number of points for the averaging |
cell_size | cell dimension |
rMAX | searching radius |
cut_sigma | for the weight's distribution cut |
Definition at line 662 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::importance_sampling | ( | const std::string | output_path, |
const std::string | model_nameA, | ||
const std::string | model_nameB, | ||
const std::vector< double > | start, | ||
const int | chain_size, | ||
const int | nwalkers, | ||
const int | burn_in, | ||
const int | thin | ||
) |
do the importance sampling for two Posterior objects, that are constructed from dataset and models
output_path | the path where the chains will be stored |
model_nameA | the name of the first model |
model_nameB | the name of the second model |
start | std::vector containing initial values for the posterior maximization |
chain_size | the size of the chains that has to be ran |
nwalkers | the number of parallal walkers for the chain |
burn_in | the minimum chain position to be written |
thin | the step used for dilution on screen |
Definition at line 713 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::initialize_chains | ( | const int | chain_size, |
const int | n_walkers | ||
) |
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
chain_size | the chain lenght |
n_walkers | the number of parallel chains |
Definition at line 779 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::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
chain_size | the chain lenght |
n_walkers | the number of parallel chains |
radius | radius of the ball in parameter space |
start | std::vector containing initial values for the posterior maximization |
max_iter | the maximum number of iterations |
tol | the tolerance in finding convergence |
epsilon | the simplex side |
Definition at line 768 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::initialize_chains | ( | const int | chain_size, |
const int | n_walkers, | ||
const std::string | input_dir, | ||
const std::string | input_file, | ||
const int | seed | ||
) |
initialize the chains reading the input values from an input file
chain_size | the chain lenght |
n_walkers | the number of parallel chains |
input_dir | input directory |
input_file | input file |
seed | the seed |
Definition at line 789 of file CombinedPosterior.cpp.
double cbl::statistics::CombinedPosterior::log | ( | std::vector< double > & | pp | ) | const |
evaluate the logarithm of the un-normalized posterior as the sum of all the logarithm of the un-normalized posteriors of the N entry objects. Considering a given point \( \vec{\theta_i} \) in the parameter space:
\[ P(\vec{\theta_i} | \vec{d}) = \sum_{j=1}^{N} \log{\mathcal{L_j}(\vec{d}|\vec{\theta_i})} + \log{Pr(\vec{\theta_i})} \]
where \(P\) is the combined posterior, \(\mathcal{L_j}(\vec{d}|\vec{\theta})\) is the j-th likelihood function and \(Pr(\vec{\theta})\) is the prior distribution.
pp | the parameters |
Definition at line 914 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::m_add_prior | ( | bool | par_is_repeated, |
std::vector< std::shared_ptr< Posterior >> | posteriors, | ||
const int | N, | ||
const int | k, | ||
std::vector< std::shared_ptr< cbl::statistics::PriorDistribution >> & | prior_distributions | ||
) |
add a prior
par_is_repeated | if true, the considered parameter is repeated |
posteriors | pointer to the Posterior objects |
N | index of the Posterior object |
k | index of the parameter of the N-th Posterior |
prior_distributions | vector containing all the Posterior objects |
Definition at line 337 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::m_check_common_repeated_par | ( | int | dummy_Nposteriors, |
std::vector< std::shared_ptr< Posterior >> | posteriors, | ||
std::vector< std::string > | repeated_par, | ||
std::vector< std::vector< std::vector< int >>> | common_repeated_par | ||
) |
check the common repeated parameters
dummy_Nposteriors | number of posterior objects |
posteriors | vector of input Posterior objects |
repeated_par | repeated parameters |
common_repeated_par | indices of the common repeated parameters |
Definition at line 309 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::m_check_repeated_par | ( | int | dummy_Nposteriors, |
std::vector< std::shared_ptr< Posterior >> | posteriors, | ||
const std::vector< std::string > | repeated_par | ||
) |
check the repeated parameters
dummy_Nposteriors | number of posterior objects |
posteriors | pointers to Posterior objects |
repeated_par | repeated parameters |
Definition at line 288 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::m_set_common_repeated_par | ( | std::vector< std::shared_ptr< Posterior >> | posteriors, |
const bool | is_in_parnames, | ||
const int | N, | ||
const int | k, | ||
const std::vector< std::string > | repeated_par, | ||
const std::vector< std::vector< std::vector< int >>> | common_repeated_par, | ||
std::vector< std::shared_ptr< cbl::statistics::PriorDistribution >> & | prior_distributions, | ||
std::vector< std::string > & | parameter_names, | ||
std::vector< std::string > & | original_names, | ||
std::vector< ParameterType > & | parameter_types | ||
) |
set a common repeated parameter
posteriors | pointer to the Posterior objects |
is_in_parnames | if true, the parameter is already in the vector of parameter names |
N | index of the Posterior object |
k | index of the parameter of the N-th Posterior |
repeated_par | repeated parameter names |
common_repeated_par | indices of the common repeated parameters |
prior_distributions | vector containing all the Posterior objects |
parameter_names | vector containing all the parameter names |
original_names | vector of original parameter names |
parameter_types | vector containing all the parameter types |
Definition at line 358 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::m_set_parameters_priors | ( | std::vector< std::shared_ptr< Posterior >> | posteriors, |
std::vector< std::string > | repeated_par = {} , |
||
const std::vector< std::vector< std::vector< int >>> | common_repeated_par = {} |
||
) |
set the parameters and the priors
posteriors | pointers to Posterior objects |
repeated_par | parameters shared by different probes, for which the user wants different posteriors for each probe. For example, if the probes A and B depend on the same parameter \(p\), whose identification string is "par", then "par" must be given in input if the user desires to have different posteriors of \(p\) for the two probes A and B. This is useful in the case of astrophysical parameters, e.g. the parameters describing galaxy cluster profiles. |
common_repeated_par | for each argument in repeated_par, a vector of vectors of integers can be defined here. Such vectors define the sets of probes for which the same priors and posteriors are provided for the same parameters. For example, let us consider the two sets of probes {A1, A2, A3} and {B1, B2, B3}, provided in the posteriors parameter. All the probes (A1, A2, A3, B1, B2, B3) depend on the parameter \(p\), whose identification string is "par", and we set repeated_par = {"par"}. If we want each pair of probes {A1, B1}, {A2, B2}, {A3, B3}, to provide a different posterior on \(p\), we must set common_repeated_par = { { {0,1}, {2,3}, {4,5} } }. This is useful when different probes in the same bin provide constraints on the same parameters. If common_repeated_par is not provided, every probe depending on \(p\) will provide a different posterior on \(p\). If two parameters are provided in repeated_par, e.g. repeated_par = {"par1", "par2"}, and only "par2" must be shared by more than one probe, then leave blank the vector of vectors corresponding to "par1" in common_repeated_par. |
Definition at line 430 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::m_set_repeated_par | ( | std::vector< std::shared_ptr< Posterior >> | posteriors, |
const bool | is_in_parnames, | ||
const int | N, | ||
const int | k, | ||
const std::vector< std::string > | repeated_par, | ||
const std::vector< std::vector< std::vector< int >>> | common_repeated_par, | ||
std::vector< std::shared_ptr< cbl::statistics::PriorDistribution >> & | prior_distributions, | ||
std::vector< std::string > & | parameter_names, | ||
std::vector< std::string > & | original_names, | ||
std::vector< ParameterType > & | parameter_types | ||
) |
set a repeated parameter
posteriors | pointer to the Posterior objects |
is_in_parnames | if true, the parameter is already in the vector of parameter names |
N | index of the Posterior object |
k | index of the parameter of the N-th Posterior |
repeated_par | repeated parameter names |
common_repeated_par | indices of the common repeated parameters |
prior_distributions | vector containing all the Posterior objects |
parameter_names | vector containing all the parameter names |
original_names | vector of original parameter names |
parameter_types | vector containing all the parameter types |
Definition at line 410 of file CombinedPosterior.cpp.
|
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)
start | std::vector containing initial values for the posterior maximization |
parameter_limits | limits for the parameters |
max_iter | the maximum number of iterations |
tol | the tolerance in finding convergence |
epsilon | the simplex side |
Definition at line 663 of file CombinedPosterior.h.
void cbl::statistics::CombinedPosterior::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
start | std::vector containing initial values for the posterior maximization |
max_iter | the maximum number of iterations |
tol | the tolerance in finding convergence |
epsilon | the simplex side |
Definition at line 946 of file CombinedPosterior.cpp.
double cbl::statistics::CombinedPosterior::operator() | ( | std::vector< double > & | pp | ) | const |
evaluate the un-normalized posterior as the product of the un-normalized entry posteriors. For N Posterior objects, at a given point \( \vec{\theta_i} \) in the parameter space:
\[ P(\vec{\theta_i} | \vec{d}) = \prod_{j=1}^{N} \mathcal{L_j}(\vec{d}|\vec{\theta_i}) \cdot Pr(\vec{\theta_i}) \]
where \(P\) is the combined posterior, \(\mathcal{L_j}(\vec{d}|\vec{\theta})\) is the j-th likelihood function and \(Pr(\vec{\theta_i})\) is the prior distribution.
pp | the parameters |
Definition at line 885 of file CombinedPosterior.cpp.
std::vector< std::vector< double > > cbl::statistics::CombinedPosterior::read | ( | const std::string | path, |
const std::string | filename | ||
) |
read an entire 2d table data (MCMC chain) of unknown dimension
path | the path in which the file is stored |
filename | the name of the file that has to be read |
Definition at line 1239 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::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)
aa | the parameter of the \(g(z)\) distribution |
parallel | false \(\rightarrow\) non-parallel sampler; true \(\rightarrow\) parallel sampler |
outputFile | output 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 |
start | the minimum chain position used to compute the median |
thin | the step used for chain dilution |
nbins | the number of bins to estimate the posterior distribution, used to assess its properties |
Definition at line 833 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::set_log_posterior | ( | const std::vector< double > | logpostA, |
const std::vector< double > | logpostB | ||
) |
set the internal values of m_log_posterior as the concatenation of the logposterior vectors of two cosmological chains
logpostA | the first logposterior distribution |
logpostB | the second logposterior distribution |
Definition at line 641 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::set_parameters | ( | const std::vector< std::vector< double >> | parametersA, |
const std::vector< std::vector< double >> | parametersB | ||
) |
set the internal values of m_parameters as the concatenation of the parameters vectors of two cosmological chains
parametersA | the first parameters vector |
parametersB | the second parameter vector |
Definition at line 626 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::set_weight | ( | const std::vector< double > | weightsA, |
const std::vector< double > | weightsB | ||
) |
set the internal values of m_weight as the concatenation of the weights vectors of two MCMC chains
weightsA | the first weights vector |
weightsB | the second weights vector |
Definition at line 651 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::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 combination on the screen.
In the case of the sum of logPosteriors:
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).
In the case of importance sampling:
The weighted average of the chains parameters are shown.
start | the minimum chain position to be written |
thin | the step used for dilution on screen |
nbins | the number of bins to estimate the posterior distribution, used to assess its properties |
show_mode | true \(\rightarrow\) show the posterior mode; false \(\rightarrow\) do not show the posterior mode |
ns | number of samples used to estimate the covariance matrix |
nb | number of data measurements, e.g. the bins of the dataset |
Definition at line 1007 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::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
output_dir | the output directory |
output_file | the output file |
start | the minimum chain position to be written |
thin | the step used for dilution |
is_FITS_format | true \(\rightarrow\) the format of the input file is FITS; false \(\rightarrow\) the format of the input file is ASCII |
prec | decimal precision |
ww | number of characters to be used as field width |
Definition at line 1059 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::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
output_dir | the output directory |
output_file | the output file |
start | the minimum chain position to be written |
thin | the step used for dilution |
prec | decimal precision |
ww | number of characters to be used as field width |
Definition at line 1071 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::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
output_dir | the output directory |
output_file | the output file |
start | the minimum chain position to be written |
thin | the step used for dilution |
Definition at line 1123 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::write_maximization_results | ( | const std::string | dir_output, |
const std::string | file | ||
) |
write maximization results on a file
dir_output | the output directory |
file | the name of the output file to be written |
Definition at line 1206 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::write_model_from_chain | ( | const std::string | output_dir, |
const std::string | output_file, | ||
const int | start, | ||
const int | thin, | ||
const std::vector< double > | xx = {} , |
||
const std::vector< double > | yy = {} |
||
) |
write the model computing 16th, 50th and 84th percentiles from the MCMC
output_dir | the output directory |
output_file | the output file |
start | the minimum chain position to be written |
thin | the step used for dilution on screen |
xx | x points where the model is computed. If not provided, the x points set for the MCMC are considered. |
yy | y points where the model is computed. If not provided, the y points set for the MCMC are considered. |
Definition at line 1174 of file CombinedPosterior.cpp.
void cbl::statistics::CombinedPosterior::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).
output_dir | the output directory |
root_file | the root of the output file to be written |
start | the minimum chain position to be written |
thin | the step used for dilution on screen |
nbins | the number of bins to estimate the posterior distribution, used to assess its properties |
fits | false \(\rightarrow\) ascii file; true \(\rightarrow\) fits file |
compute_mode | true \(\rightarrow\) compute the posterior mode; false \(\rightarrow\) do not compute the posterior mode |
ns | number of samples used to estimate the covariance matrix |
nb | number of data measurements, e.g. the bins of the dataset |
Definition at line 1021 of file CombinedPosterior.cpp.