fstats_regression Module



Contents


Variables

Type Visibility Attributes Name Initial
integer(kind=int32), public, parameter :: FS_LEVENBERG_MARQUARDT_UPDATE = 1
integer(kind=int32), public, parameter :: FS_NIELSEN_UPDATE = 3
integer(kind=int32), public, parameter :: FS_QUADRATIC_UPDATE = 2

Interfaces

interface

  • public subroutine iteration_update(iter, funvals, resid, params, step)

    Arguments

    Type IntentOptional Attributes Name
    integer(kind=int32), intent(in) :: iter
    real(kind=real64), intent(in) :: funvals(:)
    real(kind=real64), intent(in) :: resid(:)
    real(kind=real64), intent(in) :: params(:)
    real(kind=real64), intent(in) :: step(:)

interface

  • public subroutine regression_function(xdata, params, resid, stop)

    Arguments

    Type IntentOptional Attributes Name
    real(kind=real64), intent(in), dimension(:) :: xdata
    real(kind=real64), intent(in), dimension(:) :: params
    real(kind=real64), intent(out), dimension(:) :: resid
    logical, intent(out) :: stop

Derived Types

type, public ::  convergence_info

Provides information regarding convergence status.

Components

Type Visibility Attributes Name Initial
logical, public :: converge_on_gradient

True if convergence on the gradient was achieved; else, false.

logical, public :: converge_on_residual_parameter

True if convergence on the residual error parameter was achieved; else, false.

logical, public :: converge_on_solution_change

True if convergence on the change in solution was achieved; else, false.

integer(kind=int32), public :: function_evaluation_count

The function evaluation count.

real(kind=real64), public :: gradient_value

The value of the gradient test parameter.

integer(kind=int32), public :: iteration_count

The iteration count.

logical, public :: reach_function_evaluation_limit

True if the solution did not converge in the allowed number of function evaluations.

logical, public :: reach_iteration_limit

True if the solution did not converge in the allowed number of iterations.

real(kind=real64), public :: residual_value

The value of the residual error parameter.

real(kind=real64), public :: solution_change_value

The value of the change in solution parameter.

logical, public :: user_requested_stop

True if the user requested the stop; else, false.

type, public ::  iteration_controls

Provides a collection of iteration control parameters.

Components

Type Visibility Attributes Name Initial
real(kind=real64), public :: change_in_solution_tolerance

Defines a tolerance on the change in parameter values.

real(kind=real64), public :: gradient_tolerance

Defines a tolerance on the gradient of the fitted function.

real(kind=real64), public :: iteration_improvement_tolerance

Defines a tolerance to ensure adequate improvement on each iteration.

integer(kind=int32), public :: max_function_evaluations

Defines the maximum number of function evaluations allowed.

integer(kind=int32), public :: max_iteration_between_updates

Defines how many iterations can pass before a re-evaluation of the Jacobian matrix is forced.

integer(kind=int32), public :: max_iteration_count

Defines the maximum number of iterations allowed.

real(kind=real64), public :: residual_tolerance

Defines a tolerance on the metric associated with the residual error.

Type-Bound Procedures

procedure , public :: set_to_default => lm_set_default_tolerances Subroutine

type, public ::  lm_solver_options

Options to control the Levenberg-Marquardt solver.

Components

Type Visibility Attributes Name Initial
real(kind=real64), public :: damping_decrease_factor

The factor to use when decreasing the damping parameter.

real(kind=real64), public :: damping_increase_factor

The factor to use when increasing the damping parameter.

real(kind=real64), public :: finite_difference_step_size

The step size used for the finite difference calculations of the Jacobian matrix.

integer(kind=int32), public :: method

The solver method to utilize. - FS_LEVENBERG_MARQUARDT_UPDATE: - FS_QUADRATIC_UPDATE: - FS_NIELSEN_UDPATE:

Type-Bound Procedures

procedure , public :: set_to_default => lm_set_default_settings Subroutine

type, public ::  regression_statistics

A container for regression-related statistical information.

Components

Type Visibility Attributes Name Initial
real(kind=real64), public :: confidence_interval

The confidence interval for the parameter at the level determined by the regression process.

Read more…
real(kind=real64), public :: probability

The probability that the coefficient is not statistically important. A statistically important coefficient will have a low probability (p-value), typically 0.05 or lower; however, a p-value of up to ~0.2 may be acceptable dependent upon the problem. Typically any p-value larger than ~0.2 indicates the parameter is not statistically important for the model.

Read more…
real(kind=real64), public :: standard_error

The standard error for the model coefficient.

Read more…
real(kind=real64), public :: t_statistic

The T-statistic for the model coefficient.

Read more…

Functions

public function adjusted_r_squared(p, x, xm, err) result(rst)

Computes the adjusted R-squared value for a data set.

Read more…

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: p

The number of variables.

real(kind=real64), intent(in) :: x(:)

An N-element array containing the dependent variables from the data set.

real(kind=real64), intent(in) :: xm(:)

An N-element array containing the corresponding modeled values.

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 x and xm are not the same size.

Return Value real(kind=real64)

The result.

public function calculate_regression_statistics(resid, params, c, alpha, err) result(rst)

Computes statistics for the quality of fit for a regression model.

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in) :: resid(:)

An M-element array containing the model residual errors.

real(kind=real64), intent(in) :: params(:)

An N-element array containing the model parameters.

real(kind=real64), intent(in) :: c(:,:)

The N-by-N covariance matrix.

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.

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 c is not sized correctly. - FS_INVALID_INPUT_ERROR: Occurs if order is less than 1. - FS_MEMORY_ERROR: Occurs if there is a memory allocation error.

Return Value type(regression_statistics), allocatable, (:)

A regression_statistics object containing the analysis results.

public pure function correlation(x, y) result(rst)

Computes the sample correlation coefficient (an estimate to the population Pearson correlation) as follows.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in), dimension(:) :: x

The first N-element data set.

real(kind=real64), intent(in), dimension(size(x)) :: y

The second N-element data set.

Return Value real(kind=real64)

The correlation coefficient.

public function r_squared(x, xm, err) result(rst)

Computes the R-squared value for a data set.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in) :: x(:)

An N-element array containing the dependent variables from the data set.

real(kind=real64), intent(in) :: xm(:)

An N-element array containing the corresponding modeled values.

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 x and xm are not the same size.

Return Value real(kind=real64)

The result.


Subroutines

public subroutine covariance_matrix(x, c, err)

Computes the covariance matrix where and is computed by design_matrix.

Read more…

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(in) :: x(:,:)

An M-by-N matrix containing the formatted independent data matrix as computed by design_matrix.

real(kind=real64), intent(out) :: c(:,:)

The N-by-N covariance matrix.

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 matrices are not sized correctly. - FS_MEMORY_ERROR: Occurs if there is a memory allocation error.

public subroutine design_matrix(order, intercept, x, c, err)

Computes the design matrix for the linear least-squares regression problem of , where is the matrix computed here, is the vector of coefficients to be determined, and is the vector of measured dependent variables.

Read more…

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: order

The order of the equation to fit. This value must be at least one (linear equation), but can be higher as desired.

logical, intent(in) :: intercept

Set to true if the intercept is being computed as part of the regression; else, false.

real(kind=real64), intent(in) :: x(:)

An N-element array containing the independent variable measurement points.

real(kind=real64), intent(out) :: c(:,:)

An N-by-K matrix where the results will be written. K must equal order + 1 in the event intercept is true; however, if intercept is false, K must equal order.

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 c is not properly sized. - FS_INVALID_INPUT_ERROR: Occurs if order is less than 1.

public subroutine jacobian(fun, xdata, params, jac, stop, f0, f1, step, err)

Computes the Jacobian matrix for a nonlinear regression problem.

Arguments

Type IntentOptional Attributes Name
procedure(regression_function), intent(in), pointer :: fun

A pointer to the regression_function to evaluate.

real(kind=real64), intent(in) :: xdata(:)

The M-element array containing x-coordinate data.

real(kind=real64), intent(in) :: params(:)

The N-element array containing the model parameters.

real(kind=real64), intent(out) :: jac(:,:)

The M-by-N matrix where the Jacobian will be written.

logical, intent(out) :: stop

A value that the user can set in fun forcing the evaluation process to stop prior to completion.

real(kind=real64), intent(in), optional, target :: f0(:)

An optional M-element array containing the model values using the current parameters as defined in m. This input can be used to prevent the routine from performing a function evaluation at the model parameter state defined in params.

real(kind=real64), intent(out), optional, target :: f1(:)

An optional M-element workspace array used for function evaluations.

real(kind=real64), intent(in), optional :: step

The differentiation step size. The default is the square root of machine precision.

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.

public subroutine linear_least_squares(order, intercept, x, y, coeffs, ymod, resid, stats, alpha, err)

Computes a linear least-squares regression to fit a set of data.

Read more…

Arguments

Type IntentOptional Attributes Name
integer(kind=int32), intent(in) :: order

The order of the equation to fit. This value must be at least one (linear equation), but can be higher as desired, as long as there is sufficient data.

logical, intent(in) :: intercept

Set to true if the intercept is being computed as part of the regression; else, false.

real(kind=real64), intent(in) :: x(:)

An N-element array containing the independent variable measurement points.

real(kind=real64), intent(in) :: y(:)

An N-element array containing the dependent variable measurement points.

real(kind=real64), intent(out) :: coeffs(:)

An ORDER+1 element array where the coefficients will be written.

real(kind=real64), intent(out) :: ymod(:)

An N-element array where the modeled data will be written.

real(kind=real64), intent(out) :: resid(:)

An N-element array where the residual error data will be written (modeled - actual).

type(regression_statistics), intent(out), optional :: stats(:)

An M-element array of regression_statistics items where M = ORDER + 1 when intercept is set to true; however, if intercept is set to false, M = ORDER.

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.

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 approriately sized. - FS_INVALID_INPUT_ERROR: Occurs if order is less than 1. - FS_MEMORY_ERROR: Occurs if there is a memory allocation error.

public subroutine nonlinear_least_squares(fun, x, y, params, ymod, resid, weights, maxp, minp, stats, alpha, controls, settings, info, status, err)

Performs a nonlinear regression to fit a model using a version of the Levenberg-Marquardt algorithm.

Arguments

Type IntentOptional 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.