.Modelica.Math.Matrices.cholesky

Information

Syntax

H = Matrices.cholesky(A);
H = Matrices.cholesky(A, upper=true);

Description

Function cholesky computes the Cholesky factorization of a real symmetric positive definite matrix A. The optional Boolean input "upper" specifies whether the upper or the lower triangular matrix is returned, i.e.

A = H'*H   if upper is true (H is upper triangular)
A = H*H'   if upper is false (H is lower triangular)

The computation is performed by LAPACK.dpotrf.

Example

A  = [1, 0,  0;
      6, 5,  0;
      1, -2,  2];
S = A*transpose(A);

H = Matrices.cholesky(S);

results in:

H = [1.0,  6.0,  1.0;
     0.0,  5.0, -2.0;
     0.0,  0.0,  2.0]

with

transpose(H)*H = [1.0,  6.0,   1;
                  6.0, 61.0,  -4.0;
                  1.0, -4.0,   9.0] //=S

Interface

function cholesky
  extends Modelica.Icons.Function;
  import Modelica.Math.Matrices.LAPACK;
  input Real A[:, size(A, 1)] "Symmetric positive definite matrix";
  input Boolean upper = true "= true, if the right Cholesky factor (upper triangle) should be returned";
  output Real H[size(A, 1), size(A, 2)] "Cholesky factor U (upper=true) or L (upper=false) for A = U'*U or A = L*L'";
end cholesky;

Revisions


Generated at 2024-11-23T19:25:52Z by OpenModelicaOpenModelica 1.24.2 using GenerateDoc.mos