LMEDenseHessLyapunovChol#

Computes the Cholesky factor of the solution of a dense Lyapunov equation with an upper Hessenberg coefficient matrix.

Synopsis#

#include "slepclme.h" 
PetscErrorCode LMEDenseHessLyapunovChol(LME lme,PetscInt m,PetscScalar *H,PetscInt ldh,PetscInt k,PetscScalar *B,PetscInt ldb,PetscScalar *U,PetscInt ldu,PetscReal *res)

Logically Collective

Input Parameters#

  • lme - linear matrix equation solver context

  • m - number of rows and columns of H

  • H - coefficient matrix

  • ldh - leading dimension of H

  • k - number of columns of B

  • B - right-hand side matrix

  • ldb - leading dimension of B

  • ldu - leading dimension of U

Output Parameters#

  • U - Cholesky factor of the solution

  • res - (optional) residual norm, on input it should contain H(m+1,m)

Note#

The Lyapunov equation has the form HX + XH’ = -B*B’, where H is an mxm upper Hessenberg matrix, B is an mxk matrix and the solution is expressed as X = U’*U, where U is upper triangular. H is supposed to be stable.

When k=1 and the res argument is provided, the last row of X is used to compute the residual norm of a Lyapunov equation projected via Arnoldi.

See Also#

LMEDenseLyapunov(), LMESolve()

Level#

developer

Location#

src/lme/interface/lmedense.c


Index of all LME routines Table of Contents for all manual pages Index of all manual pages