Processing math: 100%
fortran-lapack
All Classes Namespaces Files Functions Variables Pages
la_qr::qr Interface Reference

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)
 

Detailed Description

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.

Parameters
[in,out]aThe input matrix of size [m,n] . If overwrite_a is true, the contents of a may be modified during computation.
[out]qThe orthonormal matrix Q , with shape [m,m] for the full problem or [m,k] for the reduced problem.
[out]rThe 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.
Returns
The QR factorization matrices Q and R are returned in the corresponding output arguments.
Note
This subroutine uses LAPACK's QR decomposition algorithm *GEQRF to compute the factorization.
Warning
If 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.

Parameters
[in]aThe input matrix of size [m,n] . This matrix is used to determine the required workspace size.
[out]lworkThe 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.
Returns
The size of the workspace array lwork that should be allocated before calling the QR factorization routine.
Note
This subroutine is useful for preallocating memory for QR factorization in large systems.
Warning
If err is not provided, the function will stop execution on error.

Member Function/Subroutine Documentation

◆ la_c_qr()

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 )
Parameters
[in,out]aInput matrix a[m,n]
[out]qOrthogonal matrix Q ([m,m], or [m,k] if reduced)
[out]rUpper 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

◆ la_d_qr()

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 )
Parameters
[in,out]aInput matrix a[m,n]
[out]qOrthogonal matrix Q ([m,m], or [m,k] if reduced)
[out]rUpper 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

◆ la_q_qr()

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 )
Parameters
[in,out]aInput matrix a[m,n]
[out]qOrthogonal matrix Q ([m,m], or [m,k] if reduced)
[out]rUpper 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

◆ la_s_qr()

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 )
Parameters
[in,out]aInput matrix a[m,n]
[out]qOrthogonal matrix Q ([m,m], or [m,k] if reduced)
[out]rUpper 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

◆ la_w_qr()

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 )
Parameters
[in,out]aInput matrix a[m,n]
[out]qOrthogonal matrix Q ([m,m], or [m,k] if reduced)
[out]rUpper 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

◆ la_z_qr()

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 )
Parameters
[in,out]aInput matrix a[m,n]
[out]qOrthogonal matrix Q ([m,m], or [m,k] if reduced)
[out]rUpper 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

The documentation for this interface was generated from the following file: