qr Interface

public interface qr

Computes the QR factorization of matrix . (Specification)

Summary

Compute the QR factorization of a real or complex matrix: , where is orthonormal and is upper-triangular. Matrix has size [m,n], with .

Description

This interface provides methods for computing the QR factorization of a matrix. Supported data types include real and complex. If a pre-allocated work space is provided, no internal memory allocations take place when using this interface.

Given k = min(m,n), one can write \cdot ). The user may want the full problem (provide shape(Q)==[m,m], shape(R)==[m,n]) or the reduced
problem only: (provide shape(Q)==[m,k], shape(R)==[k,n]).

Note

The solution is based on LAPACK's QR factorization (*GEQRF) and ordered matrix output (*ORGQR, *UNGQR).


Subroutines

private pure module subroutine stdlib_linalg_c_qr(a, q, r, overwrite_a, storage, err)

Arguments

Type IntentOptional Attributes Name
complex(kind=sp), intent(inout), target :: a(:,:)

Input matrix a[m,n]

complex(kind=sp), intent(out), contiguous, target :: q(:,:)

Orthogonal matrix Q ([m,m], or [m,k] if reduced)

complex(kind=sp), intent(out), contiguous, target :: r(:,:)

Upper triangular matrix R ([m,n], or [k,n] if reduced)

logical(kind=lk), intent(in), optional :: overwrite_a

[optional] Can A data be overwritten and destroyed?

complex(kind=sp), intent(out), optional, target :: storage(:)

[optional] Provide pre-allocated workspace, size to be checked with qr_space

type(linalg_state_type), intent(out), optional :: err

[optional] state return flag. On error if not requested, the code will stop

private pure module subroutine stdlib_linalg_d_qr(a, q, r, overwrite_a, storage, err)

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout), target :: a(:,:)

Input matrix a[m,n]

real(kind=dp), intent(out), contiguous, target :: q(:,:)

Orthogonal matrix Q ([m,m], or [m,k] if reduced)

real(kind=dp), intent(out), contiguous, target :: r(:,:)

Upper triangular matrix R ([m,n], or [k,n] if reduced)

logical(kind=lk), intent(in), optional :: overwrite_a

[optional] Can A data be overwritten and destroyed?

real(kind=dp), intent(out), optional, target :: storage(:)

[optional] Provide pre-allocated workspace, size to be checked with qr_space

type(linalg_state_type), intent(out), optional :: err

[optional] state return flag. On error if not requested, the code will stop

private pure module subroutine stdlib_linalg_s_qr(a, q, r, overwrite_a, storage, err)

Arguments

Type IntentOptional Attributes Name
real(kind=sp), intent(inout), target :: a(:,:)

Input matrix a[m,n]

real(kind=sp), intent(out), contiguous, target :: q(:,:)

Orthogonal matrix Q ([m,m], or [m,k] if reduced)

real(kind=sp), intent(out), contiguous, target :: r(:,:)

Upper triangular matrix R ([m,n], or [k,n] if reduced)

logical(kind=lk), intent(in), optional :: overwrite_a

[optional] Can A data be overwritten and destroyed?

real(kind=sp), intent(out), optional, target :: storage(:)

[optional] Provide pre-allocated workspace, size to be checked with qr_space

type(linalg_state_type), intent(out), optional :: err

[optional] state return flag. On error if not requested, the code will stop

private pure module subroutine stdlib_linalg_z_qr(a, q, r, overwrite_a, storage, err)

Arguments

Type IntentOptional Attributes Name
complex(kind=dp), intent(inout), target :: a(:,:)

Input matrix a[m,n]

complex(kind=dp), intent(out), contiguous, target :: q(:,:)

Orthogonal matrix Q ([m,m], or [m,k] if reduced)

complex(kind=dp), intent(out), contiguous, target :: r(:,:)

Upper triangular matrix R ([m,n], or [k,n] if reduced)

logical(kind=lk), intent(in), optional :: overwrite_a

[optional] Can A data be overwritten and destroyed?

complex(kind=dp), intent(out), optional, target :: storage(:)

[optional] Provide pre-allocated workspace, size to be checked with qr_space

type(linalg_state_type), intent(out), optional :: err

[optional] state return flag. On error if not requested, the code will stop