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

Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where \( A = Q R \), and \( A1 = A + U V^T \) such that \( A1 = Q1 R1 \). In the event \( V \) is complex-valued, \( V^H \) is computed instead of \( V^T \). More...

Public Member Functions

 qr_rank1_update_dbl
 
 qr_rank1_update_cmplx
 

Detailed Description

Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where \( A = Q R \), and \( A1 = A + U V^T \) such that \( A1 = Q1 R1 \). In the event \( V \) is complex-valued, \( V^H \) is computed instead of \( V^T \).

Syntax
subroutine qr_rank1_update(real(real64) q(:,:), real(real64) r(:,:), real(real64) u(:), real(real64) v(:), optional real(real64) work(:), optional class(errors) err)
subroutine qr_rank1_update(complex(real64) q(:,:), complex(real64) r(:,:), complex(real64) u(:), complex(real64) v(:), optional complex(real64) work(:), optional real(real64) rwork(:), optional class(errors) err)
Parameters
[in,out]qOn input, the original M-by-K orthogonal matrix Q. On output, the updated matrix Q1.
[in,out]rOn input, the M-by-N matrix R. On output, the updated matrix R1.
[in,out]uOn input, the M-element U update vector. On output, the original content of the array is overwritten.
[in,out]vOn input, the N-element V update vector. On output, the original content of the array is overwritten.
[out]workAn optional argument that if supplied prevents local memory allocation. If provided, the array must have at least K elements.
[out]rworkAn optional argument that if supplied prevents local memory allocation. If provided, the array must have at least K elements.
[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.
Remarks
Notice, K must either be equal to M, or equal to N. In the event that K = N, only the submatrix Qa is updated. This is appropriate as the QR factorization for an overdetermined system can be written as follows:
  A = Q * R = [Qa, Qb] * [Ra]
                         [0 ]
Note: Ra is upper triangular of dimension N-by-N.
Notes
This routine utilizes the QRUPDATE routine ZQR1UP.
See Also
Source
Usage
The following example illustrates a rank 1 update to a QR factored system. The results are compared to updating the original matrix, and then performing the factorization.
program example
use iso_fortran_env
use linalg
implicit none
! Variables
real(real64) :: a(3,3), u(3), v(3), r(3,3), tau(3), q(3,3), qu(3,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 update vectors
! | 1/2 | | 1 |
! u = | 3/2 |, v = | 5 |
! | 3 | | 2 |
u = [0.5d0, 1.5d0, 3.0d0]
v = [1.0d0, 5.0d0, 2.0d0]
! Compute the QR factorization of the original matrix
r = a ! Making a copy as the matrix will be overwritten by qr_factor
call qr_factor(r, tau)
! Form Q & R
call form_qr(r, tau, q)
! Compute the rank 1 update to the original matrix such that:
! A = A + u * v**T
call rank1_update(1.0d0, u, v, a)
! Compute the rank 1 update to the factorization. Notice, the contents
! of U & V are destroyed as part of this process.
call qr_rank1_update(q, r, u, v)
! As comparison, compute the QR factorization of the rank 1 updated matrix
call qr_factor(a, tau)
call form_qr(a, tau, qu)
! Display the matrices
print '(A)', "Updating the Factored Form:"
print '(A)', "Q = "
do i = 1, size(q, 1)
print *, q(i,:)
end do
print '(A)', "R = "
do i = 1, size(r, 1)
print *, r(i,:)
end do
print '(A)', "Updating A Directly:"
print '(A)', "Q = "
do i = 1, size(qu, 1)
print *, qu(i,:)
end do
print '(A)', "R = "
do i = 1, size(a, 1)
print *, a(i,:)
end do
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
Computes the QR factorization of an M-by-N matrix.
Definition linalg.f90:1121
Computes the rank 1 update to an M-by-N QR factored matrix A (M >= N) where , and such that ....
Definition linalg.f90:1588
Performs the rank-1 update to matrix A such that: , where is an M-by-N matrix, is a scalar,...
Definition linalg.f90:396
Provides a set of common linear algebra routines.
Definition linalg.f90:145
The above program produces the following output.
Updating the Factored Form:
Q =
-0.13031167282892092 0.98380249683206911 -0.12309149097933236
-0.47780946703937632 -0.17109608640557677 -0.86164043685532932
-0.86874448552613881 -5.3467527001743037E-002 0.49236596391733078
R =
-11.510864433221338 -26.540144032823541 -10.033998807826904
0.0000000000000000 1.0586570346345126 2.0745400476676279
0.0000000000000000 0.0000000000000000 -5.2929341121113067
Updating A Directly:
Q =
-0.13031167282892087 0.98380249683206955 -0.12309149097933178
-0.47780946703937643 -0.17109608640557616 -0.86164043685532943
-0.86874448552613903 -5.3467527001742954E-002 0.49236596391733084
R =
-11.510864433221336 -26.540144032823545 -10.033998807826906
0.0000000000000000 1.0586570346345205 2.0745400476676350
0.0000000000000000 0.0000000000000000 -5.2929341121113058

Definition at line 1588 of file linalg.f90.


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