Performs a nonlinear regression to fit a model using a version of the Levenberg-Marquardt algorithm.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
procedure(regression_function), | intent(in), | pointer | :: | fun |
A pointer to the regression_function to evaluate. |
|
real(kind=real64), | intent(in) | :: | x(:) |
The M-element array containing independent data. |
||
real(kind=real64), | intent(in) | :: | y(:) |
The M-element array containing dependent data. |
||
real(kind=real64), | intent(inout) | :: | params(:) |
On input, the N-element array containing the initial estimate of the model parameters. On output, the computed model parameters. |
||
real(kind=real64), | intent(out) | :: | ymod(:) |
An M-element array where the modeled dependent data will be written. |
||
real(kind=real64), | intent(out) | :: | resid(:) |
An M-element array where the model residuals will be written. |
||
real(kind=real64), | intent(in), | optional, | target | :: | weights(:) |
An optional M-element array allowing the weighting of individual points. |
real(kind=real64), | intent(in), | optional, | target | :: | maxp(:) |
An optional N-element array that can be used as upper limits on the parameter values. If no upper limit is requested for a particular parameter, utilize a very large value. The internal default is to utilize huge() as a value. |
real(kind=real64), | intent(in), | optional, | target | :: | minp(:) |
An optional N-element array that can be used as lower limits on the parameter values. If no lower limit is requested for a particalar parameter, utilize a very large magnitude, but negative, value. The internal default is to utilize -huge() as a value. |
type(regression_statistics), | intent(out), | optional | :: | stats(:) |
An optional N-element array that, if supplied, will be used to return statistics about the fit for each parameter. |
|
real(kind=real64), | intent(in), | optional | :: | alpha |
The significance level at which to evaluate the confidence intervals. The default value is 0.05 such that a 95% confidence interval is calculated. |
|
type(iteration_controls), | intent(in), | optional | :: | controls |
An optional input providing custom iteration controls. |
|
type(lm_solver_options), | intent(in), | optional | :: | settings |
An optional input providing custom settings for the solver. |
|
type(convergence_info), | intent(out), | optional, | target | :: | info |
An optional output that can be used to gain information about the iterative solution and the nature of the convergence. |
procedure(iteration_update), | intent(in), | optional, | pointer | :: | status |
An optional pointer to a routine that can be used to extract iteration information. |
class(errors), | intent(inout), | optional, | target | :: | err |
A mechanism for communicating errors and warnings to the caller. Possible warning and error codes are as follows. - FS_NO_ERROR: No errors encountered. - FS_ARRAY_SIZE_ERROR: Occurs if any of the arrays are not properly sized. - FS_MEMORY_ERROR: Occurs if there is a memory allocation error. - FS_UNDERDEFINED_PROBLEM_ERROR: Occurs if the problem posed is underdetetermined (M < N). - FS_TOLERANCE_TOO_SMALL_ERROR: Occurs if any supplied tolerances are too small to be practical. - FS_TOO_FEW_ITERATION_ERROR: Occurs if too few iterations are allowed. |