| Package | Description |
|---|---|
| net.dedekind.lapack |
Selected LAPACK Java bindings for Intel's MKL (Math Kernel Library) and the
F2J Fortran to Java translations from net.sourceforge.f2j as a fallback.
|
| net.frobenius.lapack |
LAPACK high-level API (sort of) which provides parameter checking and
automatic workspace allocation which should be a bit more convenient to use.
|
| Modifier and Type | Class and Description |
|---|---|
class |
LapackJ
Netlib implementation. |
class |
LapackN
Intel MKL (Math Kernel Library) native implementation. |
| Modifier and Type | Method and Description |
|---|---|
static Lapack |
Lapack.getInstance()
Returns an instance of a
Lapack implementation. |
static Lapack |
Lapack.getInstance(boolean useMKL)
Returns an instance of a specific
Lapack implementation depending
on the the value of the useMKL parameter. |
| Modifier and Type | Method and Description |
|---|---|
static void |
PlainLapack.cgeev(Lapack la,
TEigJob jobvl,
TEigJob jobvr,
int n,
float[] a,
int lda,
float[] w,
float[] vl,
int ldvl,
float[] vr,
int ldvr) |
static void |
PlainLapack.cgels(Lapack la,
TTrans trans,
int m,
int n,
int rhsCount,
float[] a,
int lda,
float[] b,
int ldb) |
static void |
PlainLapack.cgeqrf(Lapack la,
int m,
int n,
float[] a,
int lda,
float[] tau) |
static void |
PlainLapack.cgesdd(Lapack la,
TSvdJob jobz,
int m,
int n,
float[] a,
int lda,
float[] s,
float[] u,
int ldu,
float[] vt,
int ldvt) |
static void |
PlainLapack.cgesv(Lapack la,
int n,
int rhsCount,
float[] a,
int lda,
int[] indices,
float[] b,
int ldb) |
static void |
PlainLapack.cgetrf(Lapack la,
int m,
int n,
float[] a,
int lda,
int[] indices) |
static void |
PlainLapack.cungqr(Lapack la,
int m,
int n,
int k,
float[] a,
int lda,
float[] tau)
Complex counterpart of
PlainLapack.sorgqr(net.dedekind.lapack.Lapack, int, int, int, float[], int, float[]). |
static double |
PlainLapack.dgbcon(Lapack la,
TNorm norm,
int n,
int kl,
int ku,
double[] ab,
int[] indices,
double normA)
Purpose
=======
DGBCON estimates the reciprocal of the condition number of a real
general band matrix A, in either the 1-norm or the infinity-norm,
using the LU factorization computed by DGBTRF. |
static void |
PlainLapack.dgbsv(Lapack la,
int n,
int kl,
int ku,
int rhsCount,
double[] ab,
int[] indices,
double[] b,
int ldb)
Purpose
=======
DGBSV computes the solution to a real system of linear equations
A * X = B, where A is a band matrix of order N with KL subdiagonals
and KU superdiagonals, and X and B are N-by-NRHS matrices. |
static void |
PlainLapack.dgbtrf(Lapack la,
int m,
int n,
int kl,
int ku,
double[] ab,
int[] indices)
Purpose
=======
DGBTRF computes an LU factorization of a real m-by-n band matrix A
using partial pivoting with row interchanges. |
static void |
PlainLapack.dgbtrs(Lapack la,
TTrans trans,
int n,
int kl,
int ku,
int rhsCount,
double[] ab,
int[] indices,
double[] b,
int ldb)
Purpose
=======
DGBTRS solves a system of linear equations
A * X = B or A' * X = B
with a general band matrix A using the LU factorization computed
by DGBTRF. |
static double |
PlainLapack.dgecon(Lapack la,
TNorm norm,
int n,
double[] a,
int lda,
double normA)
Purpose
=======
DGECON estimates the reciprocal of the condition number of a general
real matrix A, in either the 1-norm or the infinity-norm, using
the LU factorization computed by DGETRF. |
static void |
PlainLapack.dgeev(Lapack la,
TEigJob jobvl,
TEigJob jobvr,
int n,
double[] a,
int lda,
double[] wr,
double[] wi,
double[] vl,
int ldvl,
double[] vr,
int ldvr)
Purpose
=======
DGEEV computes for an N-by-N real nonsymmetric matrix A, the
eigenvalues and, optionally, the left and/or right eigenvectors. |
static void |
PlainLapack.dgelqf(Lapack la,
int m,
int n,
double[] a,
int lda,
double[] tau)
Purpose
=======
DGELQF computes an LQ factorization of a real M-by-N matrix A:
A = L * Q. |
static void |
PlainLapack.dgels(Lapack la,
TTrans trans,
int m,
int n,
int rhsCount,
double[] a,
int lda,
double[] b,
int ldb)
Purpose
=======
DGELS solves overdetermined or underdetermined real linear systems
involving an M-by-N matrix A, or its transpose, using a QR or LQ
factorization of A. |
static void |
PlainLapack.dgeqlf(Lapack la,
int m,
int n,
double[] a,
int lda,
double[] tau)
Purpose
=======
DGEQLF computes a QL factorization of a real M-by-N matrix A:
A = Q * L. |
static void |
PlainLapack.dgeqp3(Lapack la,
int m,
int n,
double[] a,
int lda,
int[] jPivot,
double[] tau)
Purpose
=======
DGEQP3 computes a QR factorization with column pivoting of a
matrix A: A*P = Q*R using Level 3 BLAS. |
static void |
PlainLapack.dgeqrf(Lapack la,
int m,
int n,
double[] a,
int lda,
double[] tau)
Purpose
=======
DGEQRF computes a QR factorization of a real M-by-N matrix A:
A = Q * R. |
static void |
PlainLapack.dgerqf(Lapack la,
int m,
int n,
double[] a,
int lda,
double[] tau)
Purpose
=======
DGERQF computes an RQ factorization of a real M-by-N matrix A:
A = R * Q. |
static void |
PlainLapack.dgesdd(Lapack la,
TSvdJob jobz,
int m,
int n,
double[] a,
int lda,
double[] s,
double[] u,
int ldu,
double[] vt,
int ldvt)
Purpose
=======
DGESDD computes the singular value decomposition (SVD) of a real
M-by-N matrix A, optionally computing the left and right singular
vectors. |
static void |
PlainLapack.dgesv(Lapack la,
int n,
int rhsCount,
double[] a,
int lda,
int[] indices,
double[] b,
int ldb)
Purpose
=======
DGESV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N matrix and X and B are N-by-NRHS matrices. |
static void |
PlainLapack.dgetrf(Lapack la,
int m,
int n,
double[] a,
int lda,
int[] indices)
Purpose
=======
DGETRF computes an LU factorization of a general M-by-N matrix A
using partial pivoting with row interchanges. |
static void |
PlainLapack.dgetrs(Lapack la,
TTrans trans,
int n,
int rhsCount,
double[] a,
int lda,
int[] indices,
double[] b,
int ldb)
Purpose
=======
DGETRS solves a system of linear equations
A * X = B or A' * X = B
with a general N-by-N matrix A using the LU factorization computed
by DGETRF. |
static void |
PlainLapack.dgtsv(Lapack la,
int n,
int rhsCount,
double[] dl,
double[] d,
double[] du,
double[] b,
int ldb)
Purpose
=======
DGTSV solves the equation
A*X = B,
where A is an n by n tridiagonal matrix, by Gaussian elimination with
partial pivoting. |
static void |
PlainLapack.dlaswp(Lapack la,
int n,
double[] a,
int lda,
int pivFirstIdx,
int pivLastIdx,
int[] indices,
int increment)
Purpose
=======
DLASWP performs a series of row interchanges on the matrix A. |
static void |
PlainLapack.dorglq(Lapack la,
int m,
int n,
int k,
double[] a,
int lda,
double[] tau)
Purpose
=======
DORGLQ generates an M-by-N real matrix Q with orthonormal rows,
which is defined as the first M rows of a product of K elementary
reflectors of order N
Q = H(k) . . . |
static void |
PlainLapack.dorgqr(Lapack la,
int m,
int n,
int k,
double[] a,
int lda,
double[] tau)
Purpose
=======
DORGQR generates an M-by-N real matrix Q with orthonormal columns,
which is defined as the first N columns of a product of K elementary
reflectors of order M
Q = H(1) H(2) . . . |
static void |
PlainLapack.dorgrq(Lapack la,
int m,
int n,
int k,
double[] a,
int lda,
double[] tau)
Purpose
=======
DORGRQ generates an M-by-N real matrix Q with orthonormal rows,
which is defined as the last M rows of a product of K elementary
reflectors of order N
Q = H(1) H(2) . . . |
static void |
PlainLapack.dormrz(Lapack la,
TSide side,
TTrans trans,
int m,
int n,
int k,
int l,
double[] a,
int lda,
double[] tau,
double[] c,
int ldc)
Purpose
=======
DORMRZ overwrites the general real M-by-N matrix C with
SIDE = 'L' SIDE = 'R'
TRANS = 'N': Q * C C * Q
TRANS = 'T': Q**T * C C * Q**T
where Q is a real orthogonal matrix defined as the product of k
elementary reflectors
Q = H(1) H(2) . . . |
static double |
PlainLapack.dpbcon(Lapack la,
TUpLo uplo,
int n,
int diagCount,
double[] ab,
double normA)
Purpose
=======
DPBCON estimates the reciprocal of the condition number (in the
1-norm) of a real symmetric positive definite band matrix using the
Cholesky factorization A = U**T*U or A = L*L**T computed by DPBTRF. |
static void |
PlainLapack.dpbsv(Lapack la,
TUpLo uplo,
int n,
int diagCount,
int rhsCount,
double[] ab,
double[] b,
int ldb)
Purpose
=======
DPBSV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric positive definite band matrix and X
and B are N-by-NRHS matrices. |
static void |
PlainLapack.dpbtrf(Lapack la,
TUpLo uplo,
int n,
int diagCount,
double[] ab)
Purpose
=======
DPBTRF computes the Cholesky factorization of a real symmetric
positive definite band matrix A. |
static void |
PlainLapack.dpbtrs(Lapack la,
TUpLo uplo,
int n,
int diagCount,
int rhsCount,
double[] ab,
double[] b,
int ldb)
Purpose
=======
DPBTRS solves a system of linear equations A*X = B with a symmetric
positive definite band matrix A using the Cholesky factorization
A = U**T*U or A = L*L**T computed by DPBTRF. |
static double |
PlainLapack.dpocon(Lapack la,
TUpLo uplo,
int n,
double[] a,
int lda,
double normA)
Purpose
=======
DPOCON estimates the reciprocal of the condition number (in the
1-norm) of a real symmetric positive definite matrix using the
Cholesky factorization A = U**T*U or A = L*L**T computed by DPOTRF. |
static void |
PlainLapack.dposv(Lapack la,
TUpLo uplo,
int n,
int rhsCount,
double[] a,
int lda,
double[] b,
int ldb)
Purpose
=======
DPOSV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric positive definite matrix and X and B
are N-by-NRHS matrices. |
static void |
PlainLapack.dpotrf(Lapack la,
TUpLo uplo,
int n,
double[] a,
int lda)
Purpose
=======
DPOTRF computes the Cholesky factorization of a real symmetric
positive definite matrix A. |
static void |
PlainLapack.dpotrs(Lapack la,
TUpLo uplo,
int n,
int rhsCount,
double[] a,
int lda,
double[] b,
int ldb)
Purpose
=======
DPOTRS solves a system of linear equations A*X = B with a symmetric
positive definite matrix A using the Cholesky factorization
A = U**T*U or A = L*L**T computed by DPOTRF. |
static double |
PlainLapack.dppcon(Lapack la,
TUpLo uplo,
int n,
double[] ap,
double normA)
Purpose
=======
DPPCON estimates the reciprocal of the condition number (in the
1-norm) of a real symmetric positive definite packed matrix using
the Cholesky factorization A = U**T*U or A = L*L**T computed by
DPPTRF. |
static void |
PlainLapack.dppsv(Lapack la,
TUpLo uplo,
int n,
int rhsCount,
double[] ap,
double[] b,
int ldb)
Purpose
=======
DPPSV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric positive definite matrix stored in
packed format and X and B are N-by-NRHS matrices. |
static void |
PlainLapack.dpptrf(Lapack la,
TUpLo uplo,
int n,
double[] ap)
Purpose
=======
DPPTRF computes the Cholesky factorization of a real symmetric
positive definite matrix A stored in packed format. |
static void |
PlainLapack.dpptrs(Lapack la,
TUpLo uplo,
int n,
int rhsCount,
double[] ap,
double[] b,
int ldb)
Purpose
=======
DPPTRS solves a system of linear equations A*X = B with a symmetric
positive definite matrix A in packed storage using the Cholesky
factorization A = U**T*U or A = L*L**T computed by DPPTRF. |
static void |
PlainLapack.dptsv(Lapack la,
int n,
int rhsCount,
double[] d,
double[] e,
double[] b,
int ldb)
Purpose
=======
DPTSV computes the solution to a real system of linear equations
A*X = B, where A is an N-by-N symmetric positive definite tridiagonal
matrix, and X and B are N-by-NRHS matrices. |
static void |
PlainLapack.dsbevd(Lapack la,
TEigJob jobz,
TUpLo uplo,
int n,
int diagCount,
double[] ab,
double[] w,
double[] z,
int ldz)
Purpose
=======
DSBEVD computes all the eigenvalues and, optionally, eigenvectors of
a real symmetric band matrix A. |
static void |
PlainLapack.dspevd(Lapack la,
TEigJob jobz,
TUpLo uplo,
int n,
double[] ap,
double[] w,
double[] z,
int ldz)
Purpose
=======
DSPEVD computes all the eigenvalues and, optionally, eigenvectors
of a real symmetric matrix A in packed storage. |
static void |
PlainLapack.dspsv(Lapack la,
TUpLo uplo,
int n,
int rhsCount,
double[] ap,
int[] indices,
double[] b,
int ldb)
Purpose
=======
DSPSV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric matrix stored in packed format and X
and B are N-by-NRHS matrices. |
static int |
PlainLapack.dstevr(Lapack la,
TEigJob jobz,
TRange range,
int n,
double[] d,
double[] e,
double vLower,
double vUpper,
int iLower,
int iUpper,
double abstol,
double[] w,
double[] z,
int ldz,
int[] supportZ)
Purpose
=======
DSTEVR computes selected eigenvalues and, optionally, eigenvectors
of a real symmetric tridiagonal matrix T. |
static int |
PlainLapack.dsyevr(Lapack la,
TEigJob jobz,
TRange range,
TUpLo uplo,
int n,
double[] a,
int lda,
double vLower,
double vUpper,
int iLower,
int iUpper,
double abstol,
double[] w,
double[] z,
int ldz,
int[] supportZ)
Purpose
=======
DSYEVR computes selected eigenvalues and, optionally, eigenvectors
of a real symmetric matrix A. |
static void |
PlainLapack.dsygvd(Lapack la,
int type,
TEigJob jobz,
TUpLo uplo,
int n,
double[] a,
int lda,
double[] b,
int ldb,
double[] w)
Purpose
=======
DSYGVD computes all the eigenvalues, and optionally, the eigenvectors
of a real generalized symmetric-definite eigenproblem, of the form
A*x=(lambda)*B*x, A*Bx=(lambda)*x, or B*A*x=(lambda)*x. |
static void |
PlainLapack.dsysv(Lapack la,
TUpLo uplo,
int n,
int rhsCount,
double[] a,
int lda,
int[] indices,
double[] b,
int ldb)
Purpose
=======
DSYSV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N symmetric matrix and X and B are N-by-NRHS
matrices. |
static void |
PlainLapack.dtbtrs(Lapack la,
TUpLo uplo,
TTrans trans,
TDiag diag,
int n,
int diagCount,
int rhsCount,
double[] ab,
double[] b,
int ldb)
Purpose
=======
DTBTRS solves a triangular system of the form
A * X = B or A**T * X = B,
where A is a triangular band matrix of order N, and B is an
N-by NRHS matrix. |
static void |
PlainLapack.dtptrs(Lapack la,
TUpLo uplo,
TTrans trans,
TDiag diag,
int n,
int rhsCount,
double[] ap,
double[] b,
int ldb)
Purpose
=======
DTPTRS solves a triangular system of the form
A * X = B or A**T * X = B,
where A is a triangular matrix of order N stored in packed format,
and B is an N-by-NRHS matrix. |
static void |
PlainLapack.dtrtrs(Lapack la,
TUpLo uplo,
TTrans trans,
TDiag diag,
int n,
int rhsCount,
double[] a,
int lda,
double[] b,
int ldb)
Purpose
=======
DTRTRS solves a triangular system of the form
A * X = B or A**T * X = B,
where A is a triangular matrix of order N, and B is an N-by-NRHS
matrix. |
static void |
PlainLapack.sgeev(Lapack la,
TEigJob jobvl,
TEigJob jobvr,
int n,
float[] a,
int lda,
float[] wr,
float[] wi,
float[] vl,
int ldvl,
float[] vr,
int ldvr)
Purpose
=======
SGEEV computes for an N-by-N real nonsymmetric matrix A, the
eigenvalues and, optionally, the left and/or right eigenvectors. |
static void |
PlainLapack.sgels(Lapack la,
TTrans trans,
int m,
int n,
int rhsCount,
float[] a,
int lda,
float[] b,
int ldb)
Purpose
=======
SGELS solves overdetermined or underdetermined real linear systems
involving an M-by-N matrix A, or its transpose, using a QR or LQ
factorization of A. |
static void |
PlainLapack.sgeqrf(Lapack la,
int m,
int n,
float[] a,
int lda,
float[] tau)
Purpose
=======
SGEQRF computes a QR factorization of a real M-by-N matrix A:
A = Q * R. |
static void |
PlainLapack.sgesdd(Lapack la,
TSvdJob jobz,
int m,
int n,
float[] a,
int lda,
float[] s,
float[] u,
int ldu,
float[] vt,
int ldvt)
Purpose
=======
SGESDD computes the singular value decomposition (SVD) of a real
M-by-N matrix A, optionally computing the left and right singular
vectors. |
static void |
PlainLapack.sgesv(Lapack la,
int n,
int rhsCount,
float[] a,
int lda,
int[] indices,
float[] b,
int ldb)
Purpose
=======
SGESV computes the solution to a real system of linear equations
A * X = B,
where A is an N-by-N matrix and X and B are N-by-NRHS matrices. |
static void |
PlainLapack.sgetrf(Lapack la,
int m,
int n,
float[] a,
int lda,
int[] indices)
Purpose
=======
SGETRF computes an LU factorization of a general M-by-N matrix A
using partial pivoting with row interchanges. |
static void |
PlainLapack.sorgqr(Lapack la,
int m,
int n,
int k,
float[] a,
int lda,
float[] tau)
Purpose
=======
SORGQR generates an M-by-N real matrix Q with orthonormal columns,
which is defined as the first N columns of a product of K elementary
reflectors of order M
Q = H(1) H(2) . . . |
static void |
PlainLapack.zgeev(Lapack la,
TEigJob jobvl,
TEigJob jobvr,
int n,
double[] a,
int lda,
double[] w,
double[] vl,
int ldvl,
double[] vr,
int ldvr) |
static void |
PlainLapack.zgels(Lapack la,
TTrans trans,
int m,
int n,
int rhsCount,
double[] a,
int lda,
double[] b,
int ldb) |
static void |
PlainLapack.zgeqrf(Lapack la,
int m,
int n,
double[] a,
int lda,
double[] tau) |
static void |
PlainLapack.zgesdd(Lapack la,
TSvdJob jobz,
int m,
int n,
double[] a,
int lda,
double[] s,
double[] u,
int ldu,
double[] vt,
int ldvt) |
static void |
PlainLapack.zgesv(Lapack la,
int n,
int rhsCount,
double[] a,
int lda,
int[] indices,
double[] b,
int ldb) |
static void |
PlainLapack.zgetrf(Lapack la,
int m,
int n,
double[] a,
int lda,
int[] indices) |
static void |
PlainLapack.zungqr(Lapack la,
int m,
int n,
int k,
double[] a,
int lda,
double[] tau)
Complex counterpart of
PlainLapack.dorgqr(net.dedekind.lapack.Lapack, int, int, int, double[], int, double[]). |
Copyright © 2023. All rights reserved.