qr_factor Interface

public interface qr_factor

Module Procedures

private subroutine qr_factor_no_pivot(a, tau, work, olwork, err)

Computes the QR factorization of an M-by-N matrix.

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(inout), dimension(:,:) :: a

On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau, represent the orthogonal matrix Q as a product of elementary reflectors.

real(kind=real64), intent(out), dimension(:) :: tau

A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.

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

An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork.

integer(kind=int32), intent(out), optional :: olwork

An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

class(errors), intent(inout), optional, target :: err

The error object to be updated.

private subroutine qr_factor_no_pivot_cmplx(a, tau, work, olwork, err)

Computes the QR factorization of an M-by-N matrix.

Arguments

Type IntentOptional Attributes Name
complex(kind=real64), intent(inout), dimension(:,:) :: a

On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau, represent the orthogonal matrix Q as a product of elementary reflectors.

complex(kind=real64), intent(out), dimension(:) :: tau

A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.

complex(kind=real64), intent(out), optional, target, dimension(:) :: work

An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork.

integer(kind=int32), intent(out), optional :: olwork

An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

class(errors), intent(inout), optional, target :: err

The error object to be updated.

private subroutine qr_factor_pivot(a, tau, jpvt, work, olwork, err)

Computes the QR factorization of an M-by-N matrix using column pivoting such that .

Arguments

Type IntentOptional Attributes Name
real(kind=real64), intent(inout), dimension(:,:) :: a

On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau, represent the orthogonal matrix Q as a product of elementary reflectors.

real(kind=real64), intent(out), dimension(:) :: tau

A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.

integer(kind=int32), intent(inout), dimension(:) :: jpvt

On input, an N-element array that if JPVT(I) .ne. 0, the I-th column of A is permuted to the front of A * P; if JPVT(I) = 0, the I-th column of A is a free column. On output, if JPVT(I) = K, then the I-th column of A * P was the K-th column of A.

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

An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork.

integer(kind=int32), intent(out), optional :: olwork

An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

class(errors), intent(inout), optional, target :: err

The error object to be updated.

private subroutine qr_factor_pivot_cmplx(a, tau, jpvt, work, olwork, rwork, err)

Computes the QR factorization of an M-by-N matrix using column pivoting such that .

Arguments

Type IntentOptional Attributes Name
complex(kind=real64), intent(inout), dimension(:,:) :: a

On input, the M-by-N matrix to factor. On output, the elements on and above the diagonal contain the MIN(M, N)-by-N upper trapezoidal matrix R (R is upper triangular if M >= N). The elements below the diagonal, along with the array tau, represent the orthogonal matrix Q as a product of elementary reflectors.

complex(kind=real64), intent(out), dimension(:) :: tau

A MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.

integer(kind=int32), intent(inout), dimension(:) :: jpvt

On input, an N-element array that if JPVT(I) .ne. 0, the I-th column of A is permuted to the front of A * P; if JPVT(I) = 0, the I-th column of A is a free column. On output, if JPVT(I) = K, then the I-th column of A * P was the K-th column of A.

complex(kind=real64), intent(out), optional, target, dimension(:) :: work

An optional input, that if provided, prevents any local memory allocation. If not provided, the memory required is allocated within. If provided, the length of the array must be at least olwork.

integer(kind=int32), intent(out), optional :: olwork

An optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.

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

An optional input, that if provided, prevents any local allocate of real-valued memory. If not provided, the memory required is allocated within. If provided, the length of the array must be at least 2*N.

class(errors), intent(inout), optional, target :: err

The error object to be updated.