Actual source code: slepclme.h
slepc-main 2024-11-15
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: User interface for the SLEPc object for solving linear matrix equations
12: */
14: #pragma once
16: #include <slepcbv.h>
18: /* SUBMANSEC = LME */
20: SLEPC_EXTERN PetscErrorCode LMEInitializePackage(void);
21: SLEPC_EXTERN PetscErrorCode LMEFinalizePackage(void);
23: /*S
24: LME - SLEPc object that encapsulates functionality for linear matrix equations
26: Level: beginner
28: .seealso: LMECreate()
29: S*/
30: typedef struct _p_LME* LME;
32: /*J
33: LMEType - String with the name of a method for solving linear matrix equations
35: Level: beginner
37: .seealso: LMESetType(), LME
38: J*/
39: typedef const char* LMEType;
40: #define LMEKRYLOV "krylov"
42: /* Logging support */
43: SLEPC_EXTERN PetscClassId LME_CLASSID;
45: /*E
46: LMEProblemType - Determines the type of linear matrix equation
48: Level: beginner
50: .seealso: LMESetProblemType(), LMEGetProblemType()
51: E*/
52: typedef enum { LME_LYAPUNOV,
53: LME_SYLVESTER,
54: LME_GEN_LYAPUNOV,
55: LME_GEN_SYLVESTER,
56: LME_DT_LYAPUNOV ,
57: LME_STEIN} LMEProblemType;
58: SLEPC_EXTERN const char *LMEProblemTypes[];
60: SLEPC_EXTERN PetscErrorCode LMECreate(MPI_Comm,LME *);
61: SLEPC_EXTERN PetscErrorCode LMEDestroy(LME*);
62: SLEPC_EXTERN PetscErrorCode LMEReset(LME);
63: SLEPC_EXTERN PetscErrorCode LMESetType(LME,LMEType);
64: SLEPC_EXTERN PetscErrorCode LMEGetType(LME,LMEType*);
65: SLEPC_EXTERN PetscErrorCode LMESetProblemType(LME,LMEProblemType);
66: SLEPC_EXTERN PetscErrorCode LMEGetProblemType(LME,LMEProblemType*);
67: SLEPC_EXTERN PetscErrorCode LMESetCoefficients(LME,Mat,Mat,Mat,Mat);
68: SLEPC_EXTERN PetscErrorCode LMEGetCoefficients(LME,Mat*,Mat*,Mat*,Mat*);
69: SLEPC_EXTERN PetscErrorCode LMESetRHS(LME,Mat);
70: SLEPC_EXTERN PetscErrorCode LMEGetRHS(LME,Mat*);
71: SLEPC_EXTERN PetscErrorCode LMESetSolution(LME,Mat);
72: SLEPC_EXTERN PetscErrorCode LMEGetSolution(LME,Mat*);
73: SLEPC_EXTERN PetscErrorCode LMESetFromOptions(LME);
74: SLEPC_EXTERN PetscErrorCode LMESetUp(LME);
75: SLEPC_EXTERN PetscErrorCode LMESolve(LME);
76: SLEPC_EXTERN PetscErrorCode LMEView(LME,PetscViewer);
77: SLEPC_EXTERN PetscErrorCode LMEViewFromOptions(LME,PetscObject,const char[]);
78: SLEPC_EXTERN PetscErrorCode LMEConvergedReasonView(LME,PetscViewer);
79: SLEPC_EXTERN PetscErrorCode LMEConvergedReasonViewFromOptions(LME);
80: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "LMEConvergedReasonView()", ) static inline PetscErrorCode LMEReasonView(LME lme,PetscViewer v) {return LMEConvergedReasonView(lme,v);}
81: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "LMEConvergedReasonViewFromOptions()", ) static inline PetscErrorCode LMEReasonViewFromOptions(LME lme) {return LMEConvergedReasonViewFromOptions(lme);}
83: SLEPC_EXTERN PetscErrorCode LMESetBV(LME,BV);
84: SLEPC_EXTERN PetscErrorCode LMEGetBV(LME,BV*);
85: SLEPC_EXTERN PetscErrorCode LMESetTolerances(LME,PetscReal,PetscInt);
86: SLEPC_EXTERN PetscErrorCode LMEGetTolerances(LME,PetscReal*,PetscInt*);
87: SLEPC_EXTERN PetscErrorCode LMESetDimensions(LME,PetscInt);
88: SLEPC_EXTERN PetscErrorCode LMEGetDimensions(LME,PetscInt*);
89: SLEPC_EXTERN PetscErrorCode LMEGetIterationNumber(LME,PetscInt*);
91: SLEPC_EXTERN PetscErrorCode LMEGetErrorEstimate(LME,PetscReal*);
92: SLEPC_EXTERN PetscErrorCode LMEComputeError(LME,PetscReal*);
93: SLEPC_EXTERN PetscErrorCode LMESetErrorIfNotConverged(LME,PetscBool);
94: SLEPC_EXTERN PetscErrorCode LMEGetErrorIfNotConverged(LME,PetscBool*);
96: SLEPC_EXTERN PetscErrorCode LMEDenseLyapunov(LME,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt);
97: SLEPC_EXTERN PetscErrorCode LMEDenseHessLyapunovChol(LME,PetscInt,PetscScalar*,PetscInt,PetscInt,PetscScalar*,PetscInt,PetscScalar*,PetscInt,PetscReal*);
98: PETSC_DEPRECATED_FUNCTION(3, 8, 0, "LMEDenseHessLyapunovChol()", ) static inline PetscErrorCode LMEDenseLyapunovChol(LME lme,PetscScalar *H,PetscInt m,PetscInt ldh,PetscScalar *r,PetscScalar *L,PetscInt ldl,PetscReal *res) {return LMEDenseHessLyapunovChol(lme,m,H,ldh,1,r,m,L,ldl,res);}
100: SLEPC_EXTERN PetscErrorCode LMEMonitor(LME,PetscInt,PetscReal);
101: SLEPC_EXTERN PetscErrorCode LMEMonitorSet(LME,PetscErrorCode (*)(LME,PetscInt,PetscReal,void*),void*,PetscCtxDestroyFn*);
102: SLEPC_EXTERN PetscErrorCode LMEMonitorCancel(LME);
103: SLEPC_EXTERN PetscErrorCode LMEGetMonitorContext(LME,void*);
105: SLEPC_EXTERN PetscErrorCode LMEMonitorSetFromOptions(LME,const char[],const char[],void*);
106: SLEPC_EXTERN PetscErrorCode LMEMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char*[],int,int,int,int,PetscDrawLG*);
107: SLEPC_EXTERN PetscErrorCode LMEMonitorDefault(LME,PetscInt,PetscReal,PetscViewerAndFormat*);
108: SLEPC_EXTERN PetscErrorCode LMEMonitorDefaultDrawLG(LME,PetscInt,PetscReal,PetscViewerAndFormat*);
109: SLEPC_EXTERN PetscErrorCode LMEMonitorDefaultDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat**);
111: SLEPC_EXTERN PetscErrorCode LMESetOptionsPrefix(LME,const char*);
112: SLEPC_EXTERN PetscErrorCode LMEAppendOptionsPrefix(LME,const char*);
113: SLEPC_EXTERN PetscErrorCode LMEGetOptionsPrefix(LME,const char*[]);
115: /*E
116: LMEConvergedReason - reason a matrix function iteration was said to
117: have converged or diverged
119: Level: intermediate
121: .seealso: LMESolve(), LMEGetConvergedReason(), LMESetTolerances()
122: E*/
123: typedef enum {/* converged */
124: LME_CONVERGED_TOL = 1,
125: /* diverged */
126: LME_DIVERGED_ITS = -1,
127: LME_DIVERGED_BREAKDOWN = -2,
128: LME_CONVERGED_ITERATING = 0} LMEConvergedReason;
129: SLEPC_EXTERN const char *const*LMEConvergedReasons;
131: SLEPC_EXTERN PetscErrorCode LMEGetConvergedReason(LME,LMEConvergedReason *);
133: SLEPC_EXTERN PetscFunctionList LMEList;
134: SLEPC_EXTERN PetscFunctionList LMEMonitorList;
135: SLEPC_EXTERN PetscFunctionList LMEMonitorCreateList;
136: SLEPC_EXTERN PetscFunctionList LMEMonitorDestroyList;
137: SLEPC_EXTERN PetscErrorCode LMERegister(const char[],PetscErrorCode(*)(LME));
138: SLEPC_EXTERN PetscErrorCode LMEMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,PetscErrorCode(*)(LME,PetscInt,PetscReal,PetscViewerAndFormat*),PetscErrorCode(*)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**),PetscErrorCode(*)(PetscViewerAndFormat**));
140: SLEPC_EXTERN PetscErrorCode LMEAllocateSolution(LME,PetscInt);