naginterfaces.library.correg.linregm_​var_​add

naginterfaces.library.correg.linregm_var_add(q, p, x, wt=None, tol=1e-06)[source]

linregm_var_add adds a new independent variable to a general linear regression model.

For full information please refer to the NAG Library document for g02de

https://www.nag.com/numeric/nl/nagdoc_29.3/flhtml/g02/g02def.html

Parameters
qfloat, array-like, shape

If , must contain the results of the decomposition for the model with parameters as returned by linregm_fit() or a previous call to linregm_var_add.

If , the first column of should contain the values of the dependent variable, .

pfloat, array-like, shape

Contains further details of the decomposition used. The first elements of must contain the zeta values for the decomposition (see lapackeig.dgeqrf for details).

The first elements of array are provided by linregm_fit() or by previous calls to linregm_var_add.

xfloat, array-like, shape

, the new independent variable.

wtNone or float, array-like, shape , optional

If provided must contain the weights to be used with the model.

If , the th observation is not included in the model, in which case the effective number of observations is the number of observations with nonzero weights.

If is not provided the effective number of observations is .

tolfloat, optional

The value of is used to decide if the new independent variable is linearly related to independent variables already included in the model. If the new variable is linearly related then is not updated. The smaller the value of the stricter the criterion for deciding if there is a linear relationship.

Returns
qfloat, ndarray, shape

The results of the decomposition for the model with parameters:

the first column of contains the updated value of ;

the columns to are unchanged;

the first elements of column contain the new column of , while the remaining elements contain details of the matrix .

pfloat, ndarray, shape

The first elements of are unchanged and the th element contains the zeta value for .

rssfloat

The residual sum of squares for the new fitted model.

Note: this will only be valid if the model is of full rank, see Further Comments.

Raises
NagValueError
(errno )

On entry, .

Constraint: .

(errno )

On entry, and .

Constraint: .

(errno )

On entry, .

Constraint: or .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: .

(errno )

On entry, .

Constraint: , for .

Warns
NagAlgorithmicWarning
(errno )

variable is a linear combination of existing model terms.

Notes

In the NAG Library the traditional C interface for this routine uses a different algorithmic base. Please contact NAG if you have any questions about compatibility.

A linear regression model may be built up by adding new independent variables to an existing model. linregm_var_add updates the decomposition used in the computation of the linear regression model. The decomposition may come from linregm_fit() or a previous call to linregm_var_add. The general linear regression model is defined by

where

is a vector of observations on the dependent variable,

is an matrix of the independent variables of column rank ,

is a vector of length of unknown parameters,

and

is a vector of length of unknown random errors such that , where is a known diagonal matrix.

If , the identity matrix, then least squares estimation is used. If , then for a given weight matrix , weighted least squares estimation is used.

The least squares estimates, of the parameters minimize while the weighted least squares estimates, minimize .

The parameter estimates may be found by computing a decomposition of (or in the weighted case), i.e.,

where and is a upper triangular matrix and is an orthogonal matrix.

If is of full rank, then is the solution to

where (or ) and is the first elements of .

If is not of full rank a solution is obtained by means of a singular value decomposition (SVD) of .

To add a new independent variable, , and have to be updated. The matrix is found such that (or ) is upper triangular. The vector is then updated by multiplying by .

The new independent variable is tested to see if it is linearly related to the existing independent variables by checking that at least one of the values , for , is nonzero.

The new parameter estimates, , can then be obtained by a call to linregm_update().

The function can be used with , in which case and are initialized.

References

Draper, N R and Smith, H, 1985, Applied Regression Analysis, (2nd Edition), Wiley

Golub, G H and Van Loan, C F, 1996, Matrix Computations, (3rd Edition), Johns Hopkins University Press, Baltimore

Hammarling, S, 1985, The singular value decomposition in multivariate statistics, SIGNUM Newsl. (20(3)), 2–25

McCullagh, P and Nelder, J A, 1983, Generalized Linear Models, Chapman and Hall

Searle, S R, 1971, Linear Models, Wiley