![]() |
CosmoBolognaLib
Free Software C++/Python libraries for cosmological calculations
|
#include "Headers/Modelling.h"
Public Member Functions | |
void | m_set_posterior (const int seed) |
set the interal variable m_posterior More... | |
Constructors/destructors | |
Modelling ()=default | |
default constuctor | |
virtual | ~Modelling ()=default |
default destructor | |
Member functions used to get the protected members | |
std::shared_ptr< data::Data > | data () |
return the dataset More... | |
std::shared_ptr< data::Data > | data_fit () |
return the dataset More... | |
std::shared_ptr< statistics::Likelihood > | likelihood () |
return the likelihood parameters More... | |
std::shared_ptr< statistics::Posterior > | posterior () |
return the posterior parameters More... | |
std::shared_ptr< statistics::ModelParameters > | likelihood_parameters () |
return the likelihood parameters More... | |
std::shared_ptr< statistics::ModelParameters > | posterior_parameters () |
return the posterior parameters More... | |
virtual void | set_parameter_from_string (const std::string parameter, const double value) |
set the value of a parameter providing its name string More... | |
virtual double | get_parameter_from_string (const std::string parameter) const |
get the value of a parameter providing its name string More... | |
std::shared_ptr< statistics::PriorDistribution > | get_prior (const int i) |
get the internal variable m_parameter_priors More... | |
std::shared_ptr< statistics::Model > | get_response_function () |
return the response function used to compute the super-sample covariance More... | |
Member functions used to set internal parameters | |
void | reset_fit_range () |
reset the fit range More... | |
void | set_fit_range (const double xmin, const double xmax) |
set the fit range More... | |
void | set_fit_range (const double xmin, const double xmax, const double ymin, const double ymax) |
set the fit range More... | |
void | set_data (const std::shared_ptr< data::Data > dataset) |
set the dataset More... | |
void | set_likelihood (const statistics::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 function More... | |
void | set_likelihood (const cbl::statistics::Likelihood_function log_likelihood_function) |
set the likelihood function, given a user-defined log-likelihood function More... | |
Member functions used to manage likelihood/posterior | |
distributions | |
void | maximize_likelihood (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 stores them in the model More... | |
void | maximize_posterior (const std::vector< double > start, const unsigned int max_iter=10000, const double tol=1.e-6, const double epsilon=1.e-3, const int seed=666) |
function that maximizes the posterior, finds the best-fit parameters and stores them in the model More... | |
void | sample_posterior (const int chain_size, const int nwalkers, const int seed=666, const double aa=2, const bool parallel=true) |
sample the posterior, initializing the chains by drawing from the prior distributions More... | |
void | sample_posterior (const int chain_size, const int nwalkers, 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, const int seed=666, const double aa=2, const bool parallel=true) |
sample the posterior, initializing the chains in a ball around the posterior best-fit parameters values More... | |
void | sample_posterior (const int chain_size, const int nwalkers, std::vector< double > &value, const double radius, const int seed=666, const double aa=2, const bool parallel=true) |
sample the posterior, initializing the chains by drawing from the prior distributions More... | |
void | sample_posterior (const int chain_size, const std::vector< std::vector< double >> chain_value, const int seed=666, const double aa=2, const bool parallel=true) |
sample the posterior, initializing the chains with input values More... | |
void | sample_posterior (const int chain_size, const int nwalkers, const std::string input_dir, const std::string input_file, const int seed=666, const double aa=2, const bool parallel=true) |
sample the posterior, initializing the chains reading the input values from an input file More... | |
void | importance_sampling (const std::string input_dir, const std::string input_file, const int seed=666, 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) |
perform importance sampling More... | |
Member functions used for Input/Output | |
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 | read_chain (const std::string input_dir, const std::string input_file, const int nwalkers, const std::vector< size_t > columns={}, const int skip_header=1, const bool fits=false) |
read the chains More... | |
void | show_results (const int start=0, const int thin=1, const int nbins=50, const bool show_mode=false, const int ns=-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) |
write the results of the MCMC sampling to file More... | |
virtual void | write_model (const std::string output_dir, const std::string output_file, const std::vector< double > xx, const std::vector< double > parameters) |
write the model at xx for given parameters More... | |
virtual void | write_model (const std::string output_dir, const std::string output_file, const std::vector< double > xx, const std::vector< double > yy, const std::vector< double > parameters) |
write the model at xx, yy for given parameters More... | |
virtual void | write_model_at_bestfit (const std::string output_dir, const std::string output_file, const std::vector< double > xx) |
write the model at xx with best-fit parameters obtained from posterior maximization More... | |
virtual 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, yy with best-fit parameters obtained from likelihood maximization More... | |
virtual void | write_model_from_chains (const std::string output_dir, const std::string output_file, const std::vector< double > xx, const int start=0, const int thin=1) |
write the model at xx computing 16th, 50th and 84th percentiles from the chains More... | |
virtual void | write_model_from_chains (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 chains More... | |
double | reduced_chi2 (const std::vector< double > parameter={}) |
the reduced \(\chi^2\) More... | |
Protected Member Functions | |
void | m_set_prior (std::vector< statistics::PriorDistribution > prior_distribution) |
set the internal variable m_parameter_priors More... | |
void | m_isSet_response () |
check if the response function used to compute the super-sample covariance is set | |
Protected Attributes | |
std::shared_ptr< data::Data > | m_data = NULL |
input data to be modelled | |
bool | m_fit_range = false |
check if fit range has been set | |
std::shared_ptr< data::Data > | m_data_fit |
input data restricted to the range used for the fit | |
std::shared_ptr< statistics::Model > | m_model = NULL |
input model | |
std::shared_ptr< statistics::Model > | m_response_func = NULL |
response function for the computation of the super-sample covariance | |
std::shared_ptr< statistics::Likelihood > | m_likelihood = NULL |
likelihood | |
std::vector< std::shared_ptr< statistics::PriorDistribution > > | m_parameter_priors |
prior | |
std::shared_ptr< statistics::Posterior > | m_posterior = NULL |
posterior | |
The class Modelling.
This file defines the interface of the base class Modelling, used for modelling any kind of measurements
Definition at line 64 of file Modelling.h.
|
inline |
return the dataset
Definition at line 141 of file Modelling.h.
|
inline |
return the dataset
Definition at line 147 of file Modelling.h.
|
inlinevirtual |
get the value of a parameter providing its name string
parameter | parameter name to get |
Definition at line 196 of file Modelling.h.
|
inline |
get the internal variable m_parameter_priors
i | prior index |
Definition at line 206 of file Modelling.h.
|
inline |
return the response function used to compute the super-sample covariance
Definition at line 214 of file Modelling.h.
void cbl::modelling::Modelling::importance_sampling | ( | const std::string | input_dir, |
const std::string | input_file, | ||
const int | seed = 666 , |
||
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 |
||
) |
perform importance sampling
Importance sampling is a convenient technique to join independet datasets.
This function takes in input a chain and computes the likelihood or 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.
input_dir | input directory |
input_file | the input file |
seed | the seed |
column | the columns of the input file to be read |
header_lines_to_skip | the lines to be skipped in the chain file |
is_FITS_format | true \(\rightarrow\) the format of the input file is FITS; false \(\rightarrow\) the format of the input file is ASCII |
apply_to_likelihood | true \(\rightarrow\) the likelihood ratio is used as weight; false \(\rightarrow\) the posterior ratio is used as weight |
Definition at line 240 of file Modelling.cpp.
shared_ptr< statistics::Likelihood > cbl::modelling::Modelling::likelihood | ( | ) |
return the likelihood parameters
Definition at line 72 of file Modelling.cpp.
shared_ptr< cbl::statistics::ModelParameters > cbl::modelling::Modelling::likelihood_parameters | ( | ) |
return the likelihood parameters
Definition at line 98 of file Modelling.cpp.
void cbl::modelling::Modelling::m_set_posterior | ( | const int | seed | ) |
set the interal variable m_posterior
seed | the seed |
Definition at line 58 of file Modelling.cpp.
|
protected |
set the internal variable m_parameter_priors
prior_distribution | vector containing the prior distributions |
Definition at line 46 of file Modelling.cpp.
void cbl::modelling::Modelling::maximize_likelihood | ( | 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 stores 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 likelihood 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 the initial values for the likelihood maximization |
parameter_limits | limits of the parameter space in where the maximum of likelihood will be searched |
max_iter | the maximum number of iterations |
tol | the tolerance in finding convergence |
epsilon | the simplex side |
Definition at line 168 of file Modelling.cpp.
void cbl::modelling::Modelling::maximize_posterior | ( | const std::vector< double > | start, |
const unsigned int | max_iter = 10000 , |
||
const double | tol = 1.e-6 , |
||
const double | epsilon = 1.e-3 , |
||
const int | seed = 666 |
||
) |
function that maximizes the posterior, finds the best-fit parameters and stores 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 | vector containing initial values for the posterior maximization |
max_iter | the maximum number of iterations |
tol | the tolerance to find convergence |
epsilon | the relative fraction of the initial simplex size |
seed | the seed |
Definition at line 177 of file Modelling.cpp.
shared_ptr< cbl::statistics::Posterior > cbl::modelling::Modelling::posterior | ( | ) |
return the posterior parameters
Definition at line 85 of file Modelling.cpp.
shared_ptr< cbl::statistics::ModelParameters > cbl::modelling::Modelling::posterior_parameters | ( | ) |
return the posterior parameters
Definition at line 111 of file Modelling.cpp.
void cbl::modelling::Modelling::read_chain | ( | const std::string | input_dir, |
const std::string | input_file, | ||
const int | nwalkers, | ||
const std::vector< size_t > | columns = {} , |
||
const int | skip_header = 1 , |
||
const bool | fits = false |
||
) |
read the chains
input_dir | the input directory |
input_file | the input file |
nwalkers | the number of parallel chains |
columns | the columns of the input file to be read. |
skip_header | the lines to be skipped in the chain file |
fits | false \(\rightarrow\) ascii file; true \(\rightarrow\) fits file |
Definition at line 259 of file Modelling.cpp.
double cbl::modelling::Modelling::reduced_chi2 | ( | const std::vector< double > | parameter = {} | ) |
the reduced \(\chi^2\)
this function computes the reduced \(\chi^2\), i.e. \(\chi^2/(N_{\rm bins}-N_{\rm parameters})\)
parameter | vector containing the model parameters; if not provided in input, the chi2 is evaluated at best-fit parameter values |
Definition at line 293 of file Modelling.cpp.
|
inline |
reset the fit range
set m_fit_range = false, that means that the fit range is unset
Definition at line 230 of file Modelling.h.
void cbl::modelling::Modelling::sample_posterior | ( | const int | chain_size, |
const int | nwalkers, | ||
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 , |
||
const int | seed = 666 , |
||
const double | aa = 2 , |
||
const bool | parallel = true |
||
) |
sample the posterior, initializing the chains in a ball around the posterior best-fit parameters 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 |
nwalkers | 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 |
seed | the seed |
aa | the parameter of the \(g(z)\) distribution |
parallel | false \(\rightarrow\) non-parallel sampler; true \(\rightarrow\) parallel sampler |
Definition at line 197 of file Modelling.cpp.
void cbl::modelling::Modelling::sample_posterior | ( | const int | chain_size, |
const int | nwalkers, | ||
const int | seed = 666 , |
||
const double | aa = 2 , |
||
const bool | parallel = true |
||
) |
sample the posterior, initializing the chains by drawing from the prior distributions
the starting values of the chain are extracted from the (possibly different) distributions of the priors
chain_size | the chain lenght |
nwalkers | the number of parallel chains |
seed | the seed |
aa | the parameter of the \(g(z)\) distribution |
parallel | false \(\rightarrow\) non-parallel sampler; true \(\rightarrow\) parallel sampler |
Definition at line 187 of file Modelling.cpp.
void cbl::modelling::Modelling::sample_posterior | ( | const int | chain_size, |
const int | nwalkers, | ||
const std::string | input_dir, | ||
const std::string | input_file, | ||
const int | seed = 666 , |
||
const double | aa = 2 , |
||
const bool | parallel = true |
||
) |
sample the posterior, initializing the chains reading the input values 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
chain_size | the chain lenght |
nwalkers | the number of parallel chains |
input_dir | input directory |
input_file | the input file |
seed | the seed |
aa | the parameter of the \(g(z)\) distribution |
parallel | false \(\rightarrow\) non-parallel sampler; true \(\rightarrow\) parallel sampler |
Definition at line 229 of file Modelling.cpp.
void cbl::modelling::Modelling::sample_posterior | ( | const int | chain_size, |
const int | nwalkers, | ||
std::vector< double > & | value, | ||
const double | radius, | ||
const int | seed = 666 , |
||
const double | aa = 2 , |
||
const bool | parallel = true |
||
) |
sample the posterior, initializing the chains by drawing from the prior distributions
the starting values of the chain are extracted from the (possibly different) distributions of the priors
chain_size | the chain lenght |
nwalkers | the number of parallel chains |
value | input values, center of the ball in parameter space |
radius | radius of the ball in parameter space |
seed | the seed |
aa | the parameter of the \(g(z)\) distribution |
parallel | false \(\rightarrow\) non-parallel sampler; true \(\rightarrow\) parallel sampler |
Definition at line 207 of file Modelling.cpp.
void cbl::modelling::Modelling::sample_posterior | ( | const int | chain_size, |
const std::vector< std::vector< double >> | chain_value, | ||
const int | seed = 666 , |
||
const double | aa = 2 , |
||
const bool | parallel = true |
||
) |
sample the posterior, initializing the chains with input values
the starting values of the chain are the elements of the input matrix 'chain_values'
chain_size | the chain lenght |
chain_value | std::vector of size (nwalkers, nparameters) starting values of the chain |
seed | the seed |
aa | the parameter of the \(g(z)\) distribution |
parallel | false \(\rightarrow\) non-parallel sampler; true \(\rightarrow\) parallel sampler |
Definition at line 218 of file Modelling.cpp.
|
inline |
void cbl::modelling::Modelling::set_fit_range | ( | const double | xmin, |
const double | xmax | ||
) |
set the fit range
xmin | minimum x value used for the fit |
xmax | maximum x value used for the fit |
Definition at line 46 of file Modelling1D.cpp.
void cbl::modelling::Modelling::set_fit_range | ( | const double | xmin, |
const double | xmax, | ||
const double | ymin, | ||
const double | ymax | ||
) |
set the fit range
xmin | minimum x value used for the fit |
xmax | maximum x value used for the fit |
ymin | minimum y value used for the fit |
ymax | maximum y value used for the fit |
Definition at line 45 of file Modelling2D.cpp.
void cbl::modelling::Modelling::set_likelihood | ( | const cbl::statistics::Likelihood_function | log_likelihood_function | ) |
set the likelihood function, given a user-defined log-likelihood function
if the input parameter Nres is larger than 0, the inverted covariance matrix is corrected as follows (Hartlap, Simon and Schneider 2006):
\[ \hat{C}^{-1} = \left(1-\frac{n_b+1}{N_{res}-1}\right)C^{-1} \]
where \(n_b\) is the number of bins and \(N_{res}\) is the number of resamplings
log_likelihood_function | the user-defined loglikelihood function (natural logarithm) |
Definition at line 146 of file Modelling.cpp.
void cbl::modelling::Modelling::set_likelihood | ( | const statistics::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 function
if the input parameter Nres is larger than 0, the inverted covariance matrix is corrected as follows (Hartlap, Simon and Schneider 2006):
\[ \hat{C}^{-1} = \left(1-\frac{n_b+1}{N_{res}-1}\right)C^{-1} \]
where \(n_b\) is the number of bins and \(N_{res}\) is the number of resamplings
likelihood_type | the likelihood type, specified with the LikelihoodType object |
x_index | index(s) of the extra info std::vector containing the point(s) where to evaluate the model |
w_index | index of the extra info std::vector containing the data point weight |
prec | the precision required in the inversion of the covariance matrix |
Nres | \(N_{res}\), the number of catalogue resamplings used to estimate the covariance matrix; \(N_{res}=-1\) if the covariance matrix has not been estimated with resampling methods |
Definition at line 124 of file Modelling.cpp.
|
inlinevirtual |
set the value of a parameter providing its name string
parameter | parameter name to set |
value | the new value for the parameter |
Definition at line 186 of file Modelling.h.
void cbl::modelling::Modelling::show_results | ( | const int | start = 0 , |
const int | thin = 1 , |
||
const int | nbins = 50 , |
||
const bool | show_mode = false , |
||
const int | ns = -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 ( \(n_s\), the number of samples used to estimate the covariance matrix) is 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)} \,, \]
where \(n_b\) is number of data measurements (e.g. the bins of the dataset).
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).
start | the minimum chain position to be written |
thin | the step used for dilution on screen |
nbins | the number of bins |
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 |
Definition at line 269 of file Modelling.cpp.
void cbl::modelling::Modelling::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 250 of file Modelling.cpp.
|
virtual |
write the model at xx for given parameters
output_dir | the output directory |
output_file | the output file |
xx | vector of points at which the model is computed, |
parameters | vector containing the input parameters used to compute the model; if this vector is not provided, the model will be computed using the best-fit parameters |
Reimplemented in cbl::modelling::twopt::Modelling_TwoPointCorrelation_wedges.
Definition at line 56 of file Modelling1D.cpp.
|
virtual |
write the model at xx, yy for given parameters
output_dir | the output directory |
output_file | the output file |
xx | vector of points at which the model is computed, first axis |
yy | vector of points at which the model is computed, second axis |
parameters | vector containing the input parameters used to compute the model; if this vector is not provided, the model will be computed using the best-fit parameters |
Definition at line 54 of file Modelling2D.cpp.
|
virtual |
write the model at xx with best-fit parameters obtained from posterior maximization
output_dir | the output directory |
output_file | the output file |
xx | vector of points at which the model is computed |
Reimplemented in cbl::modelling::twopt::Modelling_TwoPointCorrelation_wedges.
Definition at line 68 of file Modelling1D.cpp.
|
virtual |
write the model at xx, yy with best-fit parameters obtained from likelihood maximization
output_dir | the output directory |
output_file | the output file |
xx | vector of points at which the model is computed, first axis |
yy | vector of points at which the model is computed, second axis |
Definition at line 66 of file Modelling2D.cpp.
|
virtual |
write the model at xx computing 16th, 50th and 84th percentiles from the chains
output_dir | the output directory |
output_file | the output file |
xx | vector of points at which the model is computed |
start | the starting position for each chain |
thin | the position step |
Reimplemented in cbl::modelling::twopt::Modelling_TwoPointCorrelation_wedges.
Definition at line 80 of file Modelling1D.cpp.
|
virtual |
write the model at xx, yy computing 16th, 50th and 84th percentiles from the chains
output_dir | the output directory |
output_file | the output file |
xx | vector of points at which the model is computed, first axis |
yy | vector of points at which the model is computed, second axis |
start | the starting position for each chain |
thin | the position step |
Definition at line 78 of file Modelling2D.cpp.
void cbl::modelling::Modelling::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 |
||
) |
write 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 ( \(n_s\), the number of samples used to estimate the covariance matrix) is 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)} \,, \]
where \(n_b\) is number of data measurements (e.g. the bins of the dataset).
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 director |
root_file | the root of the output files: - file_root_parameters.dat file containing the output of the MCMC sampling for each parameter - file_root_covariance.dat file containing the covariance of the parameters - file_root_chain file containing the chains: the extention can be .dat or .fits |
start | the minimum chain position to be written |
thin | the step used for dilution on screen |
nbins | the number of bins |
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 |
Definition at line 281 of file Modelling.cpp.