Processing math: 100%
fortran-lapack
All Classes Namespaces Files Functions Variables Pages
la_svd::svd Interface Reference

Compute the Singular Value Decomposition (SVD) of a matrix. More...

Public Member Functions

subroutine la_svd_s (a, s, u, vt, overwrite_a, full_matrices, err)
 SVD of matrix A = U S V^T, returning S and optionally U and V^T.
 
subroutine la_svd_d (a, s, u, vt, overwrite_a, full_matrices, err)
 SVD of matrix A = U S V^T, returning S and optionally U and V^T.
 
subroutine la_svd_q (a, s, u, vt, overwrite_a, full_matrices, err)
 SVD of matrix A = U S V^T, returning S and optionally U and V^T.
 
subroutine la_svd_c (a, s, u, vt, overwrite_a, full_matrices, err)
 SVD of matrix A = U S V^T, returning S and optionally U and V^T.
 
subroutine la_svd_z (a, s, u, vt, overwrite_a, full_matrices, err)
 SVD of matrix A = U S V^T, returning S and optionally U and V^T.
 
subroutine la_svd_w (a, s, u, vt, overwrite_a, full_matrices, err)
 SVD of matrix A = U S V^T, returning S and optionally U and V^T.
 

Detailed Description

Compute the Singular Value Decomposition (SVD) of a matrix.

This subroutine computes the Singular Value Decomposition (SVD) of a matrix A :

A = U \cdot S \cdot V^T

where A is a matrix of size [m,n] , U is an orthogonal matrix of size [m,m] , S is a diagonal matrix containing the singular values, and V^T is an orthogonal matrix of size [n,n] . The subroutine computes the singular values and optionally the matrices U and V^T .

Parameters
[in,out]aThe input matrix A of size [m,n] . If overwrite_a is true, the contents of a may be modified during computation.
[out]sThe array of singular values of size k = min(m,n) .
[out]u(Optional) The left singular vectors of matrix A , with shape [m,m] for the full problem or [m,k] for the reduced problem.
[out]vt(Optional) The right singular vectors of matrix A^T , with shape [k,n] for the reduced problem or [n,n] for the full problem.
[in]overwrite_a(Optional) A logical flag that determines whether matrix a may be overwritten. Default is false.
[in]full_matrices(Optional) If true, computes the full-sized matrices U and V^T . If false, computes the reduced matrices.
[out]err(Optional) A state return flag. If an error occurs and err is not provided, the function will stop execution.
Note
This subroutine performs the SVD computation using the LAPACK gesdd backend and can be used for both full and reduced decompositions.
Warning
If overwrite_a is enabled, the original contents of a may be lost during computation.

Compute the singular values of a matrix.

This function returns the singular values of the input matrix A . The singular values are stored in a vector s , which is an array of size k = min(m,n) .

Parameters
[in]aThe input matrix A of size [m,n] .
[out]err(Optional) A state return flag. If an error occurs and err is not provided, the function will stop execution.
Returns
The singular values of matrix A are returned in the real array s , with the same kind as the input matrix.
Note
This function returns only the singular values and does not compute the full SVD.

Member Function/Subroutine Documentation

◆ la_svd_c()

subroutine la_svd::svd::la_svd_c ( complex(sp), dimension(:,:), intent(inout), target a,
real(sp), dimension(:), intent(out) s,
complex(sp), dimension(:,:), intent(out), optional, target u,
complex(sp), dimension(:,:), intent(out), optional, target vt,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional full_matrices,
type(la_state), intent(out), optional err )

SVD of matrix A = U S V^T, returning S and optionally U and V^T.

Parameters
[in,out]aInput matrix A[m,n]
[out]sArray of singular values
[out]uThe columns of U contain the eigenvectors of A A^T
[out]vtThe rows of V^T contain the eigenvectors of A^T A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]full_matrices[optional] full matrices have shape(u)==[m,m], shape(vh)==[n,n] (default); otherwise they are shape(u)==[m,k] and shape(vh)==[k,n] with k=min(m,n)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_svd_d()

subroutine la_svd::svd::la_svd_d ( real(dp), dimension(:,:), intent(inout), target a,
real(dp), dimension(:), intent(out) s,
real(dp), dimension(:,:), intent(out), optional, target u,
real(dp), dimension(:,:), intent(out), optional, target vt,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional full_matrices,
type(la_state), intent(out), optional err )

SVD of matrix A = U S V^T, returning S and optionally U and V^T.

Parameters
[in,out]aInput matrix A[m,n]
[out]sArray of singular values
[out]uThe columns of U contain the eigenvectors of A A^T
[out]vtThe rows of V^T contain the eigenvectors of A^T A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]full_matrices[optional] full matrices have shape(u)==[m,m], shape(vh)==[n,n] (default); otherwise they are shape(u)==[m,k] and shape(vh)==[k,n] with k=min(m,n)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_svd_q()

subroutine la_svd::svd::la_svd_q ( real(qp), dimension(:,:), intent(inout), target a,
real(qp), dimension(:), intent(out) s,
real(qp), dimension(:,:), intent(out), optional, target u,
real(qp), dimension(:,:), intent(out), optional, target vt,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional full_matrices,
type(la_state), intent(out), optional err )

SVD of matrix A = U S V^T, returning S and optionally U and V^T.

Parameters
[in,out]aInput matrix A[m,n]
[out]sArray of singular values
[out]uThe columns of U contain the eigenvectors of A A^T
[out]vtThe rows of V^T contain the eigenvectors of A^T A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]full_matrices[optional] full matrices have shape(u)==[m,m], shape(vh)==[n,n] (default); otherwise they are shape(u)==[m,k] and shape(vh)==[k,n] with k=min(m,n)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_svd_s()

subroutine la_svd::svd::la_svd_s ( real(sp), dimension(:,:), intent(inout), target a,
real(sp), dimension(:), intent(out) s,
real(sp), dimension(:,:), intent(out), optional, target u,
real(sp), dimension(:,:), intent(out), optional, target vt,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional full_matrices,
type(la_state), intent(out), optional err )

SVD of matrix A = U S V^T, returning S and optionally U and V^T.

Parameters
[in,out]aInput matrix A[m,n]
[out]sArray of singular values
[out]uThe columns of U contain the eigenvectors of A A^T
[out]vtThe rows of V^T contain the eigenvectors of A^T A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]full_matrices[optional] full matrices have shape(u)==[m,m], shape(vh)==[n,n] (default); otherwise they are shape(u)==[m,k] and shape(vh)==[k,n] with k=min(m,n)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_svd_w()

subroutine la_svd::svd::la_svd_w ( complex(qp), dimension(:,:), intent(inout), target a,
real(qp), dimension(:), intent(out) s,
complex(qp), dimension(:,:), intent(out), optional, target u,
complex(qp), dimension(:,:), intent(out), optional, target vt,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional full_matrices,
type(la_state), intent(out), optional err )

SVD of matrix A = U S V^T, returning S and optionally U and V^T.

Parameters
[in,out]aInput matrix A[m,n]
[out]sArray of singular values
[out]uThe columns of U contain the eigenvectors of A A^T
[out]vtThe rows of V^T contain the eigenvectors of A^T A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]full_matrices[optional] full matrices have shape(u)==[m,m], shape(vh)==[n,n] (default); otherwise they are shape(u)==[m,k] and shape(vh)==[k,n] with k=min(m,n)
[out]err[optional] state return flag. On error if not requested, the code will stop

◆ la_svd_z()

subroutine la_svd::svd::la_svd_z ( complex(dp), dimension(:,:), intent(inout), target a,
real(dp), dimension(:), intent(out) s,
complex(dp), dimension(:,:), intent(out), optional, target u,
complex(dp), dimension(:,:), intent(out), optional, target vt,
logical(lk), intent(in), optional overwrite_a,
logical(lk), intent(in), optional full_matrices,
type(la_state), intent(out), optional err )

SVD of matrix A = U S V^T, returning S and optionally U and V^T.

Parameters
[in,out]aInput matrix A[m,n]
[out]sArray of singular values
[out]uThe columns of U contain the eigenvectors of A A^T
[out]vtThe rows of V^T contain the eigenvectors of A^T A
[in]overwrite_a[optional] Can A data be overwritten and destroyed?
[in]full_matrices[optional] full matrices have shape(u)==[m,m], shape(vh)==[n,n] (default); otherwise they are shape(u)==[m,k] and shape(vh)==[k,n] with k=min(m,n)
[out]err[optional] state return flag. On error if not requested, the code will stop

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