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_lq Interface Reference

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

Public Member Functions

 mult_lq_mtx
 
 mult_lq_mtx_cmplx
 
 mult_lq_vec
 
 mult_lq_vec_cmplx
 

Detailed Description

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

Syntax 1
Multiplies a general matrix by the orthogonal matrix \( Q \) from a LQ 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 K-by-P matrix containing the elementary reflectors output from the LQ factorization. If lside is set to true, P = M; else, if lside is set to false, P = N.
[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.
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 K-by-M matrix containing the elementary reflectors output from the LQ 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 DORMLQ (ZUNMLQ in the complex case).
Usage
The folowing example illustrates the solution of a system of equations using LQ 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 LQ factorization
call lq_factor(a, tau)
! Solve the lower triangular problem and store the solution in B.
!
! A comment about this solution noting we've factored A = L * Q.
!
! We then have to solve: L * Q * X = B for X. If we let Y = Q * X, then
! we solve the lower triangular system L * Y = B for Y.
call solve_triangular_system(.false., .false., .true., a, b)
! Now we've solved the lower triangular system L * Y = B for Y. At
! this point we solve the problem: Q * X = Y. Q is an orthogonal matrix;
! therefore, inv(Q) = Q**T. We can solve this by multiplying both
! sides by Q**T:
!
! Compute Q**T * B = X
call mult_lq(.true., a, tau, b)
! Display the results
print '(A)', "LQ Solution: X = "
print '(F8.4)', (b(i), i = 1, size(b))
end program
Computes the LQ factorization of an M-by-N matrix.
Definition linalg.f90:3570
Multiplies a general matrix by the orthogonal matrix Q from a LQ factorization.
Definition linalg.f90:3833
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.
LQ Solution: X =
0.3333
-0.6667
0.0000
See Also

Definition at line 3833 of file linalg.f90.


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