linalg  1.4.3
A linear algebra library that provides a user-friendly interface to several BLAS and LAPACK routines.
linalg_factor Module Reference

linalg_factor More...

Data Types

interface  form_lu
 Extracts the L and U matrices from the condensed [L\U] storage format used by the lu_factor. More...
 
interface  form_qr
 Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm. More...
 
interface  mult_qr
 Multiplies a general matrix by the orthogonal matrix Q from a QR factorization. More...
 
interface  mult_rz
 Multiplies a general matrix by the orthogonal matrix Z from an RZ factorization. More...
 
interface  qr_factor
 Computes the QR factorization of an M-by-N matrix. More...
 

Functions/Subroutines

subroutine, public lu_factor (a, ipvt, err)
 Computes the LU factorization of an M-by-N matrix. More...
 
subroutine form_lu_all (lu, ipvt, u, p, err)
 Extracts the L, U, and P matrices from the output of the lu_factor routine. More...
 
subroutine form_lu_only (lu, u, err)
 Extracts the L, and U matrices from the output of the lu_factor routine. More...
 
subroutine qr_factor_no_pivot (a, tau, work, olwork, err)
 Computes the QR factorization of an M-by-N matrix without pivoting. More...
 
subroutine qr_factor_pivot (a, tau, jpvt, work, olwork, err)
 Computes the QR factorization of an M-by-N matrix with column pivoting such that A * P = Q * R. More...
 
subroutine form_qr_no_pivot (r, tau, q, work, olwork, err)
 Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm. More...
 
subroutine form_qr_pivot (r, tau, pvt, q, p, work, olwork, err)
 Forms the full M-by-M orthogonal matrix Q from the elementary reflectors returned by the base QR factorization algorithm. More...
 
subroutine mult_qr_mtx (lside, trans, a, tau, c, work, olwork, err)
 Multiplies a general matrix by the orthogonal matrix Q from a QR factorization such that: C = op(Q) * C, or C = C * op(Q). More...
 
subroutine mult_qr_vec (trans, a, tau, c, work, olwork, err)
 Multiplies a vector by the orthogonal matrix Q from a QR factorization such that: C = op(Q) * C. More...
 
subroutine, public qr_rank1_update (q, r, u, v, work, err)
 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. More...
 
subroutine, public cholesky_factor (a, upper, err)
 Computes the Cholesky factorization of a symmetric, positive definite matrix. More...
 
subroutine, public cholesky_rank1_update (r, u, work, err)
 Computes the rank 1 update to a Cholesky factored matrix (upper triangular). More...
 
subroutine, public cholesky_rank1_downdate (r, u, work, err)
 Computes the rank 1 downdate to a Cholesky factored matrix (upper triangular). More...
 
subroutine, public rz_factor (a, tau, work, olwork, err)
 Factors an upper trapezoidal matrix by means of orthogonal transformations such that A = R * Z = (R 0) * Z. Z is an orthogonal matrix of dimension N-by-N, and R is an M-by-M upper triangular matrix. More...
 
subroutine mult_rz_mtx (lside, trans, l, a, tau, c, work, olwork, err)
 Multiplies a general matrix by the orthogonal matrix Z from an RZ factorization such that: C = op(Z) * C, or C = C * op(Z). More...
 
subroutine mult_rz_vec (trans, l, a, tau, c, work, olwork, err)
 Multiplies a vector by the orthogonal matrix Z from an RZ factorization such that: C = op(Z) * C. More...
 
subroutine, public svd (a, s, u, vt, work, olwork, err)
 Computes the singular value decomposition of a matrix A. The SVD is defined as: A = U * S * V**T, where U is an M-by-M orthogonal matrix, S is an M-by-N diagonal matrix, and V is an N-by-N orthogonal matrix. More...
 

Detailed Description

linalg_factor

Purpose
Provides a set of matrix factorization routines.

Function/Subroutine Documentation

◆ cholesky_factor()

subroutine, public linalg_factor::cholesky_factor ( real(real64), dimension(:,:), intent(inout)  a,
logical, intent(in), optional  upper,
class(errors), intent(inout), optional, target  err 
)

Computes the Cholesky factorization of a symmetric, positive definite matrix.

Parameters
[in,out]aOn input, the N-by-N matrix to factor. On output, the factored matrix is returned in either the upper or lower triangular portion of the matrix, dependent upon the value of upper.
[in]upperAn optional input that, if specified, provides control over whether the factorization is computed as A = U**T * U (set to true), or as A = L * L**T (set to false). The default value is true such that A = U**T * U.
[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 a is not square.
  • LA_MATRIX_FORMAT_ERROR: Occurs if a is not positive definite.
Usage
The following example illustrates the solution of a positive-definite system of equations via Cholesky factorization.
program example
use iso_fortran_env, only : real64, int32
use linalg_factor, only : cholesky_factor
implicit none
! Variables
real(real64) :: a(3, 3), b(3), bu(3)
integer(int32) :: i
! Build the 3-by-3 positive-definite matrix A.
! | 4 12 -16 |
! A = | 12 37 -43 |
! |-16 -43 98 |
a = reshape([4.0d0, 12.0d0, -16.0d0, 12.0d0, 37.0d0, -43.0d0, -16.0d0, &
-43.0d0, 98.0d0], [3, 3])
! Build the 3-element array B
! | 5 |
! b = | 1 |
! | 3 |
b = [5.0d0, 1.0d0, 3.0d0]
! Make a copy of B for later use - not necessary, but just for example to
! illustrate the long or manual method of solving a Cholesky factored system
bu = b
! Compute the Cholesky factorization of A considering only the upper
! triangular portion of A (the default configuration).
call cholesky_factor(a)
! Compute the solution
call solve_cholesky(.true., a, b)
! Display the results
print '(A)', "Cholesky Solution: X = "
print '(F8.4)', (b(i), i = 1, size(b))
! The solution could also be computed manually noting the Cholesky
! factorization causes A = U**T * U. Then U**T * U * X = B.
! Step 1 would then be to solve the problem U**T * Y = B, for Y.
call solve_triangular_system(.true., .true., .true., a, bu)
! Now, solve the problem U * X = Y, for X
call solve_triangular_system(.true., .false., .true., a, bu)
! Display the results
print '(A)', "Cholesky Solution (Manual Approach): X = "
print '(F8.4)', (bu(i), i = 1, size(bu))
end program
The above program produces the following output.
Cholesky Solution: X =
239.5833
-65.6667
10.3333
Cholesky Solution (Manual Approach): X =
239.5833
-65.6667
10.3333
Notes
This routine utilizes the LAPACK routine DPOTRF.

Definition at line 1632 of file linalg_factor.f90.

◆ cholesky_rank1_downdate()

subroutine, public linalg_factor::cholesky_rank1_downdate ( real(real64), dimension(:,:), intent(inout)  r,
real(real64), dimension(:), intent(inout)  u,
real(real64), dimension(:), intent(out), optional, target  work,
class(errors), intent(inout), optional, target  err 
)

Computes the rank 1 downdate to a Cholesky factored matrix (upper triangular).

Parameters
[in,out]rOn input, the N-by-N upper triangular matrix R. On output, the updated matrix R1.
[in,out]uOn input, the N-element update vector U. On output, the rotation sines used to transform R to R1.
[out]workAn optional argument that if supplied prevents local memory allocation. If provided, the array must have at least N elements. Additionally, this workspace array is used to contain the rotation cosines used to transform R to R1.
[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 array sizes are incorrect.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
  • LA_MATRIX_FORMAT_ERROR: Occurs if the downdated matrix is not positive definite.
  • LA_SINGULAR_MATRIX_ERROR: Occurs if r is singular.
Usage
The following example illustrates the use of the rank 1 Cholesky downdate, and compares the results to factoring the original rank 1 downdated matrix.
program example
use iso_fortran_env, only : real64, int32
use linalg_factor, only : cholesky_factor, cholesky_rank1_downdate
use linalg_core, only : rank1_update
implicit none
! Variables
real(real64) :: a(3,3), u(3), ad(3,3)
integer(int32) :: i
! Build the 3-by-3 matrix A.
! | 4.25 11.25 -15 |
! A = | 11.25 39.25 -46 |
! | -15 -46 102 |
a = reshape([4.25d0, 11.25d0, -15.0d0, 11.25d0, 39.25d0, -46.0d0, &
-15.0d0, -46.0d0, 102.0d0], [3, 3])
! The downdate vector
! | 0.5 |
! u = | -1.5 |
! | 2 |
u = [0.5d0, -1.5d0, 2.0d0]
! Compute the rank 1 downdate of A
ad = a
call rank1_update(-1.0d0, u, u, ad)
! Compute the Cholesky factorization of the original matrix
call cholesky_factor(a)
! Apply the rank 1 downdate to the factored matrix
call cholesky_rank1_downdate(a, u)
! Compute the Cholesky factorization of the downdate to the original matrix
call cholesky_factor(ad)
! Display the matrices
print '(A)', "Downdating the Factored Form:"
do i = 1, size(a, 1)
print *, a(i,:)
end do
print '(A)', "Downdating A Directly:"
do i = 1, size(ad, 1)
print *, ad(i,:)
end do
end program
The above program produces the following output.
Downdating the Factored Form:
2.0000000000000000 6.0000000000000000 -8.0000000000000000
0.0000000000000000 1.0000000000000000 4.9999999999999973
0.0000000000000000 0.0000000000000000 3.0000000000000049
Downdating A Directly:
2.0000000000000000 6.0000000000000000 -8.0000000000000000
0.0000000000000000 1.0000000000000000 5.0000000000000000
0.0000000000000000 0.0000000000000000 3.0000000000000000
Notes
This routine utilizes the QRUPDATE routine DCH1DN.
See Also
Source

Definition at line 1944 of file linalg_factor.f90.

◆ cholesky_rank1_update()

subroutine, public linalg_factor::cholesky_rank1_update ( real(real64), dimension(:,:), intent(inout)  r,
real(real64), dimension(:), intent(inout)  u,
real(real64), dimension(:), intent(out), optional, target  work,
class(errors), intent(inout), optional, target  err 
)

Computes the rank 1 update to a Cholesky factored matrix (upper triangular).

Parameters
[in,out]rOn input, the N-by-N upper triangular matrix R. On output, the updated matrix R1.
[in,out]uOn input, the N-element update vector U. On output, the rotation sines used to transform R to R1.
[out]workAn optional argument that if supplied prevents local memory allocation. If provided, the array must have at least N elements. Additionally, this workspace array is used to contain the rotation cosines used to transform R to R1.
[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 array sizes are incorrect.
  • LA_OUT_OF_MEMORY_ERROR: Occurs if local memory must be allocated, and there is insufficient memory available.
Usage
The following example illustrates the use of the rank 1 Cholesky update, and compares the results to factoring the original rank 1 updated matrix.
program example
use iso_fortran_env, only : real64, int32
use linalg_factor, only : cholesky_factor, cholesky_rank1_update
use linalg_core, only : rank1_update
implicit none
! Variables
real(real64) :: a(3,3), u(3), au(3,3)
integer(int32) :: i
! Build the 3-by-3 positive-definite matrix A.
! | 4 12 -16 |
! A = | 12 37 -43 |
! |-16 -43 98 |
a = reshape([4.0d0, 12.0d0, -16.0d0, 12.0d0, 37.0d0, -43.0d0, -16.0d0, &
-43.0d0, 98.0d0], [3, 3])
! Build the update vector U
u = [0.5d0, -1.5d0, 2.0d0]
! Compute the rank 1 update of A
au = a
call rank1_update(1.0d0, u, u, au)
! Compute the Cholesky factorization of the original matrix
call cholesky_factor(a)
! Apply the rank 1 update to the factored matrix
call cholesky_rank1_update(a, u)
! Compute the Cholesky factorization of the update of the original matrix
call cholesky_factor(au)
! Display the matrices
print '(A)', "Updating the Factored Form:"
do i = 1, size(a, 1)
print *, a(i,:)
end do
print '(A)', "Updating A Directly:"
do i = 1, size(au, 1)
print *, au(i,:)
end do
end program
The above program produces the following output.
Updating the Factored Form:
2.0615528128088303 5.4570515633174921 -7.2760687510899889
0.0000000000000000 3.0774320845949008 -2.0452498947307731
0.0000000000000000 0.0000000000000000 6.6989384530323566
Updating A Directly:
2.0615528128088303 5.4570515633174921 -7.2760687510899889
0.0000000000000000 3.0774320845949008 -2.0452498947307736
0.0000000000000000 0.0000000000000000 6.6989384530323557
Notes
This routine utilizes the QRUPDATE routine DCH1UP.
See Also
Source

Definition at line 1785 of file linalg_factor.f90.

◆ form_lu_all()

subroutine linalg_factor::form_lu_all ( real(real64), dimension(:,:), intent(inout)  lu,
integer(int32), dimension(:), intent(in)  ipvt,
real(real64), dimension(:,:), intent(out)  u,
real(real64), dimension(:,:), intent(out)  p,
class(errors), intent(inout), optional, target  err 
)
private

Extracts the L, U, and P matrices from the output of the lu_factor routine.

Parameters
[in,out]luOn input, the N-by-N matrix as output by lu_factor. On output, the N-by-N lower triangular matrix L.
[in]ipvtThe N-element pivot array as output by lu_factor.
[out]uAn N-by-N matrix where the U matrix will be written.
[out]pAn N-by-N matrix where the row permutation matrix will be written.
[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 array sizes are incorrect.
Remarks
This routine allows extraction of the actual "L", "U", and "P" matrices of the decomposition. To use these matrices to solve the system A*X = B, the following approach is used.
  1. First, solve the linear system: L*Y = P*B for Y.
  2. Second, solve the linear system: U*X = Y for X.

Notice, as both L and U are triangular in structure, the above equations can be solved by forward and backward substitution.

See Also

Definition at line 498 of file linalg_factor.f90.

◆ form_lu_only()

subroutine linalg_factor::form_lu_only ( real(real64), dimension(:,:), intent(inout)  lu,
real(real64), dimension(:,:), intent(out)  u,
class(errors), intent(inout), optional, target  err 
)
private

Extracts the L, and U matrices from the output of the lu_factor routine.

Parameters
[in,out]luOn input, the N-by-N matrix as output by lu_factor. On output, the N-by-N lower triangular matrix L.
[out]uAn N-by-N matrix where the U matrix will be written.
[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 array sizes are incorrect.

Definition at line 575 of file linalg_factor.f90.

◆ form_qr_no_pivot()

subroutine linalg_factor::form_qr_no_pivot ( real(real64), dimension(:,:), intent(inout)  r,
real(real64), dimension(:), intent(in)  tau,
real(real64), dimension(:,:), intent(out)  q,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

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

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.
[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.

Definition at line 882 of file linalg_factor.f90.

◆ form_qr_pivot()

subroutine linalg_factor::form_qr_pivot ( real(real64), dimension(:,:), intent(inout)  r,
real(real64), dimension(:), intent(in)  tau,
integer(int32), dimension(:), intent(in)  pvt,
real(real64), dimension(:,:), intent(out)  q,
real(real64), dimension(:,:), intent(out)  p,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

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

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.
[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.

Definition at line 1013 of file linalg_factor.f90.

◆ lu_factor()

subroutine, public linalg_factor::lu_factor ( real(real64), dimension(:,:), intent(inout)  a,
integer(int32), dimension(:), intent(out)  ipvt,
class(errors), intent(inout), optional, target  err 
)

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

Parameters
[in,out]aOn input, the M-by-N matrix on which to operate. On output, the LU factored matrix in the form [L\U] where the unit diagonal elements of L are not stored.
[out]ipvtAn MIN(M, N)-element array used to track row-pivot operations. The array stored pivot information such that row I is interchanged with row IPVT(I).
[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 ipvt is not sized appropriately.
  • LA_SINGULAR_MATRIX_ERROR: Occurs as a warning if a is found to be singular.
Usage
To solve a system of 3 equations of 3 unknowns using LU factorization, the following code will suffice.
program example
use iso_fortran_env
use linalg_factor, only : lu_factor
use linalg_solve, only : solve_lu
implicit none
! Local Variables
real(real64) :: a(3,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 LU factorization
call lu_factor(a, pvt)
! Compute the solution. The results overwrite b.
call solve_lu(a, pvt, b)
! Display the results.
print '(A)', "LU Solution: X = "
print '(F8.4)', (b(i), i = 1, size(b))
end program
The program generates the following output.
LU Solution: X =
0.3333
-0.6667
0.0000
Notes
This routine utilizes the LAPACK routine DGETRF.
See Also

Definition at line 414 of file linalg_factor.f90.

◆ mult_qr_mtx()

subroutine linalg_factor::mult_qr_mtx ( logical, intent(in)  lside,
logical, intent(in)  trans,
real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(in)  tau,
real(real64), dimension(:,:), intent(inout)  c,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

Multiplies a general matrix by the orthogonal matrix Q from a QR factorization such that: C = op(Q) * C, or C = C * op(Q).

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]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.
[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]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.

Definition at line 1115 of file linalg_factor.f90.

◆ mult_qr_vec()

subroutine linalg_factor::mult_qr_vec ( logical, intent(in)  trans,
real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(in)  tau,
real(real64), dimension(:), intent(inout)  c,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

Multiplies a vector by the orthogonal matrix Q from a QR factorization such that: C = op(Q) * C.

Parameters
[in]transSet to true to apply Q**T; else, set to false.
[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.
[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 is based upon the LAPACK routine DORM2R.

Definition at line 1242 of file linalg_factor.f90.

◆ mult_rz_mtx()

subroutine linalg_factor::mult_rz_mtx ( logical, intent(in)  lside,
logical, intent(in)  trans,
integer(int32), intent(in)  l,
real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(in)  tau,
real(real64), dimension(:,:), intent(inout)  c,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

Multiplies a general matrix by the orthogonal matrix Z from an RZ factorization such that: C = op(Z) * C, or C = C * op(Z).

Parameters
[in]lsideSet to true to apply Z or Z**T from the left; else, set to false to apply Z or Z**T from the right.
[in]transSet to true to apply Z**T; else, set to false.
[in]lThe number of columns in matrix a containing the meaningful part of the Householder vectors. If lside is true, M >= L >= 0; else, if lside is false, N >= L >= 0.
[in,out]aOn input the K-by-LTA matrix Z, where LTA = M if lside is true; else, LTA = N if lside is false. The I-th row must contain the Householder vector in the last k rows. Notice, the contents of this matrix are restored on exit.
[in]tauA K-element array containing the scalar factors of the elementary reflectors, where M >= K >= 0 if lside is true; else, N >= K >= 0 if lside is false.
[in,out]cOn input, the M-by-N matrix C. On output, the product of the orthogonal matrix Z 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.
[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]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 DORMRZ.

Definition at line 2195 of file linalg_factor.f90.

◆ mult_rz_vec()

subroutine linalg_factor::mult_rz_vec ( logical, intent(in)  trans,
integer(int32), intent(in)  l,
real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(in)  tau,
real(real64), dimension(:), intent(inout)  c,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

Multiplies a vector by the orthogonal matrix Z from an RZ factorization such that: C = op(Z) * C.

Parameters
[in]transSet to true to apply Z**T; else, set to false.
[in]lThe number of columns in matrix a containing the meaningful part of the Householder vectors. If lside is true, M >= L >= 0; else, if lside is false, N >= L >= 0.
[in,out]aOn input the K-by-LTA matrix Z, where LTA = M if lside is true; else, LTA = N if lside is false. The I-th row must contain the Householder vector in the last k rows. Notice, the contents of this matrix are restored on exit.
[in]tauA K-element array containing the scalar factors of the elementary reflectors, where M >= K >= 0 if lside is true; else, N >= K >= 0 if lside is false.
[in,out]cOn input, the M-element array C. On output, the product of the orthogonal matrix Z and the original array 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.
[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 DORMRZ.

Definition at line 2337 of file linalg_factor.f90.

◆ qr_factor_no_pivot()

subroutine linalg_factor::qr_factor_no_pivot ( real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(out)  tau,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

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

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.
[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.

Definition at line 664 of file linalg_factor.f90.

◆ qr_factor_pivot()

subroutine linalg_factor::qr_factor_pivot ( real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(out)  tau,
integer(int32), dimension(:), intent(inout)  jpvt,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)
private

Computes the QR factorization of an M-by-N matrix with column pivoting such that A * P = Q * R.

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]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.

Definition at line 767 of file linalg_factor.f90.

◆ qr_rank1_update()

subroutine, public linalg_factor::qr_rank1_update ( real(real64), dimension(:,:), intent(inout)  q,
real(real64), dimension(:,:), intent(inout)  r,
real(real64), dimension(:), intent(inout)  u,
real(real64), dimension(:), intent(inout)  v,
real(real64), dimension(:), intent(out), optional, target  work,
class(errors), intent(inout), optional, target  err 
)

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.

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 2*K elements.
[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 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.
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
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
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
Notes
This routine utilizes the QRUPDATE routine DQR1UP.
See Also
Source

Definition at line 1462 of file linalg_factor.f90.

◆ rz_factor()

subroutine, public linalg_factor::rz_factor ( real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(out)  tau,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)

Factors an upper trapezoidal matrix by means of orthogonal transformations such that A = R * Z = (R 0) * Z. Z is an orthogonal matrix of dimension N-by-N, and R is an M-by-M upper triangular matrix.

Parameters
[in,out]aOn input, the M-by-N upper trapezoidal matrix to factor. On output, the leading M-by-M upper triangular part of the matrix contains the upper triangular matrix R, and elements N-L+1 to N of the first M rows of A, with the array tau, represent the orthogonal matrix Z as a product of M elementary reflectors.
[out]tauAn M-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.
[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.
Further Details
  The factorization is obtained by Householder's method.  The kth
  transformation matrix, Z( k ), which is used to introduce zeros into
  the ( m - k + 1 )th row of A, is given in the form

     Z( k ) = ( I     0   ),
              ( 0  T( k ) )

  where

     T( k ) = I - tau*u( k )*u( k )**T,   u( k ) = (   1    ),
                                                   (   0    )
                                                   ( z( k ) )

  tau is a scalar and z( k ) is an l element vector. tau and z( k )
  are chosen to annihilate the elements of the kth row of A2.

  The scalar tau is returned in the kth element of TAU and the vector
  u( k ) in the kth row of A2, such that the elements of z( k ) are
  in  a( k, l + 1 ), ..., a( k, n ). The elements of R are returned in
  the upper triangular part of A1.

  Z is given by

     Z =  Z( 1 ) * Z( 2 ) * ... * Z( m ).
Notes
This routine is based upon the LAPACK routine DTZRZF.
See Also

Definition at line 2083 of file linalg_factor.f90.

◆ svd()

subroutine, public linalg_factor::svd ( real(real64), dimension(:,:), intent(inout)  a,
real(real64), dimension(:), intent(out)  s,
real(real64), dimension(:,:), intent(out), optional  u,
real(real64), dimension(:,:), intent(out), optional  vt,
real(real64), dimension(:), intent(out), optional, target  work,
integer(int32), intent(out), optional  olwork,
class(errors), intent(inout), optional, target  err 
)

Computes the singular value decomposition of a matrix A. The SVD is defined as: A = U * S * V**T, where U is an M-by-M orthogonal matrix, S is an M-by-N diagonal matrix, and V is an N-by-N orthogonal matrix.

Parameters
[in,out]aOn input, the M-by-N matrix to factor. The matrix is overwritten on output.
[out]sA MIN(M, N)-element array containing the singular values of a sorted in descending order.
[out]uAn optional argument, that if supplied, is used to contain the orthogonal matrix U from the decomposition. The matrix U contains the left singular vectors, and can be either M-by-M (all left singular vectors are computed), or M-by-MIN(M,N) (only the first MIN(M, N) left singular vectors are computed).
[out]vtAn optional argument, that if supplied, is used to contain the transpose of the N-by-N orthogonal matrix V. The matrix V contains the right singular vectors.
[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]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.
  • LA_CONVERGENCE_ERROR: Occurs as a warning if the QR iteration process could not converge to a zero value.
Usage
! Decompose matrix the M-by-N matrix A such that A = U * S * V**T with
! M >= N.
! Variables
real(real64), dimension(m, n) :: a
real(real64), dimension(m, m) :: u
real(real64), dimension(n, n) :: vt
real(real64), dimension(n) :: s
! Initialize A...
! Compute the SVD of A. On output, S contains the MIN(M,N) singular
! values of A in descending order, U contains the left singular vectors
! (one per column), and VT contains the right singular vectors (one per
! row).
call svd(a, s, u, vt)
! Note: If M > N, then we can make U M-by-N, and compute the N
! left singular vectors of A, as there are at most N singular values
! of A. Also, if M < N, then there are at most M singular values of A,
! and as such, the length of the array s should be m.
Notes
This routine utilizes the LAPACK routine DGESVD.
See Also

Definition at line 2496 of file linalg_factor.f90.