linalg 1.8.2
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
Loading...
Searching...
No Matches
linalg::qr_factor Interface Reference

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

Public Member Functions

 qr_factor_no_pivot
 
 qr_factor_no_pivot_cmplx
 
 qr_factor_pivot
 
 qr_factor_pivot_cmplx
 

Detailed Description

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

Syntax 1
subroutine qr_factor(real(real64) a(:,:), real(real64) tau(:), optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine qr_factor(complex(real64) a(:,:), complex(real64) tau(:), optional complex(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
Parameters
[in,out]aOn 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.
[out]tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
[out]workAn 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.
[out]olworkAn optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.
[in,out]errAn optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
  • LA_ARRAY_SIZE_ERROR: Occurs if tau or work are not sized appropriately.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
Remarks
QR factorization without pivoting is best suited to solving an overdetermined system in least-squares terms, or to solve a normally defined system. To solve an underdetermined system, it is recommended to use either LQ factorization, or a column-pivoting based QR factorization.
Notes
This routine utilizes the LAPACK routine DGEQRF (ZGEQRF for the complex case).
Syntax 2
Computes the QR factorization of an M-by-N matrix with column pivoting such that \( A P = Q R \).
subroutine qr_factor(real(real64) a(:,:), real(real64) tau(:), integer(int32) jpvt(:), optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine qr_factor(complex(real64) a(:,:), complex(real64) tau(:), integer(int32) jpvt(:), optional complex(real64) work(:), optional integer(int32) olwork, optional real(real64) rwork(:), optional class(errors) err)
Parameters
[in,out]aOn 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.
[out]tauA MIN(M, N)-element array used to store the scalar factors of the elementary reflectors.
[in,out]jpvtOn 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.
[out]workAn 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.
[out]olworkAn optional output used to determine workspace size. If supplied, the routine determines the optimal size for work, and returns without performing any actual calculations.
[out]rworkAn 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.
[in,out]errAn optional errors-based object that if provided can be used to retrieve information relating to any errors encountered during execution. If not provided, a default implementation of the errors class is used internally to provide error handling. Possible errors and warning messages that may be encountered are as follows.
  • LA_ARRAY_SIZE_ERROR: Occurs if any of the input arrays are not sized appropriately.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
Notes
This routine utilizes the LAPACK routine DGEQP3 (ZGEQP3 for the complex case).
Usage
The following example illustrates the solution of a system of equations using QR factorization.
program example
use iso_fortran_env, only : real64, int32
use linalg
implicit none
! Local Variables
real(real64) :: a(3,3), tau(3), b(3)
integer(int32) :: i, pvt(3)
! Build the 3-by-3 matrix A.
! | 1 2 3 |
! A = | 4 5 6 |
! | 7 8 0 |
a = reshape( &
[1.0d0, 4.0d0, 7.0d0, 2.0d0, 5.0d0, 8.0d0, 3.0d0, 6.0d0, 0.0d0], &
[3, 3])
! Build the right-hand-side vector B.
! | -1 |
! b = | -2 |
! | -3 |
b = [-1.0d0, -2.0d0, -3.0d0]
! The solution is:
! | 1/3 |
! x = | -2/3 |
! | 0 |
! Compute the QR factorization, using pivoting
pvt = 0 ! Zero every entry in order not to lock any column in place
call qr_factor(a, tau, pvt)
! Compute the solution. The results overwrite b.
call solve_qr(a, tau, pvt, b)
! Display the results.
print '(A)', "QR Solution: X = "
print '(F8.4)', (b(i), i = 1, size(b))
! Notice, QR factorization without pivoting could be accomplished in the
! same manner. The only difference is to omit the PVT array (column pivot
! tracking array).
end program
Computes the QR factorization of an M-by-N matrix.
Definition linalg.f90:1121
Solves a system of M QR-factored equations of N unknowns.
Definition linalg.f90:2557
Provides a set of common linear algebra routines.
Definition linalg.f90:145
The above program produces the following output.
QR Solution: X =
0.3333
-0.6667
0.0000
See Also

Definition at line 1121 of file linalg.f90.


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