slepc-main 2024-11-22
Report Typos and Errors

LMEDenseLyapunov

Computes the solution of a dense continuous-time Lyapunov equation.

Synopsis

#include "slepclme.h" 
PetscErrorCode LMEDenseLyapunov(LME lme,PetscInt m,PetscScalar *A,PetscInt lda,PetscScalar *B,PetscInt ldb,PetscScalar *X,PetscInt ldx)
Logically Collective

Input Parameters

lme  - linear matrix equation solver context
m  - number of rows and columns of A
A  - coefficient matrix
lda  - leading dimension of A
B  - right-hand side matrix
ldb  - leading dimension of B
ldx  - leading dimension of X

Output Parameter

X  - the solution

Note

The Lyapunov equation has the form A*X + X*A' = -B, where all are mxm matrices, and B is symmetric.

See Also

LMEDenseHessLyapunovChol(), 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