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::form_qr Interface Reference

Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm. More...

Public Member Functions

 form_qr_no_pivot
 
 form_qr_no_pivot_cmplx
 
 form_qr_pivot
 
 form_qr_pivot_cmplx
 

Detailed Description

Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm.

Syntax 1
subroutine form_qr(real(real64) r(:,:), real(real64) tau(:), real(real64) q(:,:), optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine form_qr(complex(real64) r(:,:), complex(real64) tau(:), complex(real64) q(:,:), optional complex(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
Parameters
[in,out]rOn input, an M-by-N matrix where the elements below the diagonal contain the elementary reflectors generated from the QR factorization. On and above the diagonal, the matrix contains the matrix R. On output, the elements below the diagonal are zeroed such that the remaining matrix is simply the M-by-N matrix R.
[in]tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
[out]qAn M-by-M matrix where the full orthogonal matrix Q will be written. In the event that M > N, Q may be supplied as M-by-N, and therefore only return the useful submatrix Q1 (Q = [Q1, Q2]) as the factorization can be written as Q * R = [Q1, Q2] * [R1; 0].
[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 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 DORGQR (ZUNQR in the complex case).
Syntax 2
subroutine form_qr(real(real64) r(:,:), real(real64) tau(:), integer(int32) pvt(:), real(real64) q(:,:), real(real64) p(:,:), optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine form_qr(complex(real64) r(:,:), complex(real64) tau(:), integer(int32) pvt(:), complex(real64) q(:,:), complex(real64) p(:,:), optional complex(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
Parameters
[in,out]rOn input, an M-by-N matrix where the elements below the diagonal contain the elementary reflectors generated from the QR factorization. On and above the diagonal, the matrix contains the matrix R. On output, the elements below the diagonal are zeroed such that the remaining matrix is simply the M-by-N matrix R.
[in]tauA MIN(M, N)-element array containing the scalar factors of each elementary reflector defined in r.
[in]pvtAn N-element column pivot array as returned by the QR factorization.
[out]qAn M-by-M matrix where the full orthogonal matrix Q will be written. In the event that M > N, Q may be supplied as M-by-N, and therefore only return the useful submatrix Q1 (Q = [Q1, Q2]) as the factorization can be written as Q * R = [Q1, Q2] * [R1; 0].
[out]pAn N-by-N matrix where the pivot matrix will be written.
[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 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 DORGQR (ZUNQR in the complex case).
Usage
The following example illustrates how to explicitly form the Q and R matrices from the output of qr_factor, and then use the resulting matrices to solve a system of linear equations.
program example
use iso_fortran_env, only : real64, int32
use linalg
implicit none
! Variables
real(real64) :: a(3,3), b(3), q(3,3), tau(3)
integer(int32) :: i
! 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 without column pivoting
call qr_factor(a, tau)
! Build Q and R. A is overwritten with R
call form_qr(a, tau, q)
! As this system is square, matrix R is upper triangular. Also, Q is
! always orthogonal such that it's inverse and transpose are equal. As the
! system is now factored, its form is: Q * R * X = B. Solving this system
! is then as simple as solving the upper triangular system:
! R * X = Q**T * B.
! Compute Q**T * B, and store the results in B
b = matmul(transpose(q), b)
! Solve the upper triangular system R * X = Q**T * B for X
call solve_triangular_system(.true., .false., .true., a, b)
! Display the results
print '(A)', "QR Solution: X = "
print '(F8.4)', (b(i), i = 1, size(b))
! Notice, QR factorization with column pivoting could be accomplished via
! a similar approach, but the column pivoting would need to be accounted
! for by noting that Q * R = A * P, where P is an N-by-N matrix describing
! the column pivoting operations.
end program
Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR fact...
Definition linalg.f90:1281
Performs sparse matrix multiplication C = A * B.
Definition linalg.f90:5382
Computes the QR factorization of an M-by-N matrix.
Definition linalg.f90:1121
Solves a triangular system of equations.
Definition linalg.f90:2316
Provides the transpose of a sparse matrix.
Definition linalg.f90:5468
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 1281 of file linalg.f90.


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