libkpl
6.0
A Library for Graphical Presentation of Data Sets and Functions
|
Parameter fit class. More...
#include <lmfit.h>
Signals | |
void | updateMessage (const QString &message) |
Emitted when a message should be displayed. | |
Public Member Functions | |
LMFit (const QList< const ArrayItem * > *aList, QList< FunItem * > *fList, const bool *bFit, double *par, QObject *parent=nullptr, const double *sigma=nullptr) | |
Constructor. | |
~LMFit () | |
Destructor. | |
bool | linFit (int n, int nPar, double par[], double resVector[], double fJac[], int iPvt[]) |
Performs general linear least square parameter fit. | |
void | lmdif (int n, int nPar, double par[], double resVector[], double tol, int itMax, int *info, double fJac[], int iPvt[]) |
Performs nonlinear Levenberg-Marquardt parameter fit. | |
void | setSigma (const double *sigma) |
Sets data error array. | |
void | setUserBreak (bool userBreak=true) |
Sets condition to stop calculations. | |
bool | userBreak () const |
Returns current break condition for calculations. | |
Static Public Member Functions | |
static double | enorm (int n, const double x[]) |
Calculates the euclidean norm of x. | |
static double | igamc (double a, double x) |
Calculates complemented incomplete gamma integral. | |
static int | minv (double a[], double x[], int n) |
Finds the inverse of a square matrix. | |
Protected Attributes | |
bool | m_userBreak |
const bool * | m_bFit |
double * | m_par |
const double * | m_sigma |
const QList< const ArrayItem * > * | m_aList |
QList< FunItem * > * | m_fList |
Parameter fit class.
Provides methods for nonlinear Levenberg-Marquardt parameter fits and general linear least square parameter fits.
The functions lmdif, lmpar, fdjac2, qrsolv, qrfac, and enorm use code from the public domain Fortran version of Argonne National Laboratories MINPACK, C translation by Steve Moshier. Modifications by Werner Stille.
The functions minv, igamc, igam, simq, and mtransp are from the Cephes Math Library Release 2.7, Copyright by Stephen L. Moshier.
The function lgam (original name: gsl_sf_lngamma_impl, Author: G. Jungman) is from the GNU Scientific Library. Modifications by Werner Stille.
LMFit::LMFit | ( | const QList< const ArrayItem * > * | aList, |
QList< FunItem * > * | fList, | ||
const bool * | bFit, | ||
double * | par, | ||
QObject * | parent = nullptr, | ||
const double * | sigma = nullptr ) |
Constructor.
aList | list of array items. |
fList | list of function items. |
bFit | array containing true values for parameters to be fitted. |
par | array of parameters. |
sigma | array for data errors. |
parent | pointer to parent widget. |
|
static |
Calculates the euclidean norm of x.
n | length of x. |
x | input array. |
|
static |
Calculates complemented incomplete gamma integral.
a | argument of gamma function. |
x | upper boundary of integral. |
bool LMFit::linFit | ( | int | n, |
int | nPar, | ||
double | par[], | ||
double | resVector[], | ||
double | fJac[], | ||
int | iPvt[] ) |
Performs general linear least square parameter fit.
n | number of data points. |
nPar | number of parameters to be fitted. |
par | array of parameters to be fitted. |
resVector | array for weighted residuals. |
fJac | array for calculating the Jacobian matrix. |
iPvt | array defining a permutation matrix for calculating the Jacobian matrix. |
void LMFit::lmdif | ( | int | n, |
int | nPar, | ||
double | par[], | ||
double | resVector[], | ||
double | tol, | ||
int | itMax, | ||
int * | info, | ||
double | fJac[], | ||
int | iPvt[] ) |
Performs nonlinear Levenberg-Marquardt parameter fit.
n | number of data points. |
nPar | number of parameters to be fitted. |
par | array of parameters to be fitted. |
resVector | array for weighted residuals. |
tol | tolerance value for stopping nonlinear fit iterations. |
itMax | maximum number of nonlinear fit iterations. |
info | error code. |
fJac | array for calculating the Jacobian matrix. |
iPvt | array defining a permutation matrix for calculating the Jacobian matrix. |
|
static |
Finds the inverse of a square matrix.
a | input matrix. |
x | output matrix. |
n | numbers of rows and columns of input matrix. |
|
inline |
Sets data error array.
sigma | array for data errors. |
|
inline |
Sets condition to stop calculations.
userBreak | true to break calculations. |
|
signal |
Emitted when a message should be displayed.
message | message text |