Operator.Linalg
val inv : Symbol.Shape.Type.arr -> Symbol.Shape.Type.arr
inv a
computes the inverse of the matrix a
. Returns a new array representing the inverse matrix.
val logdet : Symbol.Shape.Type.arr -> Symbol.Shape.Type.elt
logdet a
computes the natural logarithm of the determinant of the matrix a
. Returns the logarithm of the determinant as a scalar.
val chol : ?upper:bool -> Symbol.Shape.Type.arr -> Symbol.Shape.Type.arr
chol ?upper a
performs the Cholesky decomposition of the positive-definite matrix a
.
upper
specifies whether to return the upper or lower triangular matrix. If upper
is true, returns the upper triangular matrix, otherwise the lower triangular matrix. Returns a new array representing the Cholesky factor.val qr : Symbol.Shape.Type.arr -> Symbol.Shape.Type.arr * Symbol.Shape.Type.arr
qr a
performs the QR decomposition of the matrix a
. Returns a tuple of two arrays (Q, R), where Q
is an orthogonal matrix and R
is an upper triangular matrix.
val lq : Symbol.Shape.Type.arr -> Symbol.Shape.Type.arr * Symbol.Shape.Type.arr
lq a
performs the LQ decomposition of the matrix a
. Returns a tuple of two arrays (L, Q), where L
is a lower triangular matrix and Q
is an orthogonal matrix.
val svd :
?thin:bool ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr * Symbol.Shape.Type.arr * Symbol.Shape.Type.arr
svd ?thin a
performs the Singular Value Decomposition (SVD) of the matrix a
.
thin
specifies whether to return the reduced form of the SVD. Returns a tuple of three arrays (U, S, V), where U
and V
are orthogonal matrices, and S
is a diagonal matrix containing the singular values.val sylvester :
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr
sylvester a b c
solves the Sylvester equation A*X + X*B = C for the unknown matrix X. Returns a new array representing the solution matrix X.
val lyapunov :
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr
lyapunov a q
solves the continuous Lyapunov equation A*X + X*A^T = Q for the unknown matrix X. Returns a new array representing the solution matrix X.
val discrete_lyapunov :
?solver:[ `default | `bilinear | `direct ] ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr
discrete_lyapunov ?solver a q
solves the discrete Lyapunov equation A*X*A^T - X + Q = 0 for the unknown matrix X.
solver
specifies the method to use: `default`, `bilinear`, or `direct`. Returns a new array representing the solution matrix X.val linsolve :
?trans:bool ->
?typ:[ `n | `u | `l ] ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr
linsolve ?trans ?typ a b
solves the linear system A*X = B for the unknown matrix X.
trans
specifies whether to transpose the matrix A.typ
specifies the type of matrix A: `n` for normal, `u` for upper triangular, and `l` for lower triangular. Returns a new array representing the solution matrix X.val care :
?diag_r:bool ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr
care ?diag_r a b q r
solves the Continuous-time Algebraic Riccati Equation (CARE) A*X + X*A^T - X*B*R^-1*B^T*X + Q = 0 for the unknown matrix X.
diag_r
if true, R
is assumed to be diagonal. Returns a new array representing the solution matrix X.val dare :
?diag_r:bool ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr ->
Symbol.Shape.Type.arr
dare ?diag_r a b q r
solves the Discrete-time Algebraic Riccati Equation (DARE) A*X*A^T - X - (A*X*B^T)*inv(B*X*B^T + R)*(A*X*B^T)^T + Q = 0 for the unknown matrix X.
diag_r
if true, R
is assumed to be diagonal. Returns a new array representing the solution matrix X.