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

Multiplies a general matrix by the orthogonal matrix Q from a QR factorization. More...

Public Member Functions

 mult_qr_mtx
 
 mult_qr_mtx_cmplx
 
 mult_qr_vec
 
 mult_qr_vec_cmplx
 

Detailed Description

Multiplies a general matrix by the orthogonal matrix Q from a QR factorization.

Syntax 1
Multiplies a general matrix by the orthogonal matrix \( Q \) from a QR factorization such that: \( C = op(Q) C \), or \( C = C op(Q) \).
subroutine mult_qr(logical lside, logical trans, real(real64) a(:,:), real(real64) tau(:), real(real64) c(:,:), optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine mult_qr(logical lside, logical trans, complex(real64) a(:,:), complex(real64) tau(:), complex(real64) c(:,:), optional complex(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
Parameters
[in]lsideSet to true to apply \( Q \) or \( Q^T \) from the left; else, set to false to apply \( Q \) or \( Q^T \) from the right.
[in]transSet to true to apply \( Q^T \); else, set to false. In the event \( Q \) is complex-valued, \( Q^H \) is computed instead of \( Q^T \).
[in]aOn input, an LDA-by-K matrix containing the elementary reflectors output from the QR factorization. If lside is set to true, LDA = M, and M >= K >= 0; else, if lside is set to false, LDA = N, and N >= K >= 0. Notice, the contents of this matrix are restored on exit.
[in]tauA K-element array containing the scalar factors of each elementary reflector defined in a.
[in,out]cOn input, the M-by-N matrix C. On output, the product of the orthogonal matrix Q and the original matrix C.
[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.
[in,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 DORMQR (ZUNMQR in the complex case).
Syntax 2
Multiplies a vector by the orthogonal matrix \( Q \) from a QR factorization such that: \( C = op(Q) C\).
subroutine mult_qr(logical trans, real(real64) a(:,:), real(real64) tau(:), real(real64) c(:), optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine mult_qr(logical trans, complex(real64) a(:,:), complex(real64) tau(:), complex(real64) c(:), optional complex(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
Parameters
[in]transSet to true to apply \( Q^T \); else, set to false. In the event \( Q \) is complex-valued, \( Q^H \) is computed instead of \( Q^T \).
[in]aOn input, an M-by-K matrix containing the elementary reflectors output from the QR factorization. Notice, the contents of this matrix are restored on exit.
[in]tauA K-element array containing the scalar factors of each elementary reflector defined in a.
[in,out]cOn input, the M-element vector C. On output, the product of the orthogonal matrix Q and the original vector C.
[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 DORMQR (ZUNMQR in the complex case).
Usage
The following example illustrates how to perform the multiplication \( Q^T B \) when solving a system of QR factored equations without explicitly forming the matrix \( Q \).
program example
use iso_fortran_env, only : real64, int32
use linalg
implicit none
! Variables
real(real64) :: a(3,3), b(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)
! 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. Notice, using mult_qr
! avoids direct construction of the full Q and R matrices.
call mult_qr(.true., a, tau, 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
Multiplies a general matrix by the orthogonal matrix Q from a QR factorization.
Definition linalg.f90:1438
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 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

Definition at line 1438 of file linalg.f90.


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