Top |
glong | burnin | Read / Write / Construct |
char * | filename | Read / Write / Construct |
int | m2lnp-var | Read / Write / Construct |
NcmMSet * | mset | Read / Write / Construct Only |
GStrv | nadd-val-names | Read / Write / Construct Only |
GStrv | nadd-val-symbols | Read / Write / Construct Only |
guint | nadd-vals | Read / Write / Construct Only |
guint | nchains | Read / Write / Construct Only |
gboolean | read-only | Read / Write / Construct Only |
NcmRNG * | rng | Read / Write |
char * | run-type-string | Read / Write / Construct |
NcmMSetCatalogSync | smode | Read / Write / Construct |
double | sync-interval | Read / Write / Construct |
NcmMSetCatalogTauMethod | tau-method | Read / Write / Construct |
gboolean | weighted | Read / Write / Construct Only |
#define | NCM_TYPE_MSET_CATALOG |
enum | NcmMSetCatalogSync |
enum | NcmMSetCatalogTrimType |
enum | NcmMSetCatalogPostNormMethod |
enum | NcmMSetCatalogTauMethod |
#define | NCM_MSET_CATALOG_EXTNAME |
#define | NCM_MSET_CATALOG_M2LNL_COLNAME |
#define | NCM_MSET_CATALOG_M2LNL_SYMBOL |
#define | NCM_MSET_CATALOG_FIRST_ID_LABEL |
#define | NCM_MSET_CATALOG_M2LNP_ID_LABEL |
#define | NCM_MSET_CATALOG_RNG_ALGO_LABEL |
#define | NCM_MSET_CATALOG_RNG_SEED_LABEL |
#define | NCM_MSET_CATALOG_RNG_STAT_LABEL |
#define | NCM_MSET_CATALOG_RNG_INIS_LABEL |
#define | NCM_MSET_CATALOG_NROWS_LABEL |
#define | NCM_MSET_CATALOG_RTYPE_LABEL |
#define | NCM_MSET_CATALOG_NCHAINS_LABEL |
#define | NCM_MSET_CATALOG_NADDVAL_LABEL |
#define | NCM_MSET_CATALOG_WEIGHTED_LABEL |
#define | NCM_MSET_CATALOG_RTYPE_BSTRAP_MEAN |
#define | NCM_MSET_CATALOG_RTYPE_UNDEFINED |
#define | NCM_MSET_CATALOG_FSYMB_LABEL |
#define | NCM_MSET_CATALOG_ASYMB_LABEL |
#define | NCM_MSET_CATALOG_DIST_EST_SD_SCALE |
NcmMSetCatalog |
GEnum ├── NcmMSetCatalogPostNormMethod ├── NcmMSetCatalogSync ╰── NcmMSetCatalogTauMethod GFlags ╰── NcmMSetCatalogTrimType GObject ╰── NcmMSetCatalog
This class defines a catalog type object. This object can automatically synchronize with a fits file (thought cfitsio).
For Mote Carlo studies, like resampling from a fiducial model or bootstrap, it is used to save the best-fitting values of each realization. Since the order of the resampling is important, due to the fact that we use the same pseudo-random number generator for all resamplings, this object also guarantees the order of the samples added.
For Markov Chain Monte Carlo (MCMC) this object saves the value of the same likelihood in different points of the parameter space.
For both applications this object keeps an interactive mean and variance of the parameters added, this allows a sample by sample analyses of the convergence. Some MCMC convergence diagnostic functions are also implemented here.
NcmMSetCatalog * ncm_mset_catalog_new (NcmMSet *mset
,guint nadd_vals
,guint nchains
,gboolean weighted
,...
);
Creates a new NcmMSetCatalog based on the NcmFit object fit
. The catalog assumes that
the fit
object will remain with the same set of free parameters during its whole lifetime.
If nchains
is larger than one, the catalog will keep track of the statistics of each chain
separately.
mset |
a NcmMSet |
|
nadd_vals |
number of additional values |
|
nchains |
number of different chains in the catalog (>=1) |
|
weighted |
set to TRUE whenever the catalog is weighted |
|
... |
additional values name/symbol pairs |
NcmMSetCatalog * ncm_mset_catalog_new_array (NcmMSet *mset
,guint nadd_vals
,guint nchains
,gboolean weighted
,gchar **names
,gchar **symbols
);
Creates a new NcmMSetCatalog based on the NcmFit object fit
. The catalog assumes that
the fit
object will remain with the same set of free parameters during its whole lifetime.
If nchains
is larger than one, the catalog will keep track of the statistics of each chain
separately.
mset |
a NcmMSet |
|
nadd_vals |
number of additional values |
|
nchains |
number of different chains in the catalog (>=1) |
|
weighted |
set to TRUE whenever the catalog is weighted |
|
names |
additional values name NULL-terminated array. |
[array zero-terminated=1] |
symbols |
additional values symbol NULL-terminated array. |
[array zero-terminated=1] |
NcmMSetCatalog * ncm_mset_catalog_new_from_file (const gchar *filename
,glong burnin
);
Creates a new NcmMSetCatalog from the catalog in the file file
.
It will use also the mset file (same name but with .mset extension).
NcmMSetCatalog * ncm_mset_catalog_new_from_file_ro (const gchar *filename
,glong burnin
);
Creates a new NcmMSetCatalog from the catalog in the file file
.
The file
is opened in a read-only fashion.
It will use also the mset file (same name but with .mset extension).
NcmMSetCatalog *
ncm_mset_catalog_ref (NcmMSetCatalog *mcat
);
Increases the reference count of mcat
atomically.
void
ncm_mset_catalog_free (NcmMSetCatalog *mcat
);
Decreases the reference count of mcat
atomically.
void
ncm_mset_catalog_clear (NcmMSetCatalog **mcat
);
Decrese the reference count of *mcat
atomically and sets the pointer *mcat
to null.
void ncm_mset_catalog_set_file (NcmMSetCatalog *mcat
,const gchar *filename
);
Sets the data filename to be used to sync/save data.
void ncm_mset_catalog_set_sync_mode (NcmMSetCatalog *mcat
,NcmMSetCatalogSync smode
);
Sets the sync mode to smode
.
void ncm_mset_catalog_set_sync_interval (NcmMSetCatalog *mcat
,gdouble interval
);
Sets the minimum time interval between syncs.
void ncm_mset_catalog_set_first_id (NcmMSetCatalog *mcat
,gint first_id
);
Sets the first id of the catalog, mainly used to inform in which realization the catalog starts.
void ncm_mset_catalog_set_run_type (NcmMSetCatalog *mcat
,const gchar *rtype_str
);
Sets the run type string.
void ncm_mset_catalog_set_rng (NcmMSetCatalog *mcat
,NcmRNG *rng
);
Sets the random number generator.
void ncm_mset_catalog_sync (NcmMSetCatalog *mcat
,gboolean check
);
Synchronize memory and data file. If no file was defined, it simply returns.
void ncm_mset_catalog_timed_sync (NcmMSetCatalog *mcat
,gboolean check
);
Synchronize memory and data file if enough time was passed after
the last sync, see ncm_mset_catalog_set_sync_interval()
. If no
file was defined, it simply returns.
void
ncm_mset_catalog_reset_stats (NcmMSetCatalog *mcat
);
Reset catalog statistical quantities.
void
ncm_mset_catalog_reset (NcmMSetCatalog *mcat
);
Clean all catalog data from memory and file. Otherwise it does not change any object's parameter.
void
ncm_mset_catalog_erase_data (NcmMSetCatalog *mcat
);
Erases all data from the fits file associated with the catalog.
void ncm_mset_catalog_set_m2lnp_var (NcmMSetCatalog *mcat
,const gint p
);
Sets p
as the $-2\ln(p)$ parameter index.
gint
ncm_mset_catalog_get_m2lnp_var (NcmMSetCatalog *mcat
);
Gets the $-2\ln(p)$ parameter index.
const gchar *
ncm_mset_catalog_peek_filename (NcmMSetCatalog *mcat
);
Gets the filename associated with mcat
.
NcmRNG *
ncm_mset_catalog_get_rng (NcmMSetCatalog *mcat
);
This function checks if any pseudo random number generator (RNG) is registred in the catalog. If so, it returns it or NULL.
NcmRNG *
ncm_mset_catalog_peek_rng (NcmMSetCatalog *mcat
);
This function checks if any pseudo random number generator (RNG) is registred in the catalog. If so, it returns it or NULL.
gdouble
ncm_mset_catalog_largest_error (NcmMSetCatalog *mcat
);
This function calculates the largest proportional error of the parameters included, i.e., $\text{lre} = \sigma_{\hat{p}}/(|\hat{p}|\sqrt{n})$ where $n$ represents the number of samples in the catalog, $\hat{p}$ is the estimated mean of the parameter $p$ and $\sigma_{\hat{p}}$ its standard deviation.
It tries to guess when $p = 0$. In this case $\sigma_{\hat{p}} \approx |\hat{p}|\sqrt{n}$. Therefore, for $n > 10$, it tests if $\text{lre} \approx 1$ and, if it is the case, it returns $\text{lre} = \sigma_{\hat{p}}/\sqrt{n}$ instead.
guint
ncm_mset_catalog_len (NcmMSetCatalog *mcat
);
Number of itens in the catalog.
guint
ncm_mset_catalog_max_time (NcmMSetCatalog *mcat
);
Number of itens in the catalog divided by the number of chains.
guint
ncm_mset_catalog_nchains (NcmMSetCatalog *mcat
);
Number of chains in the catalog.
guint
ncm_mset_catalog_nadd_vals (NcmMSetCatalog *mcat
);
Number of additional variables in the catalog.
gboolean
ncm_mset_catalog_weighted (NcmMSetCatalog *mcat
);
Whether the catalog has weights.
guint ncm_mset_catalog_get_row_from_time (NcmMSetCatalog *mcat
,gint t
);
const gchar * ncm_mset_catalog_col_name (NcmMSetCatalog *mcat
,guint i
);
const gchar * ncm_mset_catalog_col_full_name (NcmMSetCatalog *mcat
,guint i
);
const gchar * ncm_mset_catalog_col_symb (NcmMSetCatalog *mcat
,guint i
);
gboolean ncm_mset_catalog_col_by_name (NcmMSetCatalog *mcat
,const gchar *name
,guint *col_index
);
Finds the column name
in the catalog mcat
.
void ncm_mset_catalog_set_burnin (NcmMSetCatalog *mcat
,glong burnin
);
Sets the number of elements to ignore when reading from a catalogue, it must be set before loading data from a file.
It will not affect a catalogue in any other context, only when reading data from a file. It is recommended to be used only when analysing a catalogue.
glong
ncm_mset_catalog_get_burnin (NcmMSetCatalog *mcat
);
Gets the burn-in size, see ncm_mset_catalog_set_burnin()
.
void ncm_mset_catalog_set_tau_method (NcmMSetCatalog *mcat
,NcmMSetCatalogTauMethod tau_method
);
Sets the autocorrelation time method to tau_method
.
NcmMSetCatalogTauMethod
ncm_mset_catalog_get_tau_method (NcmMSetCatalog *mcat
);
void ncm_mset_catalog_add_from_mset (NcmMSetCatalog *mcat
,NcmMSet *mset
,...
);
This function adds a new element to the catalog using the parameters from mset
.
It assumes that mset
is compatible with the catalog and expect the
right number of additional values.
void ncm_mset_catalog_add_from_mset_array (NcmMSetCatalog *mcat
,NcmMSet *mset
,gdouble *ax
);
This funtion adds a new element to the catalog using the parameters from mset
.
It assumes that mset
is compatible with the catalog and expect the
right number of additional values in the array ax
.
void ncm_mset_catalog_add_from_vector (NcmMSetCatalog *mcat
,NcmVector *vals
);
Adds a new element to the catalog using the values from the vector
vals
.
void ncm_mset_catalog_add_from_vector_array (NcmMSetCatalog *mcat
,NcmVector *vals
,gdouble *ax
);
Adds a new element to the catalog using the parameter values from the
vector vals
and additional parameters from array ax
.
void
ncm_mset_catalog_log_current_stats (NcmMSetCatalog *mcat
);
Logs the current means and standard deviations of the catalog's parameters.
void
ncm_mset_catalog_log_current_chain_stats
(NcmMSetCatalog *mcat
);
Logs the current means and standard deviations of the catalog's parameters for each chain.
NcmMSet *
ncm_mset_catalog_get_mset (NcmMSetCatalog *mcat
);
Gets the NcmMSet catalog from mcat
.
NcmMSet *
ncm_mset_catalog_peek_mset (NcmMSetCatalog *mcat
);
Gets the NcmMSet catalog from mcat
.
const gchar *
ncm_mset_catalog_get_run_type (NcmMSetCatalog *mcat
);
Gets the run type string from mcat
.
NcmStatsVec *
ncm_mset_catalog_peek_pstats (NcmMSetCatalog *mcat
);
Peeks the parameters NcmStatsVec object.
NcmStatsVec *
ncm_mset_catalog_peek_e_mean_stats (NcmMSetCatalog *mcat
);
Peeks the ensemble NcmStatsVec object.
NcmStatsVec * ncm_mset_catalog_peek_chain_pstats (NcmMSetCatalog *mcat
,const guint i
);
Peeks the chain NcmStatsVec object.
GArray *
ncm_mset_catalog_peek_accept_ratio_array
(NcmMSetCatalog *mcat
);
Peeks the acceptance ratio array.
NcmVector * ncm_mset_catalog_peek_row (NcmMSetCatalog *mcat
,guint i
);
Gets the i
-th row.
NcmVector *
ncm_mset_catalog_peek_current_row (NcmMSetCatalog *mcat
);
Gets the last added row.
NcmVector *
ncm_mset_catalog_peek_current_e_mean (NcmMSetCatalog *mcat
);
Gets the last ensemble mean.
NcmVector *
ncm_mset_catalog_peek_current_e_var (NcmMSetCatalog *mcat
);
Gets the last ensemble variance.
NcmVector * ncm_mset_catalog_peek_e_mean_t (NcmMSetCatalog *mcat
,guint t
);
Gets the mean of the t
-th ensemble.
NcmVector * ncm_mset_catalog_peek_e_var_t (NcmMSetCatalog *mcat
,guint t
);
Gets the variance of the t
-th ensemble.
gdouble ncm_mset_catalog_get_post_lnnorm (NcmMSetCatalog *mcat
,gdouble *post_lnnorm_sd
);
Computes, if necessary, the posterior normalization.
gdouble ncm_mset_catalog_get_post_lnvol (NcmMSetCatalog *mcat
,const gdouble level
,gdouble *glnvol
);
Computes the volume of the level
region of the posterior.
Sets into glnvol
the log volume of the Gaussian approximation
of the posterior.
gdouble ncm_mset_catalog_get_nth_m2lnL_percentile (NcmMSetCatalog *mcat
,const gdouble p
,guint *nth
);
Computes the p
percentile of the likelihood.
gdouble
ncm_mset_catalog_get_bestfit_m2lnL (NcmMSetCatalog *mcat
);
NcmVector *
ncm_mset_catalog_get_bestfit_row (NcmMSetCatalog *mcat
);
void ncm_mset_catalog_get_mean (NcmMSetCatalog *mcat
,NcmVector **mean
);
Gets the current mean vector.
void ncm_mset_catalog_get_covar (NcmMSetCatalog *mcat
,NcmMatrix **cov
);
Gets the current covariance matrix.
void ncm_mset_catalog_get_full_covar (NcmMSetCatalog *mcat
,NcmMatrix **cov
);
Gets the current full (including additional values) covariance matrix.
void
ncm_mset_catalog_log_full_covar (NcmMSetCatalog *mcat
);
Logs the current full (including additional values) covariance matrix.
void ncm_mset_catalog_estimate_autocorrelation_tau (NcmMSetCatalog *mcat
,gboolean force_single_chain
);
Updates the internal estimates of the integrate autocorrelation time.
NcmVector *
ncm_mset_catalog_peek_autocorrelation_tau
(NcmMSetCatalog *mcat
);
Gets the last estimate of the autocorrelation tau calculated
in the last call of ncm_mset_catalog_estimate_autocorrelation_tau()
.
gdouble ncm_mset_catalog_get_param_shrink_factor (NcmMSetCatalog *mcat
,guint p
);
Gets the current shrink factor of parameter p
.
gdouble
ncm_mset_catalog_get_shrink_factor (NcmMSetCatalog *mcat
);
Gets the current shrink factor which is the multivatiate potential scale reduction factor (MPSRF), namely, $$\hat{R}^p = \sqrt{\frac{n - 1}{n} + \left( \frac{m + 1}{m} \right) \lambda_1},$$ where $n$ is the number of points of one chain, $m$ is the number of chains and $\lambda_1$ is the largest eigenvalue of the positive definite matrix $W^{-1}B/n$.
$W$ is the within-chain covariance: $$W = $$ arithmetical mean of the covariance matrices of each chain.
$B$ is the between-chain covariance: $$B = $$ covariance between the means of each chain.
Refined version: $$\hat{R}^p = \sqrt{\frac{\hat{d} + 3}{\hat{d} + 1} \left(\frac{n - 1}{n} + \left( \frac{m + 1}{m} \right) \lambda_1\right)},$$ where $\hat{d} = 2 \hat{V}^2 / \widehat{Var}(\hat{V})$, $$\hat{V} = \frac{n -1}{n}W + \frac{m + 1}{m} \frac{B}{n}.$$
Some references for this MCMC convergence diagnostic: Brooks and Gelman (1998), Gelman and Rubin (1992), SAS/STAT.
void ncm_mset_catalog_param_pdf (NcmMSetCatalog *mcat
,guint i
);
Bins and calculates the pdf associated with the parameter i
.
gdouble ncm_mset_catalog_param_pdf_pvalue (NcmMSetCatalog *mcat
,gdouble pvalue
,gboolean both
);
Calculates the p-value associated with the parameter value pvalue
.
NcmMatrix * ncm_mset_catalog_calc_ci_direct (NcmMSetCatalog *mcat
,NcmMSetFunc *func
,NcmVector *x_v
,GArray *p_val
);
Calculates the mean and the confidence interval (CI) for the value of func
for each p-value in p_val
. It stores the results in a NcmVector, where the
first element contains the mean and the following contain the lower and
upper bounds for each p-value in p_val
.
This function calculates the quantiles directly using:
gsl_stats_quantile_from_sorted_data for this reason it must allocates the
catalog size times the number of elements in x
, for a less memory intensive
version use ncm_mset_catalog_calc_ci_interp()
.
The NcmMSetFunc func
must be of dimension one.
If p_val
contains two values ($1\sigma$) 0.6827 and ($\sigma$) 0.9545,
the first element will contain the mean, the second and third, the lower
and upper bounds, respectively. Then, the fourth and fifth elements the
lower and upper bounds of $2\sigma$ CI.
mcat |
||
func |
a NcmMSetFunc of type n-n |
|
x_v |
NcmVector of arguments of |
|
p_val |
p-values for the confidence intervals. |
[element-type double] |
a NcmVector containing the mean and lower/upper bound of the confidence interval for func
.
[transfer full]
NcmMatrix * ncm_mset_catalog_calc_ci_interp (NcmMSetCatalog *mcat
,NcmMSetFunc *func
,NcmVector *x_v
,GArray *p_val
,guint nodes
,NcmFitRunMsgs mtype
);
Calculates the mean and the confidence interval (CI) for the value of func
for each p-value in p_val
. It stores the results in a NcmMatrix, where the
first element contains the mean and the following contain the lower and
upper bounds for each p-value in p_val
.
This function creates an approximation of the distribution for each value of
the function func
and calculates the quantiles from this approximation.
The NcmMSetFunc func
must be of dimension one.
If p_val
contains two values ($1\sigma$) 0.6827 and ($\sigma$) 0.9545,
the first element will contain the mean, the second and third, the lower
and upper bounds, respectively. Then, the fourth and fifth elements the
lower and upper bounds of $2\sigma$ CI.
mcat |
||
func |
a NcmMSetFunc of type n-n |
|
x_v |
NcmVector of arguments of |
|
p_val |
p-values for the confidence intervals. |
[element-type double] |
nodes |
number of nodes in the distribution approximations |
|
mtype |
NcmFitRunMsgs log level |
a NcmMatrix containing the mean and lower/upper bound of the confidence interval for func
.
[transfer full]
NcmMatrix * ncm_mset_catalog_calc_pvalue (NcmMSetCatalog *mcat
,NcmMSetFunc *func
,NcmVector *x_v
,GArray *lim
,guint nodes
,NcmFitRunMsgs mtype
);
Calculates the p-values for the value of func
for each limit in lim
, integrating the probability distribution function from
the left tail to lim
. It stores the results in a NcmMatrix, where the
first element contains the p-value with respect to the first lim
, and so on.
The NcmMSetFunc func
must be of dimension one.
mcat |
||
func |
a NcmMSetFunc of type n-n |
|
x_v |
NcmVector of arguments of |
|
lim |
integration limits to compute the p-value. |
[element-type double] |
nodes |
number of nodes in the distribution approximations |
|
mtype |
NcmFitRunMsgs log level |
NcmStatsDist1d * ncm_mset_catalog_calc_distrib (NcmMSetCatalog *mcat
,NcmMSetFunc *func
,NcmFitRunMsgs mtype
);
Calculates the distribution of func
.
This function creates an approximation of the distribution for each value of
the function func
calculated in each model in mcat
.
NcmStatsDist1d * ncm_mset_catalog_calc_param_distrib (NcmMSetCatalog *mcat
,const NcmMSetPIndex *pi
,NcmFitRunMsgs mtype
);
Calculates the distribution of parameter pi
.
This function creates an approximation of the distribution for each value of
the parameter pi
in mcat
.
NcmStatsDist1d * ncm_mset_catalog_calc_add_param_distrib (NcmMSetCatalog *mcat
,guint add_param
,NcmFitRunMsgs mtype
);
Calculates the distribution of parameter pi
.
This function creates an approximation of the distribution for each value of
the parameter pi
in mcat
.
void ncm_mset_catalog_calc_param_ensemble_evol (NcmMSetCatalog *mcat
,const NcmMSetPIndex *pi
,guint nsteps
,NcmFitRunMsgs mtype
,NcmVector **pval
,NcmMatrix **t_evol
);
Calculates the time evolution of the parameter pi
distribution.
mcat |
||
pi |
||
nsteps |
number of steps to calculate the distribution |
|
mtype |
NcmFitRunMsgs log level |
|
pval |
output NcmVector containing parameter values. |
[out callee-allocates] |
t_evol |
output NcmMatrix containing probability distribution evolution. |
[out callee-allocates] |
void ncm_mset_catalog_calc_add_param_ensemble_evol (NcmMSetCatalog *mcat
,guint add_param
,guint nsteps
,NcmFitRunMsgs mtype
,NcmVector **pval
,NcmMatrix **t_evol
);
Calculates the time evolution of the parameter pi
distribution.
mcat |
||
add_param |
additional parameter index |
|
nsteps |
number of steps to calculate the distribution |
|
mtype |
NcmFitRunMsgs log level |
|
pval |
output NcmVector containing parameter values. |
[out callee-allocates] |
t_evol |
output NcmMatrix containing probability distribution evolution. |
[out callee-allocates] |
void ncm_mset_catalog_trim (NcmMSetCatalog *mcat
,const guint tc
,const guint thin
);
Drops all points in the catalog such that $t < t_c$ and skips
every thin
-1 rows, creating a thinner catalog.
This function trims the first $t_c \times n_\mathrm{chains}$
points from the catalog. Creates a backup of the original file.
void ncm_mset_catalog_trim_p (NcmMSetCatalog *mcat
,const gdouble p
);
Drops all points in the catalog such that the first p
percent of the catalog is dropped.
guint ncm_mset_catalog_trim_oob (NcmMSetCatalog *mcat
,const gchar *out_file
);
Remove all points that are outside the bounds defined by
the catalog mset file. The catalog will always have a
single chain after the trimming. The result is saved to out_file
.
void
ncm_mset_catalog_remove_last_ensemble (NcmMSetCatalog *mcat
);
Removes the last ensemble point in the catalog, i.e., removes the last 'number of chains' points of the catalog. Creates a backup of the original file.
guint ncm_mset_catalog_calc_max_ess_time (NcmMSetCatalog *mcat
,const guint ntests
,gdouble *max_ess
,NcmFitRunMsgs mtype
);
Calculates the time $t_m$ that maximizes the ESS for all
elements of the catalog. If the number of chains in the catalog is larger
than one, it considers the whole catalog otherwise it considers the ensemble
means. The variable ntests
control the number of divisions where the ESS
will be calculated, if it is zero the default 10 tests will be used.
mcat |
||
ntests |
number of tests |
|
max_ess |
the maximum effective sample size (ESS). |
[out] |
mtype |
NcmFitRunMsgs log level |
guint ncm_mset_catalog_calc_heidel_diag (NcmMSetCatalog *mcat
,const guint ntests
,const gdouble pvalue
,NcmFitRunMsgs mtype
);
Applies the Heidelberger and Welch's convergence diagnostic to the catalog,
see ncm_stats_vec_heidel_diag()
for mode details. If the number of chains in
the catalog is larger than one, it considers the whole catalog otherwise it
considers the ensemble means. The variable ntests
control the number of
divisions where the test will be applied, if it is zero the default 10 tests
will be used.
mcat |
||
ntests |
number of tests |
|
pvalue |
the required Schruben test p-value |
|
mtype |
NcmFitRunMsgs log level |
guint ncm_mset_catalog_calc_const_break (NcmMSetCatalog *mcat
,guint p
,NcmFitRunMsgs mtype
);
Fits the model: $$f(t) = c_0 + \theta_r(t-t_0)\left[c_1(t-t_0) + c_2\frac{(t-t_0)^2}{2}\right].$$ to estimate the time $t_0$ where the chain stops evolving.
void ncm_mset_catalog_trim_by_type (NcmMSetCatalog *mcat
,guint ntests
,NcmMSetCatalogTrimType trim_type
,NcmFitRunMsgs mtype
);
Calculates the time $t_m$ that satisfies all trimming options
in trim_type
. Then drops all elements of the catalog and drops
all points $t < t_m$.
mcat |
||
ntests |
number of tests |
|
trim_type |
the trimming type to apply NcmMSetCatalogTrimType |
|
mtype |
NcmFitRunMsgs log level |
guint ncm_mset_catalog_max_ess_time_by_chain (NcmMSetCatalog *mcat
,const guint ntests
,gdouble *max_ess
,NcmFitRunMsgs mtype
);
Calculates the time $t_m$ that maximizes the ESS for each chain of the catalog.
The variable ntests
control the number of divisions where the ESS
will be calculated, if it is zero the default 10 tests will be used.
mcat |
||
ntests |
number of tests |
|
max_ess |
the maximum effective sample size (ESS). |
[out] |
mtype |
NcmFitRunMsgs log level |
guint ncm_mset_catalog_heidel_diag_by_chain (NcmMSetCatalog *mcat
,const guint ntests
,const gdouble pvalue
,gdouble *wp_pvalue
,NcmFitRunMsgs mtype
);
Calculates the lowest time $t_m$ where all chains satisfy the Heidelberger
and Welch's convergence diagnostic. The variable ntests
control the number
of divisions where the test will be calculated, if it is zero the default
10 tests will be used.
mcat |
||
ntests |
number of tests |
|
pvalue |
required p-value |
|
wp_pvalue |
worst parameter p-value. |
[out] |
mtype |
NcmFitRunMsgs log level |
Catalog sync modes.
Catalog will be synchronized only when closing the file or with an explicit call of |
||
Catalog will be synchronized in every catalog addition. |
||
Catalog will be synchronized with a minimum time interval between syncs. |
See ncm_mset_catalog_calc_max_ess_time()
and ncm_mset_catalog_calc_heidel_diag()
.
See ncm_mset_catalog_calc_max_ess_time()
and ncm_mset_catalog_calc_heidel_diag()
.
“burnin”
property “burnin” glong
Burn-in size.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct
Allowed values: >= 0
Default value: 0
“filename”
property “filename” char *
Catalog filename.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct
Default value: NULL
“m2lnp-var”
property “m2lnp-var” int
Index of the variable representing m2lnp.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct
Allowed values: >= -1
Default value: -1
“mset”
property“mset” NcmMSet *
NcmMSet object.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct Only
“nadd-val-names”
property “nadd-val-names” GStrv
Additional value names.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct Only
“nadd-val-symbols”
property “nadd-val-symbols” GStrv
Additional value symbols.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct Only
“nadd-vals”
property “nadd-vals” guint
Number of additional values.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct Only
Default value: 0
“nchains”
property “nchains” guint
Number of different chains in the catalog.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct Only
Allowed values: >= 1
Default value: 1
“read-only”
property “read-only” gboolean
If the fits catalogue must be open in the readonly mode.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct Only
Default value: FALSE
“rng”
property“rng” NcmRNG *
Random number generator object.
Owner: NcmMSetCatalog
Flags: Read / Write
“run-type-string”
property “run-type-string” char *
Run type string.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct
Default value: "undefined-run"
“smode”
property“smode” NcmMSetCatalogSync
Catalog sync mode.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct
Default value: NCM_MSET_CATALOG_SYNC_AUTO
“sync-interval”
property “sync-interval” double
Data sync interval.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct
Allowed values: [0,1000]
Default value: 10
“tau-method”
property“tau-method” NcmMSetCatalogTauMethod
Method used to calculate the autocorrelation time.
Owner: NcmMSetCatalog
Flags: Read / Write / Construct
Default value: NCM_MSET_CATALOG_TAU_METHOD_AR_MODEL