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

Solves the overdetermined or underdetermined system \( A X = B \) of M equations of N unknowns using a singular value decomposition of matrix A. More...

Public Member Functions

 solve_least_squares_mtx_svd
 
 solve_least_squares_vec_svd
 

Detailed Description

Solves the overdetermined or underdetermined system \( A X = B \) of M equations of N unknowns using a singular value decomposition of matrix A.

Syntax
subroutine solve_least_squares_svd(real(real64) a(:,:), real(real64) b(:,:), optional real(real64) s(:), optional integer(int32) arnk, optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine solve_least_squares_svd(complex(real64) a(:,:), complex(real64) b(:,:), optional real(real64) s(:), optional integer(int32) arnk, optional complex(real64) work(:), optional integer(int32) olwork, optional real(real64) rwork(:), optional class(errors) err)
subroutine solve_least_squares_svd(real(real64) a(:,:), real(real64) b(:), optional real(real64) s(:), optional integer(int32) arnk, optional real(real64) work(:), optional integer(int32) olwork, optional class(errors) err)
subroutine solve_least_squares_svd(complex(real64) a(:,:), complex(real64) b(:), optional real(real64) s(:), optional integer(int32) arnk, optional complex(real64) work(:), optional integer(int32) olwork, optional real(real64) rwork(:), optional class(errors) err)
Parameters
[in,out]aOn input, the M-by-N matrix A. On output, the matrix is overwritten by the details of its complete orthogonal factorization.
[in,out]bIf M >= N, the M-by-NRHS matrix B. On output, the first N rows contain the N-by-NRHS solution matrix X. If M < N, an N-by-NRHS matrix with the first M rows containing the matrix B. On output, the N-by-NRHS solution matrix X.
[out]arnkAn optional output, that if provided, will return the rank of a.
[out]sAn optional MIN(M, N)-element array that on output contains the singular values of a in descending order. Notice, the condition number of a can be determined by S(1) / S(MIN(M, N)).
[out]arnkAn optional output, that if provided, will return the rank of a.
[out]workAn optional input, that if provided, prevents any local memory allocation for complex-valued workspaces. 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]rworkAn optional input, that if provided, prevents any local memory allocation for real-valued workspaces. If not provided, the memory required is allocated within. If provided, the length of the array must be at least 5 * MIN(M, N).
[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.
  • LA_CONVERGENCE_ERROR: Occurs as a warning if the QR iteration process could not converge to a zero value.
Notes
This routine utilizes the LAPACK routine DGELSS (ZGELSS in the complex case).
Usage
The following example illustrates the least squares solution of an overdetermined system of linear equations.
program example
use iso_fortran_env, only : real64, int32
use linalg
implicit none
! Local Variables
real(real64) :: a(3,2), b(3)
integer(int32) :: i
! Build the 3-by-2 matrix A
! | 2 1 |
! A = |-3 1 |
! |-1 1 |
a = reshape([2.0d0, -3.0d0, -1.0d0, 1.0d0, 1.0d0, 1.0d0], [3, 2])
! Build the right-hand-side vector B.
! |-1 |
! b = |-2 |
! | 1 |
b = [-1.0d0, -2.0d0, 1.0d0]
! The solution is:
! x = [0.13158, -0.57895]**T
! Compute the solution via a least-squares approach. The results overwrite
! the first 2 elements in b.
! Display the results
print '(A)', "Least Squares Solution: X = "
print '(F9.5)', (b(i), i = 1, size(a, 2))
end program
Solves the overdetermined or underdetermined system of M equations of N unknowns using a singular va...
Definition linalg.f90:2956
Provides a set of common linear algebra routines.
Definition linalg.f90:145
The above program produces the following output.
Least Squares Solution: X =
0.13158
-0.57895

Definition at line 2956 of file linalg.f90.


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