fortran-lapack
|
Compute the QR factorization of a matrix. More...
Public Member Functions | |
pure subroutine | la_s_qr (a, q, r, overwrite_a, storage, err) |
pure subroutine | la_d_qr (a, q, r, overwrite_a, storage, err) |
pure subroutine | la_q_qr (a, q, r, overwrite_a, storage, err) |
pure subroutine | la_c_qr (a, q, r, overwrite_a, storage, err) |
pure subroutine | la_z_qr (a, q, r, overwrite_a, storage, err) |
pure subroutine | la_w_qr (a, q, r, overwrite_a, storage, err) |
Compute the QR factorization of a matrix.
This subroutine computes the QR factorization of a real or complex matrix:
A = Q \cdot R
where A is a matrix of size [m,n] with m \ge n , Q is an orthonormal matrix, and R is an upper-triangular matrix. The result is returned in the output matrices Q and R , which are of the same type and kind as A . The function can solve both full and reduced problems depending on the sizes of the matrices.
[in,out] | a | The input matrix of size [m,n] . If overwrite_a is true, the contents of a may be modified during computation. |
[out] | q | The orthonormal matrix Q , with shape [m,m] for the full problem or [m,k] for the reduced problem. |
[out] | r | The upper-triangular matrix R , with shape [m,n] for the full problem or [k,n] for the reduced problem. |
[in] | overwrite_a | (Optional) A logical flag that determines whether matrix a may be overwritten. Default is false. |
[in,out] | storage | (Optional) Pre-allocated workspace array for intermediate calculations. Its size is checked with qr_space . |
[out] | err | (Optional) A state return flag. If an error occurs and err is not provided, the function will stop execution. |
*GEQRF
to compute the factorization.overwrite_a
is enabled, the original contents of a
may be lost during computation.Get the required workspace size for QR factorization.
This subroutine returns the minimum workspace size required for performing QR factorization. It computes the size of the workspace array that is necessary for both QR factorization and solving the reduced problem, depending on the input matrix.
[in] | a | The input matrix of size [m,n] . This matrix is used to determine the required workspace size. |
[out] | lwork | The minimum size of the workspace array needed for QR operations. |
[out] | err | (Optional) A state return flag. If an error occurs and err is not provided, the function will stop execution. |
lwork
that should be allocated before calling the QR factorization routine.err
is not provided, the function will stop execution on error. pure subroutine la_qr::qr::la_c_qr | ( | complex(sp), dimension(:,:), intent(inout), target | a, |
complex(sp), dimension(:,:), intent(out), target, contiguous | q, | ||
complex(sp), dimension(:,:), intent(out), target, contiguous | r, | ||
logical(lk), intent(in), optional | overwrite_a, | ||
complex(sp), dimension(:), intent(inout), optional, target | storage, | ||
type(la_state), intent(out), optional | err ) |
[in,out] | a | Input matrix a[m,n] |
[out] | q | Orthogonal matrix Q ([m,m], or [m,k] if reduced) |
[out] | r | Upper triangular matrix R ([m,n], or [k,n] if reduced) |
[in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
[in,out] | storage | [optional] Provide pre-allocated workspace, size to be checked with qr_space |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
pure subroutine la_qr::qr::la_d_qr | ( | real(dp), dimension(:,:), intent(inout), target | a, |
real(dp), dimension(:,:), intent(out), target, contiguous | q, | ||
real(dp), dimension(:,:), intent(out), target, contiguous | r, | ||
logical(lk), intent(in), optional | overwrite_a, | ||
real(dp), dimension(:), intent(inout), optional, target | storage, | ||
type(la_state), intent(out), optional | err ) |
[in,out] | a | Input matrix a[m,n] |
[out] | q | Orthogonal matrix Q ([m,m], or [m,k] if reduced) |
[out] | r | Upper triangular matrix R ([m,n], or [k,n] if reduced) |
[in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
[in,out] | storage | [optional] Provide pre-allocated workspace, size to be checked with qr_space |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
pure subroutine la_qr::qr::la_q_qr | ( | real(qp), dimension(:,:), intent(inout), target | a, |
real(qp), dimension(:,:), intent(out), target, contiguous | q, | ||
real(qp), dimension(:,:), intent(out), target, contiguous | r, | ||
logical(lk), intent(in), optional | overwrite_a, | ||
real(qp), dimension(:), intent(inout), optional, target | storage, | ||
type(la_state), intent(out), optional | err ) |
[in,out] | a | Input matrix a[m,n] |
[out] | q | Orthogonal matrix Q ([m,m], or [m,k] if reduced) |
[out] | r | Upper triangular matrix R ([m,n], or [k,n] if reduced) |
[in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
[in,out] | storage | [optional] Provide pre-allocated workspace, size to be checked with qr_space |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
pure subroutine la_qr::qr::la_s_qr | ( | real(sp), dimension(:,:), intent(inout), target | a, |
real(sp), dimension(:,:), intent(out), target, contiguous | q, | ||
real(sp), dimension(:,:), intent(out), target, contiguous | r, | ||
logical(lk), intent(in), optional | overwrite_a, | ||
real(sp), dimension(:), intent(inout), optional, target | storage, | ||
type(la_state), intent(out), optional | err ) |
[in,out] | a | Input matrix a[m,n] |
[out] | q | Orthogonal matrix Q ([m,m], or [m,k] if reduced) |
[out] | r | Upper triangular matrix R ([m,n], or [k,n] if reduced) |
[in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
[in,out] | storage | [optional] Provide pre-allocated workspace, size to be checked with qr_space |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
pure subroutine la_qr::qr::la_w_qr | ( | complex(qp), dimension(:,:), intent(inout), target | a, |
complex(qp), dimension(:,:), intent(out), target, contiguous | q, | ||
complex(qp), dimension(:,:), intent(out), target, contiguous | r, | ||
logical(lk), intent(in), optional | overwrite_a, | ||
complex(qp), dimension(:), intent(inout), optional, target | storage, | ||
type(la_state), intent(out), optional | err ) |
[in,out] | a | Input matrix a[m,n] |
[out] | q | Orthogonal matrix Q ([m,m], or [m,k] if reduced) |
[out] | r | Upper triangular matrix R ([m,n], or [k,n] if reduced) |
[in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
[in,out] | storage | [optional] Provide pre-allocated workspace, size to be checked with qr_space |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |
pure subroutine la_qr::qr::la_z_qr | ( | complex(dp), dimension(:,:), intent(inout), target | a, |
complex(dp), dimension(:,:), intent(out), target, contiguous | q, | ||
complex(dp), dimension(:,:), intent(out), target, contiguous | r, | ||
logical(lk), intent(in), optional | overwrite_a, | ||
complex(dp), dimension(:), intent(inout), optional, target | storage, | ||
type(la_state), intent(out), optional | err ) |
[in,out] | a | Input matrix a[m,n] |
[out] | q | Orthogonal matrix Q ([m,m], or [m,k] if reduced) |
[out] | r | Upper triangular matrix R ([m,n], or [k,n] if reduced) |
[in] | overwrite_a | [optional] Can A data be overwritten and destroyed? |
[in,out] | storage | [optional] Provide pre-allocated workspace, size to be checked with qr_space |
[out] | err | [optional] state return flag. On error if not requested, the code will stop |