Actual source code: lmesetup.c
slepc-3.18.1 2022-11-02
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: LME routines related to problem setup
12: */
14: #include <slepc/private/lmeimpl.h>
16: static inline PetscErrorCode LMESetUp_Lyapunov(LME lme)
17: {
18: Mat C1,C2,X1,X2;
19: Vec dc,dx;
21: MatLRCGetMats(lme->C,NULL,&C1,&dc,&C2);
24: if (lme->X) {
25: MatLRCGetMats(lme->X,NULL,&X1,&dx,&X2);
28: }
29: return 0;
30: }
32: /*@
33: LMESetUp - Sets up all the internal data structures necessary for the
34: execution of the linear matrix equation solver.
36: Collective on lme
38: Input Parameter:
39: . lme - linear matrix equation solver context
41: Notes:
42: This function need not be called explicitly in most cases, since LMESolve()
43: calls it. It can be useful when one wants to measure the set-up time
44: separately from the solve time.
46: Level: developer
48: .seealso: LMECreate(), LMESolve(), LMEDestroy()
49: @*/
50: PetscErrorCode LMESetUp(LME lme)
51: {
52: PetscInt N;
56: /* reset the convergence flag from the previous solves */
57: lme->reason = LME_CONVERGED_ITERATING;
59: if (lme->setupcalled) return 0;
60: PetscLogEventBegin(LME_SetUp,lme,0,0,0);
62: /* Set default solver type (LMESetFromOptions was not called) */
63: if (!((PetscObject)lme)->type_name) LMESetType(lme,LMEKRYLOV);
65: /* Check problem dimensions */
67: MatGetSize(lme->A,&N,NULL);
68: if (lme->ncv > N) lme->ncv = N;
70: /* setup options for the particular equation type */
71: switch (lme->problem_type) {
72: case LME_LYAPUNOV:
73: LMESetUp_Lyapunov(lme);
74: break;
75: case LME_SYLVESTER:
76: LMECheckCoeff(lme,lme->B,"B","Sylvester");
77: break;
78: case LME_GEN_LYAPUNOV:
79: LMECheckCoeff(lme,lme->D,"D","Generalized Lyapunov");
80: break;
81: case LME_GEN_SYLVESTER:
82: LMECheckCoeff(lme,lme->B,"B","Generalized Sylvester");
83: LMECheckCoeff(lme,lme->D,"D","Generalized Sylvester");
84: LMECheckCoeff(lme,lme->E,"E","Generalized Sylvester");
85: break;
86: case LME_DT_LYAPUNOV:
87: break;
88: case LME_STEIN:
89: LMECheckCoeff(lme,lme->D,"D","Stein");
90: break;
91: }
94: /* call specific solver setup */
95: PetscUseTypeMethod(lme,setup);
97: /* set tolerance if not yet set */
98: if (lme->tol==PETSC_DEFAULT) lme->tol = SLEPC_DEFAULT_TOL;
100: PetscLogEventEnd(LME_SetUp,lme,0,0,0);
101: lme->setupcalled = 1;
102: return 0;
103: }
105: static inline PetscErrorCode LMESetCoefficients_Private(LME lme,Mat A,Mat *lmeA)
106: {
107: PetscInt m,n;
109: MatGetSize(A,&m,&n);
111: if (!lme->setupcalled) MatDestroy(lmeA);
112: PetscObjectReference((PetscObject)A);
113: *lmeA = A;
114: return 0;
115: }
117: /*@
118: LMESetCoefficients - Sets the coefficient matrices that define the linear matrix
119: equation to be solved.
121: Collective on lme
123: Input Parameters:
124: + lme - the matrix function context
125: . A - first coefficient matrix
126: . B - second coefficient matrix
127: . D - third coefficient matrix
128: - E - fourth coefficient matrix
130: Notes:
131: The matrix equation takes the general form A*X*E+D*X*B=C, where matrix C is not
132: provided here but with LMESetRHS(). Not all four matrices must be passed, some
133: can be NULL instead, see LMESetProblemType() for details.
135: It must be called before LMESetUp(). If it is called again after LMESetUp() then
136: the LME object is reset.
138: In order to delete a previously set matrix, pass a NULL in the corresponding
139: argument.
141: Level: beginner
143: .seealso: LMESolve(), LMESetUp(), LMESetRHS(), LMESetProblemType()
144: @*/
145: PetscErrorCode LMESetCoefficients(LME lme,Mat A,Mat B,Mat D,Mat E)
146: {
150: if (B) {
153: }
154: if (D) {
157: }
158: if (E) {
161: }
163: if (lme->setupcalled) LMEReset(lme);
165: LMESetCoefficients_Private(lme,A,&lme->A);
166: if (B) LMESetCoefficients_Private(lme,B,&lme->B);
167: else if (!lme->setupcalled) MatDestroy(&lme->B);
168: if (D) LMESetCoefficients_Private(lme,D,&lme->D);
169: else if (!lme->setupcalled) MatDestroy(&lme->D);
170: if (E) LMESetCoefficients_Private(lme,E,&lme->E);
171: else if (!lme->setupcalled) MatDestroy(&lme->E);
173: lme->setupcalled = 0;
174: return 0;
175: }
177: /*@
178: LMEGetCoefficients - Gets the coefficient matrices of the matrix equation.
180: Collective on lme
182: Input Parameter:
183: . lme - the LME context
185: Output Parameters:
186: + A - first coefficient matrix
187: . B - second coefficient matrix
188: . D - third coefficient matrix
189: - E - fourth coefficient matrix
191: Level: intermediate
193: .seealso: LMESolve(), LMESetCoefficients()
194: @*/
195: PetscErrorCode LMEGetCoefficients(LME lme,Mat *A,Mat *B,Mat *D,Mat *E)
196: {
198: if (A) *A = lme->A;
199: if (B) *B = lme->B;
200: if (D) *D = lme->D;
201: if (E) *E = lme->E;
202: return 0;
203: }
205: /*@
206: LMESetRHS - Sets the right-hand side of the matrix equation, as a low-rank
207: matrix.
209: Collective on lme
211: Input Parameters:
212: + lme - the matrix function context
213: - C - the right-hand side matrix
215: Notes:
216: The matrix equation takes the general form A*X*E+D*X*B=C, where matrix C is
217: given with this function. C must be a low-rank matrix of type MATLRC, that is,
218: C = U*D*V' where D is diagonal of order k, and U, V are dense tall-skinny
219: matrices with k columns. No sparse matrix must be provided when creating the
220: MATLRC matrix.
222: In equation types that require C to be symmetric, such as Lyapunov, C must be
223: created with V=U (or V=NULL).
225: Level: beginner
227: .seealso: LMESetSolution(), LMESetProblemType()
228: @*/
229: PetscErrorCode LMESetRHS(LME lme,Mat C)
230: {
231: Mat A;
238: MatLRCGetMats(C,&A,NULL,NULL,NULL);
241: PetscObjectReference((PetscObject)C);
242: MatDestroy(&lme->C);
243: lme->C = C;
244: return 0;
245: }
247: /*@
248: LMEGetRHS - Gets the right-hand side of the matrix equation.
250: Collective on lme
252: Input Parameter:
253: . lme - the LME context
255: Output Parameters:
256: . C - the low-rank matrix
258: Level: intermediate
260: .seealso: LMESolve(), LMESetRHS()
261: @*/
262: PetscErrorCode LMEGetRHS(LME lme,Mat *C)
263: {
266: *C = lme->C;
267: return 0;
268: }
270: /*@
271: LMESetSolution - Sets the placeholder for the solution of the matrix
272: equation, as a low-rank matrix.
274: Collective on lme
276: Input Parameters:
277: + lme - the matrix function context
278: - X - the solution matrix
280: Notes:
281: The matrix equation takes the general form A*X*E+D*X*B=C, where the solution
282: matrix is of low rank and is written in factored form X = U*D*V'. This function
283: provides a Mat object of type MATLRC that stores U, V and (optionally) D.
284: These factors will be computed during LMESolve().
286: In equation types whose solution X is symmetric, such as Lyapunov, X must be
287: created with V=U (or V=NULL).
289: If the user provides X with this function, then the solver will
290: return a solution with rank at most the number of columns of U. Alternatively,
291: it is possible to let the solver choose the rank of the solution, by
292: setting X to NULL and then calling LMEGetSolution() after LMESolve().
294: Level: intermediate
296: .seealso: LMEGetSolution(), LMESetRHS(), LMESetProblemType(), LMESolve()
297: @*/
298: PetscErrorCode LMESetSolution(LME lme,Mat X)
299: {
300: Mat A;
303: if (X) {
307: MatLRCGetMats(X,&A,NULL,NULL,NULL);
309: PetscObjectReference((PetscObject)X);
310: }
311: MatDestroy(&lme->X);
312: lme->X = X;
313: return 0;
314: }
316: /*@
317: LMEGetSolution - Gets the solution of the matrix equation.
319: Collective on lme
321: Input Parameter:
322: . lme - the LME context
324: Output Parameters:
325: . X - the low-rank matrix
327: Level: intermediate
329: .seealso: LMESolve(), LMESetSolution()
330: @*/
331: PetscErrorCode LMEGetSolution(LME lme,Mat *X)
332: {
335: *X = lme->X;
336: return 0;
337: }
339: /*@
340: LMEAllocateSolution - Allocate memory storage for common variables such
341: as the basis vectors.
343: Collective on lme
345: Input Parameters:
346: + lme - linear matrix equation solver context
347: - extra - number of additional positions, used for methods that require a
348: working basis slightly larger than ncv
350: Developer Notes:
351: This is SLEPC_EXTERN because it may be required by user plugin LME
352: implementations.
354: Level: developer
356: .seealso: LMESetUp()
357: @*/
358: PetscErrorCode LMEAllocateSolution(LME lme,PetscInt extra)
359: {
360: PetscInt oldsize,requested;
361: Vec t;
363: requested = lme->ncv + extra;
365: /* oldsize is zero if this is the first time setup is called */
366: BVGetSizes(lme->V,NULL,NULL,&oldsize);
368: /* allocate basis vectors */
369: if (!lme->V) LMEGetBV(lme,&lme->V);
370: if (!oldsize) {
371: if (!((PetscObject)(lme->V))->type_name) BVSetType(lme->V,BVSVEC);
372: MatCreateVecsEmpty(lme->A,&t,NULL);
373: BVSetSizesFromVec(lme->V,t,requested);
374: VecDestroy(&t);
375: } else BVResize(lme->V,requested,PETSC_FALSE);
376: return 0;
377: }