fortran-lapack
Loading...
Searching...
No Matches
la_least_squares::lstsq Interface Reference

Compute a least squares solution to system Ax=b, i.e. such that the 2-norm abs(b-Ax) is minimized. More...

Public Member Functions

real(sp) function, dimension(:), allocatable, target la_slstsq_one (a, b, cond, overwrite_a, rank, err)
 
real(dp) function, dimension(:), allocatable, target la_dlstsq_one (a, b, cond, overwrite_a, rank, err)
 
real(qp) function, dimension(:), allocatable, target la_qlstsq_one (a, b, cond, overwrite_a, rank, err)
 
complex(sp) function, dimension(:), allocatable, target la_clstsq_one (a, b, cond, overwrite_a, rank, err)
 
complex(dp) function, dimension(:), allocatable, target la_zlstsq_one (a, b, cond, overwrite_a, rank, err)
 
complex(qp) function, dimension(:), allocatable, target la_wlstsq_one (a, b, cond, overwrite_a, rank, err)
 
real(sp) function, dimension(:,:), allocatable, target la_slstsq_multiple (a, b, cond, overwrite_a, rank, err)
 
real(dp) function, dimension(:,:), allocatable, target la_dlstsq_multiple (a, b, cond, overwrite_a, rank, err)
 
real(qp) function, dimension(:,:), allocatable, target la_qlstsq_multiple (a, b, cond, overwrite_a, rank, err)
 
complex(sp) function, dimension(:,:), allocatable, target la_clstsq_multiple (a, b, cond, overwrite_a, rank, err)
 
complex(dp) function, dimension(:,:), allocatable, target la_zlstsq_multiple (a, b, cond, overwrite_a, rank, err)
 
complex(qp) function, dimension(:,:), allocatable, target la_wlstsq_multiple (a, b, cond, overwrite_a, rank, err)
 

Detailed Description

Compute a least squares solution to system Ax=b, i.e. such that the 2-norm abs(b-Ax) is minimized.

Member Function/Subroutine Documentation

◆ la_clstsq_multiple()

complex(sp) function, dimension(:,:), allocatable, target la_least_squares::lstsq::la_clstsq_multiple ( complex(sp), dimension(:,:), intent(inout), target  a,
complex(sp), dimension(:,:), intent(in)  b,
real(sp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_clstsq_one()

complex(sp) function, dimension(:), allocatable, target la_least_squares::lstsq::la_clstsq_one ( complex(sp), dimension(:,:), intent(inout), target  a,
complex(sp), dimension(:), intent(in)  b,
real(sp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_dlstsq_multiple()

real(dp) function, dimension(:,:), allocatable, target la_least_squares::lstsq::la_dlstsq_multiple ( real(dp), dimension(:,:), intent(inout), target  a,
real(dp), dimension(:,:), intent(in)  b,
real(dp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_dlstsq_one()

real(dp) function, dimension(:), allocatable, target la_least_squares::lstsq::la_dlstsq_one ( real(dp), dimension(:,:), intent(inout), target  a,
real(dp), dimension(:), intent(in)  b,
real(dp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_qlstsq_multiple()

real(qp) function, dimension(:,:), allocatable, target la_least_squares::lstsq::la_qlstsq_multiple ( real(qp), dimension(:,:), intent(inout), target  a,
real(qp), dimension(:,:), intent(in)  b,
real(qp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_qlstsq_one()

real(qp) function, dimension(:), allocatable, target la_least_squares::lstsq::la_qlstsq_one ( real(qp), dimension(:,:), intent(inout), target  a,
real(qp), dimension(:), intent(in)  b,
real(qp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_slstsq_multiple()

real(sp) function, dimension(:,:), allocatable, target la_least_squares::lstsq::la_slstsq_multiple ( real(sp), dimension(:,:), intent(inout), target  a,
real(sp), dimension(:,:), intent(in)  b,
real(sp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_slstsq_one()

real(sp) function, dimension(:), allocatable, target la_least_squares::lstsq::la_slstsq_one ( real(sp), dimension(:,:), intent(inout), target  a,
real(sp), dimension(:), intent(in)  b,
real(sp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_wlstsq_multiple()

complex(qp) function, dimension(:,:), allocatable, target la_least_squares::lstsq::la_wlstsq_multiple ( complex(qp), dimension(:,:), intent(inout), target  a,
complex(qp), dimension(:,:), intent(in)  b,
real(qp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_wlstsq_one()

complex(qp) function, dimension(:), allocatable, target la_least_squares::lstsq::la_wlstsq_one ( complex(qp), dimension(:,:), intent(inout), target  a,
complex(qp), dimension(:), intent(in)  b,
real(qp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_zlstsq_multiple()

complex(dp) function, dimension(:,:), allocatable, target la_least_squares::lstsq::la_zlstsq_multiple ( complex(dp), dimension(:,:), intent(inout), target  a,
complex(dp), dimension(:,:), intent(in)  b,
real(dp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

◆ la_zlstsq_one()

complex(dp) function, dimension(:), allocatable, target la_least_squares::lstsq::la_zlstsq_one ( complex(dp), dimension(:,:), intent(inout), target  a,
complex(dp), dimension(:), intent(in)  b,
real(dp), intent(in), optional  cond,
logical(lk), intent(in), optional  overwrite_a,
integer(ilp), intent(out), optional  rank,
type(la_state), intent(out), optional  err 
)
Parameters
[in,out]aInput matrix a[n,n]
[in]bRight hand side vector or array, b[n] or b[n,nrhs]
[in]cond[optional] cutoff for rank evaluation: singular values s(i)<=cond*maxval(s) are considered 0.
[in]overwrite_a[optional] Can A,b data be overwritten and destroyed?
[out]rank[optional] Return rank of A
[out]err[optional] state return flag. On error if not requested, the code will stop
Returns
Result array/matrix x[n] or x[n,nrhs]

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