SVD

   more...
Computes the singular value decomposition of a matrix, A = USVT.
Function Return Value
× n diagonal matrix of singular values, S, from the A = USVT decomposition. (Output)
Required Argument
A — Array of size × n to be decomposed. Must be rank-2 or rank-3 array of type single, double, complex, or double complex. (Input)
Optional Arguments, Packaged Options
U — Array of size × m containing the singular vectors U. (Output)
V — Array of size × n containing the singular vectors V. (Output)
The option and derived type names are given in the following tables:
Option Names for SVD
Option Value
Options_for_lin_svd
1
Options_for_lin_sol_svd
2
skip_error_processing
5
Name of Unallocated Option Array to Use for Setting Options
Use
Derived Type
?_svd_options(:)
Use when setting options for calls hereafter.
?_options
?_svd_options_once(:)
Use when setting options for next call only.
?_options
For a description on how to use these options, see Matrix Optional Data Changes. See LIN_SVD and LIN_SOL_SVD in Chapter 1, “Linear Systems” for the specific options for these routines.
FORTRAN 90 Interface
SVD (A [])
Description
Computes the singular value decomposition of a rank-2 or rank-3 array, A = USVT.
This function uses one of the routines LIN_SVD and LIN_SOL_SVD. If a complete decomposition is required, LIN_SVD is used. If singular values only, or singular values and one of the right and left singular vectors are required, then LIN_SOL_SVD is called.
Examples
Example 1: operator_ex14.f90
 
use linear_operators
implicit none
! This is the equivalent of Example 2 for LIN_SOL_SVD using operators
! and functions.
integer, parameter :: n=32
real(kind(1d0)) :: one=1d0, zero=0d0
real(kind(1d0)) A(n,n), P(n,n), Q(n,n), &
S_D(n), U_D(n,n), V_D(n,n)
! Generate a random matrix.
A = rand(A)
! Compute the singular value decomposition.
S_D = SVD(A, U=U_D, V=V_D)
! Compute the (left) orthogonal factor.
P = U_D .xt. V_D
! Compute the (right) self-adjoint factor.
Q = V_D .x. diag(S_D) .xt. V_D
! Check the results.
if (norm( EYE(n) - (P .xt. P)) &
<= sqrt(epsilon(one))) then
if (norm(A - (P .x. Q))/norm(A) &
<= sqrt(epsilon(one))) then
write (*,*) 'Example 2 for LIN_SOL_SVD (operators) is correct.'
end if
end if
end
Parallel Example (parallel_ex14.f90)
Systems of least-squares problems are solved, but now using the SVD() function. A box data type is used. This is an example which uses optional arguments and a generic function overloaded for parallel execution of a box data type. Any number of nodes can be used.
 
use linear_operators
use mpi_setup_int
implicit none
 
! This is the equivalent of Parallel Example 14
! for SVD, .tx. , .x. and NORM.
integer, parameter :: m=128, n=32, nr=4
real(kind(1d0)) :: one=1d0, err(nr)
real(kind(1d0)) A(m,n,nr), b(m,1,nr), x(n,1,nr), U(m,m,nr), &
V(n,n,nr), S(n,nr), g(m,1,nr)
 
! Setup for MPI:
mp_nprocs=mp_setup()
 
if(mp_rank == 0) then
! Generate a random matrix and right-hand side.
A = rand(A); b = rand(b)
endif
 
! Compute the least-squares solution matrix of Ax=b.
S = SVD(A, U = U, V = V)
g = U .tx. b
x = V .x. (diag(one/S) .x. g(1:n,:,:))
 
! Check the results.
err = norm(A .tx. (b - (A .x. x)))/(norm(A)+norm(x))
if (ALL(err <= sqrt(epsilon(one)))) then
if(mp_rank == 0) &
write (*,*) 'Parallel Example 14 is correct.'
end if
 
! See to any error messages and quit MPI
mp_nprocs = mp_setup('Final')
 
end
Published date: 03/19/2020
Last modified date: 03/19/2020