![]() |
libSUFR
a LIBrary of Some Useful Fortran Routines
|
Procedures to fit functions to data. More...
Functions/Subroutines | |
subroutine, public | linear_fit_yerr (ndat, xdat, ydat, yerr, ncoef, fcoef, freecoef, ncov, covmat, chisq, myfunc) |
Fit the model coefficients using a chi^2 method. | |
subroutine, public | nonlin_fit_yerr (ny, ndat, xdat, ydat, yerr, ncoef, fcoef, freecoef, nmat, covmat, curvmat, chisq, iterstat, myfunc) |
Levenberg-Marquardt method to reduce the value of chi-squared of a fit of a set of data points with errors in Y to a non-linear function. | |
subroutine, public | nonlin_fit_example_myfunc (ny, xdat, ncoef, fcoef, ydat, dyda) |
Dummy example function myFunc for nonlin_fit_yerr(): return the value and partial derivatives. | |
pure subroutine, public | basefunc_polynomial (xval, ncoef, basefunc) |
User-provided base function to fit a polynomial to. | |
Procedures to fit functions to data.
pure subroutine, public sufr_fitting::basefunc_polynomial | ( | real(double), intent(in) | xval, |
integer, intent(in) | ncoef, | ||
real(double), dimension(ncoef), intent(out) | basefunc ) |
User-provided base function to fit a polynomial to.
xVal | X value to evaluate yVal = basefunc_polynimial(xVal) for |
nCoef | Number of coefficients for the base function |
baseFunc | Base function i for coefficent i. For a polynomial, BF(i) = x**(i-1) (output) |
Definition at line 585 of file fitting.f90.
References sufr_kinds::double.
subroutine, public sufr_fitting::linear_fit_yerr | ( | integer, intent(in) | ndat, |
real(double), dimension(ndat), intent(in) | xdat, | ||
real(double), dimension(ndat), intent(in) | ydat, | ||
real(double), dimension(ndat), intent(in) | yerr, | ||
integer, intent(in) | ncoef, | ||
real(double), dimension(ncoef), intent(inout) | fcoef, | ||
integer, dimension(ncoef), intent(in) | freecoef, | ||
integer, intent(in) | ncov, | ||
real(double), dimension(ncov,ncov), intent(out) | covmat, | ||
real(double), intent(out) | chisq, | ||
external | myfunc ) |
Fit the model coefficients using a chi^2 method.
Do a chi^2-minimisation fit of nDat data points xDat and yDat, the latter with errors yErr, to function func. The fit uses the nCoef coefficients stored in fCoef() that should be varied if the corresponding entry in freeCoef=1, fixed if =0.
The fitting function must be provided as the subroutine myFunc which returns a value y = myFunc(x, basefunc, nCoef), where y is a linear combination of x, i.e. y = SUM(i=1,nCoef) fCoef(i) * basefunc(i)(x). For example, for a polynomial, basefunc(i) = x**(i-1).
The routine returns the fit coefficients in fCoef(), the variance-covariance matrix in covMat sith size nCov, and the chi^2 value in chiSq.
nDat | Number of data points |
xDat | X values of the data (nDat) |
yDat | Y values of the data (nDat) |
yErr | Errors (standard deviations) of the y values (nDat) |
nCoef | Number of coefficients of the fitting function |
fCoef | Coefficients of the fitting function (nCoef) (output) |
freeCoef | Fix coefficient coef(i) if freeCoef(i)=0, otherwise let if vary freely and fit it (nCoef) |
nCov | Size of both dimensions of covMat() |
covMat | Covariance matrix (nCov,nCov) (output) |
chiSq | Chi^2 - chi squared (output) |
myFunc | External subroutine that describes the model value of Y for given value X |
Definition at line 69 of file fitting.f90.
References sufr_kinds::double, sufr_system::quit_program_error(), and sufr_system::warn().
subroutine, public sufr_fitting::nonlin_fit_example_myfunc | ( | integer, intent(in) | ny, |
real(double), intent(in) | xdat, | ||
integer, intent(in) | ncoef, | ||
real(double), dimension(ncoef), intent(in) | fcoef, | ||
real(double), dimension(ny), intent(out) | ydat, | ||
real(double), dimension(ny,ncoef), intent(out) | dyda ) |
Dummy example function myFunc for nonlin_fit_yerr(): return the value and partial derivatives.
nY | Number of Y values for each X value (normally 1, but e.g. 2 for a sky position with two coordinates) |
xDat | Input X value for the current data point |
nCoef | Number of coefficients |
fCoef | Vector of coefficients that describe the function |
yDat | Y values for the current data point (output) |
dyda | Partial derivatives for yDat: 1: dy/dfCoef(1), ..., n: dy/dfCoef(n) (output) |
Definition at line 406 of file fitting.f90.
References sufr_kinds::double, sufr_system::quit_program_error(), sufr_system::swapdbl(), and sufr_system::warn().
subroutine, public sufr_fitting::nonlin_fit_yerr | ( | integer, intent(in) | ny, |
integer, intent(in) | ndat, | ||
real(double), dimension(ndat), intent(in) | xdat, | ||
real(double), dimension(ny,ndat), intent(in) | ydat, | ||
real(double), dimension(ndat), intent(in) | yerr, | ||
integer, intent(in) | ncoef, | ||
real(double), dimension(ncoef), intent(inout) | fcoef, | ||
integer, dimension(ncoef), intent(in) | freecoef, | ||
integer, intent(in) | nmat, | ||
real(double), dimension(nmat,nmat), intent(out) | covmat, | ||
real(double), dimension(nmat,nmat), intent(out) | curvmat, | ||
real(double), intent(out) | chisq, | ||
real(double), intent(inout) | iterstat, | ||
external | myfunc ) |
Levenberg-Marquardt method to reduce the value of chi-squared of a fit of a set of data points with errors in Y to a non-linear function.
You will need to call this routine repeatedly, with different values for iterStat, until chiSq no longer decreases (significantly)
nY | Number of Y values for each X value (normally 1, but e.g. 2 for a sky position with two coordinates) |
nDat | Number of data points in X and Y |
xDat | X data points to fit |
yDat | Y data points to fit |
yErr | Errors (standard deviations) for yDat |
nCoef | Number of coefficients used to describe the non-linear function myFunc |
fCoef | Coefficients for the non-linear function myFunc, updated after each call |
freeCoef | Determines which coefficients for the non-linear function myFunc should be fitted for. Set freeCoef(i) = 0 in order to keep fCoef(i) fixed |
nMat | Size of the matrices; > nCoef |
covMat | Variance-covariance matrix - returned on last call with iterStat. Fixed parameters return zero (co)variances |
curvMat | Hessian/curvature matrix - double partial derivative of chi squared to two coefficients fCoef |
chiSq | Chi squared |
iterStat | Iteration status; Set to <0 for initialisation in the first call, and let it vary in subsequent calls, until the fit converges. After that, set to 0 to return the variance-covariance and curvature matrices on the final call. IterStat decreases 10x if chiSq becomes smaller, and increases 10x otherwise. |
myFunc | External subroutine that describes the model value of Y and partial derivatives dY/dXi for given value X and function coefficients fCoef |
Definition at line 191 of file fitting.f90.
References sufr_numerics::deq0(), and sufr_kinds::double.