Top |
NcmObjArray * | equality-constraints | Read / Write |
NcmVector * | equality-constraints-tot | Read / Write |
NcmFitGradType | grad-type | Read / Write / Construct |
NcmObjArray * | inequality-constraints | Read / Write |
NcmVector * | inequality-constraints-tot | Read / Write |
NcmLikelihood * | likelihood | Read / Write / Construct Only |
double | m2lnL-abstol | Read / Write / Construct |
double | m2lnL-reltol | Read / Write / Construct |
guint | maxiter | Read / Write / Construct |
NcmMSet * | mset | Read / Write / Construct Only |
double | params-reltol | Read / Write / Construct |
NcmFitState * | state | Read / Write / Construct Only |
NcmFit * | sub-fit | Read / Write |
#define | NCM_TYPE_FIT |
enum | NcmFitType |
enum | NcmFitGradType |
enum | NcmFitRunMsgs |
struct | NcmFitClass |
#define | NCM_FIT_DEFAULT_M2LNL_RELTOL |
#define | NCM_FIT_DEFAULT_M2LNL_ABSTOL |
#define | NCM_FIT_DEFAULT_PARAMS_RELTOL |
#define | NCM_FIT_DEFAULT_MAXITER |
NcmFit |
GEnum ├── NcmFitGradType ├── NcmFitRunMsgs ╰── NcmFitType GObject ╰── NcmFit ├── NcmFitGSLLS ├── NcmFitGSLMM ├── NcmFitGSLMMS ├── NcmFitLevmar ╰── NcmFitNLOpt
NcmFit * ncm_fit_factory (NcmFitType ftype
,gchar *algo_name
,NcmLikelihood *lh
,NcmMSet *mset
,NcmFitGradType gtype
);
Creates a new NcmFit object.
NcmFit * ncm_fit_copy_new (NcmFit *fit
,NcmLikelihood *lh
,NcmMSet *mset
,NcmFitGradType gtype
);
Duplicates the NcmFit object with new references for its contents.
NcmFit * ncm_fit_dup (NcmFit *fit
,NcmSerialize *ser
);
Duplicates the NcmFit object duplicating all its contents.
void
ncm_fit_free (NcmFit *fit
);
Atomically decrements the reference count of fit
by one. If the reference count drops to 0,
all memory allocated by fit
is released.
void
ncm_fit_clear (NcmFit **fit
);
The reference count of fit
is decreased and the pointer is set to NULL.
void ncm_fit_set_sub_fit (NcmFit *fit
,NcmFit *sub_fit
);
Sets a NcmFit object to be used as subsidiary fit.
NcmFit *
ncm_fit_get_sub_fit (NcmFit *fit
);
Gets a NcmFit object to be used as subsidiary fit.
gboolean
ncm_fit_has_sub_fit (NcmFit *fit
);
Checks if fit
has a subsidiary fit.
void ncm_fit_set_grad_type (NcmFit *fit
,NcmFitGradType gtype
);
Sets the differentiation method to be used.
void ncm_fit_set_maxiter (NcmFit *fit
,guint maxiter
);
Sets the maximum number of iterations.
void ncm_fit_set_m2lnL_reltol (NcmFit *fit
,gdouble tol
);
Sets the relative tolerance for the m2lnL.
void ncm_fit_set_m2lnL_abstol (NcmFit *fit
,gdouble tol
);
Sets the absolute tolerance for the m2lnL.
void ncm_fit_set_params_reltol (NcmFit *fit
,gdouble tol
);
Sets the relative tolerance for the fitted parameters.
void ncm_fit_set_messages (NcmFit *fit
,NcmFitRunMsgs mtype
);
Sets the log level for the messages to be printed during the fit.
NcmFitGradType
ncm_fit_get_grad_type (NcmFit *fit
);
Gets the differentiation method to be used.
guint
ncm_fit_get_maxiter (NcmFit *fit
);
Gets the maximum number of iterations.
gdouble
ncm_fit_get_m2lnL_reltol (NcmFit *fit
);
Gets the relative tolerance for the m2lnL.
gdouble
ncm_fit_get_m2lnL_abstol (NcmFit *fit
);
Gets the absolute tolerance for the m2lnL.
gdouble
ncm_fit_get_params_reltol (NcmFit *fit
);
Gets the relative tolerance for the fitted parameters.
NcmFitRunMsgs
ncm_fit_get_messages (NcmFit *fit
);
Gets the log level for the messages to be printed during the fit.
gboolean
ncm_fit_is_least_squares (NcmFit *fit
);
Indicates if the least squares fitting is being used (TRUE) or not (FALSE).
NcmLikelihood *
ncm_fit_peek_likelihood (NcmFit *fit
);
Peeks the NcmLikelihood object.
void ncm_fit_params_set (NcmFit *fit
,guint i
,const gdouble x
);
Sets the parameters vector.
void ncm_fit_params_set_vector (NcmFit *fit
,NcmVector *x
);
Sets the parameters vector.
void ncm_fit_params_set_vector_offset (NcmFit *fit
,NcmVector *x
,guint offset
);
Sets the parameters from vector x
starting at offset
.
void ncm_fit_params_set_array (NcmFit *fit
,const gdouble *x
);
Sets the parameters from array x
.
void ncm_fit_params_set_gsl_vector (NcmFit *fit
,const gsl_vector *x
);
Sets the parameters from a gsl_vector x
.
[skip]
void ncm_fit_add_equality_constraint (NcmFit *fit
,NcmMSetFunc *func
,const gdouble tot
);
Adds an equality constraint with the function func
and the tolerance tot
.
void ncm_fit_add_inequality_constraint (NcmFit *fit
,NcmMSetFunc *func
,const gdouble tot
);
Adds an inequality constraint with the function func
and the tolerance tot
.
void
ncm_fit_remove_equality_constraints (NcmFit *fit
);
Removes all equality constraints.
void
ncm_fit_remove_inequality_constraints (NcmFit *fit
);
Removes all inequality constraints.
guint
ncm_fit_equality_constraints_len (NcmFit *fit
);
Gets the number of equality constraints.
guint
ncm_fit_inequality_constraints_len (NcmFit *fit
);
Gets the number of inequality constraints.
void ncm_fit_get_equality_constraint (NcmFit *fit
,guint i
,NcmMSetFunc **func
,gdouble *tot
);
Gets the equality constraint at index i
.
fit |
a NcmFit |
|
i |
index of the equality constraint |
|
func |
a NcmMSetFunc. |
[out][transfer none] |
tot |
a double. |
[out] |
void ncm_fit_get_inequality_constraint (NcmFit *fit
,guint i
,NcmMSetFunc **func
,gdouble *tot
);
Gets the inequality constraint at index i
.
fit |
a NcmFit |
|
i |
index of the inequality constraint |
|
func |
a NcmMSetFunc. |
[out][transfer none] |
tot |
a double. |
[out] |
void ncm_fit_set_logger (NcmFit *fit
,NcmFitWriter writer
,NcmFitUpdater updater
,NcmFitUpdateChange start_update
,NcmFitUpdateChange end_update
);
Sets the logger functions. The writer
function is called to write the
messages to the log. The updater
function is called to update the
parameters. The start_update
is called before the minimization starts and
the end_update
is called after the minimization ends.
fit |
a NcmFit |
|
writer |
a NcmFitWriter. |
[scope notified] |
updater |
[scope notified] | |
start_update |
[scope notified][nullable] | |
end_update |
[scope notified][nullable] |
void
ncm_fit_log_info (NcmFit *fit
);
Prints to the log the data set and the model set.
void
ncm_fit_log_covar (NcmFit *fit
);
Prints to the log file the names and indices of the fitted parameters, their best-fit values, standard deviations and correlation matrix.
void
ncm_fit_log_start (NcmFit *fit
);
This function prints in the log the initial state.
void
ncm_fit_log_state (NcmFit *fit
);
This function prints in the log the current state.
void
ncm_fit_log_step (NcmFit *fit
);
This function prints in the log one step of the minimization.
void ncm_fit_log_step_error (NcmFit *fit
,const gchar *strerror
,...
);
This function prints in the log the error message.
void
ncm_fit_log_end (NcmFit *fit
);
This function prints in the log the precision with which the best-fit was found.
void ncm_fit_data_m2lnL_val (NcmFit *fit
,gdouble *data_m2lnL
);
This function computes minus two times the logarithm base e of the likelihood
using only the data set and not considering any prior. The result is set
on data_m2lnL
.
void ncm_fit_priors_m2lnL_val (NcmFit *fit
,gdouble *priors_m2lnL
);
This function computes minus two times the logarithm base e of the likelihood
using the data set and taking into account the assumed priors. The result is
set on priors_m2lnL
.
void ncm_fit_m2lnL_val (NcmFit *fit
,gdouble *m2lnL
);
Computes minus two times the logarithm base e of the likelihood.
void ncm_fit_m2lnL_grad (NcmFit *fit
,NcmVector *df
);
Computes the gradient of the minus two times the logarithm base e of the likelihood.
void ncm_fit_m2lnL_val_grad (NcmFit *fit
,gdouble *result
,NcmVector *df
);
Computes the minus two times the logarithm base e of the likelihood and its gradient.
void ncm_fit_ls_J (NcmFit *fit
,NcmMatrix *J
);
Computes the Jacobian matrix for the least squares problem.
void ncm_fit_ls_f_J (NcmFit *fit
,NcmVector *f
,NcmMatrix *J
);
Computes the residuals vector and the Jacobian matrix for the least squares problem.
void
ncm_fit_obs_fisher (NcmFit *fit
);
Computes the covariance matrix using the inverse of the Hessian matrix
$\partial_i\partial_j -\ln(L)$, where the derivatives are taken with respect to the
free parameters. This function does not use the gradient defined in the fit
object,
it always uses the accurate numerical differentiation methods implemented in the
NcmDiff object.
It sets both the covariance matrix and the Hessian matrix in the NcmFitState object
associated to the fit
object.
void
ncm_fit_ls_fisher (NcmFit *fit
);
Computes the covariance matrix using the jacobian matrix and the least squares
problem, see ncm_fit_ls_J()
. Note that this function uses the gradient defined
in the fit
object using ncm_fit_set_grad_type()
to compute the jacobian matrix.
This function is similar to ncm_fit_fisher()
but it computes the covariance
using the full jacobian matrix instead adding the individual contributions
of each NcmData in the data set.
It sets both the covariance matrix in the NcmFitState object associated to the
fit
object.
void
ncm_fit_fisher (NcmFit *fit
);
Calculates the covariance from the Fisher matrix, see ncm_dataset_fisher_matrix()
.
Note that this function does not use the gradient defined in the fit
object, it
always uses the accurate numerical differentiation methods implemented in the
NcmDiff object.
It sets both the covariance matrix in the NcmFitState object associated to the
fit
object.
NcmVector * ncm_fit_fisher_bias (NcmFit *fit
,NcmVector *f_true
);
Calculates the covariance from the Fisher matrix and the bias vector, see
ncm_dataset_fisher_matrix_bias()
. The bias vector is calculated using the
the theory vector f_true
as the true model expectation values.
Note that this function does not use the gradient defined in the fit
object, it
always uses the accurate numerical differentiation methods implemented in the
NcmDiff object.
It sets the covariance matrix in the NcmFitState object associated to the
fit
object. Moreover, it computes the final bias vector, that is, the inverse
of the Fisher matrix times the bias vector returns it.
void
ncm_fit_numdiff_m2lnL_covar (NcmFit *fit
);
ncm_fit_numdiff_m2lnL_covar
has been deprecated since version 0.18.2 and should not be used in newly-written code.
Use ncm_fit_obs_fisher()
instead.
Calcualtes the covariance matrix using the inverse of the Hessian matrix $\partial_i\partial_j -\ln(L)$, where the derivatives are taken with respect to the free parameters.
gdouble
ncm_fit_numdiff_m2lnL_lndet_covar (NcmFit *fit
);
Calculates the logarithm of the determinant of the covariance matrix using the inverse of the Hessian matrix $\partial_i\partial_j -\ln(L)$, where the derivatives are taken with respect to the free parameters.
It sets the covariance matrix in the NcmFitState object associated to the fit
object.
NcmMatrix *
ncm_fit_get_covar (NcmFit *fit
);
Returns a copy of the covariance matrix (pre-calculated by, e.g, ncm_fit_numdiff_m2lnL_covar()
).
gdouble ncm_fit_covar_var (NcmFit *fit
,NcmModelID mid
,guint pid
);
Computes the variance of the fitted parameter pid
of the model mid
.
gdouble ncm_fit_covar_sd (NcmFit *fit
,NcmModelID mid
,guint pid
);
Computes the standard deviation of the fitted parameter pid
of the model mid
.
gdouble ncm_fit_covar_cov (NcmFit *fit
,NcmModelID mid1
,guint pid1
,NcmModelID mid2
,guint pid2
);
Computes the covariance between the parameters pid1
and pid2
of the models
mid1
and mid2
, respectively.
fit |
a NcmFit |
|
mid1 |
||
pid1 |
the parameter's index of the model |
|
mid2 |
||
pid2 |
the parameter's index of the model |
gdouble ncm_fit_covar_cor (NcmFit *fit
,NcmModelID mid1
,guint pid1
,NcmModelID mid2
,guint pid2
);
Computes the correlation between the parameters pid1
and pid2
of the models
mid1
and mid2
, respectively.
fit |
a NcmFit |
|
mid1 |
||
pid1 |
the parameter's index of the model |
|
mid2 |
||
pid2 |
the parameter's index of the model |
gdouble ncm_fit_covar_fparam_var (NcmFit *fit
,guint fpi
);
Computes the variance of the fitted parameter fpi
.
This index refers to the list of all FREE parameters set in the MSet.
See also the similar function ncm_fit_covar_var()
to which one has to provide
the respective model of the parameter.
gdouble ncm_fit_covar_fparam_sd (NcmFit *fit
,guint fpi
);
Computes the standard deviation of the fitted parameter fpi
.
This index refers to the list of all FREE parameters set in the MSet.
See also the similar function ncm_fit_covar_sd()
to which one has to provide
the respective model of the parameter.
gdouble ncm_fit_covar_fparam_cov (NcmFit *fit
,guint fpi1
,guint fpi2
);
Computes the covariance between the fitted parameters fpi1
and fpi2
.
These indices refers to the list of all FREE parameters set in the MSet.
See also the similar function ncm_fit_covar_cov()
to which one has to provide
the respective models of the parameters.
gdouble ncm_fit_covar_fparam_cor (NcmFit *fit
,guint fpi1
,guint fpi2
);
Computes the correlation between the fitted parameters fpi1
and fpi2
.
These indices refers to the list of all FREE parameters set in the MSet.
See also the similar function ncm_fit_covar_cor()
to which one has to provide
the respective models of the parameters.
gboolean ncm_fit_run_restart (NcmFit *fit
,NcmFitRunMsgs mtype
,const gdouble abstol
,const gdouble reltol
,NcmMSet *save_mset
,const gchar *mset_file
);
Re-runs the fit until the difference between fits are less than the required tolerance, i.e., $$ m2lnL_{i-1} - m2lnL_i < \mathrm{abstol} + \mathrm{reltol}\vert m2lnL_{i-1}\vert. $$
NcmMatrix * ncm_fit_lr_test_range (NcmFit *fit
,NcmModelID mid
,guint pid
,gdouble start
,gdouble stop
,guint nsteps
);
Computes the likelihood ratio test for the parameter pid
of the model mid
in the interval [start
, stop
] subdivided by nsteps
. The function returns a
NcmMatrix with the following columns:
Parameter value.
The difference in -2 times the natural logarithm of the likelihood
between the full model and the model with the parameter pid
fixed to the
value in the first column.
The probability density of the difference in -2 times the natural logarithm
of the likelihood between the full model and the model with the parameter
pid
fixed to the value in the first column, assuming a chi-squared distribution
with one degree of freedom.
Cumulative probability (two-sides) of the difference (Column 3).
fit |
a NcmFit |
|
mid |
a NcmModelID. |
|
pid |
the parameter id |
|
start |
starting value |
|
stop |
ending value |
|
nsteps |
step size |
gdouble ncm_fit_lr_test (NcmFit *fit
,NcmModelID mid
,guint pid
,gdouble val
,gint dof
);
Computes the likelihood ratio test for the parameter pid
of the model mid
with the value val
. The function returns the probability of the null hypothesis
assuming a chi-squared distribution with dof
degrees of freedom. That is,
it computes the left tail of the chi-squared distribution with dof
degrees of
freedom.
fit |
a NcmFit. |
|
mid |
a NcmModelID. |
|
pid |
the parameter id |
|
val |
the parameter value |
|
dof |
degrees of freedom |
void ncm_fit_function_error (NcmFit *fit
,NcmMSetFunc *func
,gdouble *x
,gdouble *f
,gdouble *sigma_f
);
Propagates the errors in the free parameters to the function func
.
If no free parameters are set in the NcmMSet object associated to the
fit
object, the computed error is 0.
Note that the covariance matrix must be calculated before calling this function. Otherwise, it will raise an error.
fit |
a NcmFit |
|
func |
||
x |
the argument for the function |
|
f |
the function value. |
[out] |
sigma_f |
the error in the function value. |
[out] |
“equality-constraints”
property“equality-constraints” NcmObjArray *
Equality constraints array.
Owner: NcmFit
Flags: Read / Write
“equality-constraints-tot”
property“equality-constraints-tot” NcmVector *
Equality constraints tolerance.
Owner: NcmFit
Flags: Read / Write
“grad-type”
property“grad-type” NcmFitGradType
Differentiation method.
Owner: NcmFit
Flags: Read / Write / Construct
Default value: NCM_FIT_GRAD_NUMDIFF_FORWARD
“inequality-constraints”
property“inequality-constraints” NcmObjArray *
Inequality constraints array.
Owner: NcmFit
Flags: Read / Write
“inequality-constraints-tot”
property“inequality-constraints-tot” NcmVector *
Inequality constraints tolerance.
Owner: NcmFit
Flags: Read / Write
“likelihood”
property“likelihood” NcmLikelihood *
Likelihood object.
Owner: NcmFit
Flags: Read / Write / Construct Only
“m2lnL-abstol”
property “m2lnL-abstol” double
Absolute tolarence in m2lnL.
Owner: NcmFit
Flags: Read / Write / Construct
Allowed values: >= 0
Default value: 0
“m2lnL-reltol”
property “m2lnL-reltol” double
Relative tolarence in m2lnL.
Owner: NcmFit
Flags: Read / Write / Construct
Allowed values: >= 0
Default value: 1e-08
“maxiter”
property “maxiter” guint
Maximum number of interations.
Owner: NcmFit
Flags: Read / Write / Construct
Default value: 100000
“mset”
property“mset” NcmMSet *
Model set object.
Owner: NcmFit
Flags: Read / Write / Construct Only
“params-reltol”
property “params-reltol” double
Relative tolarence in fitted parameters.
Owner: NcmFit
Flags: Read / Write / Construct
Allowed values: >= 0
Default value: 1e-05
“state”
property“state” NcmFitState *
Fit state object.
Owner: NcmFit
Flags: Read / Write / Construct Only