1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: ST interface routines related to the KSP object associated with it
12: */
14: #include <slepc/private/stimpl.h> /*I "slepcst.h" I*/
16: /*
17: This is used to set a default type for the KSP and PC objects.
18: It is called at STSetFromOptions (before KSPSetFromOptions)
19: and also at STSetUp (in case STSetFromOptions was not called).
20: */
21: PetscErrorCode STSetDefaultKSP(ST st) 22: {
28: if (!st->ksp) { STGetKSP(st,&st->ksp); }
29: if (st->ops->setdefaultksp) { (*st->ops->setdefaultksp)(st); }
30: return(0);
31: }
33: /*
34: Checks whether the ST matrix is symmetric
35: */
36: PETSC_STATIC_INLINE PetscErrorCode STMatIsSymmetricKnown(ST st,PetscBool *symm) 37: {
39: PetscInt i;
40: PetscBool set,sbaij=PETSC_FALSE,asymm;
43: *symm = PETSC_FALSE;
44: if (!st->nmat) return(0); /* STSetMatrices() not called yet */
45: /* check if problem matrices are all sbaij */
46: for (i=0;i<st->nmat;i++) {
47: PetscObjectTypeCompareAny((PetscObject)st->A[i],&sbaij,MATSEQSBAIJ,MATMPISBAIJ,"");
48: if (!sbaij) break;
49: }
50: if (sbaij) *symm = PETSC_TRUE;
51: else {
52: /* check if user has set the symmetric flag */
53: for (i=0;i<st->nmat;i++) {
54: MatIsSymmetricKnown(st->A[i],&set,&asymm);
55: if (!set) asymm = PETSC_FALSE;
56: if (!asymm) return(0);
57: if (PetscRealPart(st->sigma)==0.0) break;
58: }
59: *symm = PETSC_TRUE;
60: }
61: return(0);
62: }
64: /*
65: This is done by all ST types except PRECOND.
66: The default is an LU direct solver, or GMRES+Jacobi if matmode=shell.
67: */
68: PetscErrorCode STSetDefaultKSP_Default(ST st) 69: {
71: PC pc;
72: PCType pctype;
73: KSPType ksptype;
74: PetscBool asymm;
77: KSPGetPC(st->ksp,&pc);
78: KSPGetType(st->ksp,&ksptype);
79: PCGetType(pc,&pctype);
80: if (!pctype && !ksptype) {
81: if (st->matmode == ST_MATMODE_SHELL) {
82: KSPSetType(st->ksp,KSPGMRES);
83: PCSetType(pc,PCJACOBI);
84: } else {
85: STMatIsSymmetricKnown(st,&asymm);
86: KSPSetType(st->ksp,KSPPREONLY);
87: PCSetType(pc,asymm?PCCHOLESKY:PCLU);
88: }
89: }
90: KSPSetErrorIfNotConverged(st->ksp,PETSC_TRUE);
91: return(0);
92: }
94: /*@
95: STMatMult - Computes the matrix-vector product y = T[k] x, where T[k] is
96: the k-th matrix of the spectral transformation.
98: Neighbor-wise Collective on st
100: Input Parameters:
101: + st - the spectral transformation context
102: . k - index of matrix to use
103: - x - the vector to be multiplied
105: Output Parameter:
106: . y - the result
108: Level: developer
110: .seealso: STMatMultTranspose()
111: @*/
112: PetscErrorCode STMatMult(ST st,PetscInt k,Vec x,Vec y)113: {
121: STCheckMatrices(st,1);
122: if (k<0 || k>=PetscMax(2,st->nmat)) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat);
123: if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
124: VecSetErrorIfLocked(y,3);
126: if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
127: VecLockReadPush(x);
128: PetscLogEventBegin(ST_MatMult,st,x,y,0);
129: if (!st->T[k]) {
130: /* T[k]=NULL means identity matrix */
131: VecCopy(x,y);
132: } else {
133: MatMult(st->T[k],x,y);
134: }
135: PetscLogEventEnd(ST_MatMult,st,x,y,0);
136: VecLockReadPop(x);
137: return(0);
138: }
140: /*@
141: STMatMultTranspose - Computes the matrix-vector product y = T[k]' x, where T[k] is
142: the k-th matrix of the spectral transformation.
144: Neighbor-wise Collective on st
146: Input Parameters:
147: + st - the spectral transformation context
148: . k - index of matrix to use
149: - x - the vector to be multiplied
151: Output Parameter:
152: . y - the result
154: Level: developer
156: .seealso: STMatMult()
157: @*/
158: PetscErrorCode STMatMultTranspose(ST st,PetscInt k,Vec x,Vec y)159: {
167: STCheckMatrices(st,1);
168: if (k<0 || k>=PetscMax(2,st->nmat)) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat);
169: if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
170: VecSetErrorIfLocked(y,3);
172: if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
173: VecLockReadPush(x);
174: PetscLogEventBegin(ST_MatMultTranspose,st,x,y,0);
175: if (!st->T[k]) {
176: /* T[k]=NULL means identity matrix */
177: VecCopy(x,y);
178: } else {
179: MatMultTranspose(st->T[k],x,y);
180: }
181: PetscLogEventEnd(ST_MatMultTranspose,st,x,y,0);
182: VecLockReadPop(x);
183: return(0);
184: }
186: /*@
187: STMatSolve - Solves P x = b, where P is the preconditioner matrix of
188: the spectral transformation, using a KSP object stored internally.
190: Collective on st
192: Input Parameters:
193: + st - the spectral transformation context
194: - b - right hand side vector
196: Output Parameter:
197: . x - computed solution
199: Level: developer
201: .seealso: STMatSolveTranspose()
202: @*/
203: PetscErrorCode STMatSolve(ST st,Vec b,Vec x)204: {
206: PetscInt its;
207: PetscBool flg;
213: STCheckMatrices(st,1);
214: if (x == b) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
215: VecSetErrorIfLocked(x,3);
217: if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
218: VecLockReadPush(b);
219: PetscLogEventBegin(ST_MatSolve,st,b,x,0);
220: PetscObjectTypeCompareAny((PetscObject)st,&flg,STPRECOND,STSHELL,"");
221: if (!flg && !st->P) {
222: /* P=NULL means identity matrix */
223: VecCopy(b,x);
224: } else {
225: if (!st->ksp) { STGetKSP(st,&st->ksp); }
226: KSPSolve(st->ksp,b,x);
227: KSPGetIterationNumber(st->ksp,&its);
228: PetscInfo1(st,"Linear solve iterations=%D\n",its);
229: }
230: PetscLogEventEnd(ST_MatSolve,st,b,x,0);
231: VecLockReadPop(b);
232: return(0);
233: }
235: /*@
236: STMatSolveTranspose - Solves P' x = b, where P is the preconditioner matrix of
237: the spectral transformation, using a KSP object stored internally.
239: Collective on st
241: Input Parameters:
242: . st - the spectral transformation context
243: . b - right hand side vector
245: Output Parameter:
246: . x - computed solution
248: Level: developer
250: .seealso: STMatSolve()
251: @*/
252: PetscErrorCode STMatSolveTranspose(ST st,Vec b,Vec x)253: {
255: PetscInt its;
256: PetscBool flg;
262: STCheckMatrices(st,1);
263: if (x == b) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
264: VecSetErrorIfLocked(x,3);
266: if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
267: VecLockReadPush(b);
268: PetscLogEventBegin(ST_MatSolveTranspose,st,b,x,0);
269: PetscObjectTypeCompareAny((PetscObject)st,&flg,STPRECOND,STSHELL,"");
270: if (!flg && !st->P) {
271: /* P=NULL means identity matrix */
272: VecCopy(b,x);
273: } else {
274: if (!st->ksp) { STGetKSP(st,&st->ksp); }
275: KSPSolveTranspose(st->ksp,b,x);
276: KSPGetIterationNumber(st->ksp,&its);
277: PetscInfo1(st,"Linear solve iterations=%D\n",its);
278: }
279: PetscLogEventEnd(ST_MatSolveTranspose,st,b,x,0);
280: VecLockReadPop(b);
281: return(0);
282: }
284: /*
285: STMatSetHermitian - Sets the symmetric and/or Hermitian flag to the ST matrix.
287: Input Parameters:
288: . st - the spectral transformation context
289: . M - matrix
290: */
291: PetscErrorCode STMatSetHermitian(ST st,Mat M)292: {
294: PetscBool set,aherm,mherm;
295: PetscInt i;
298: if (!st->nmat) return(0); /* STSetMatrices() not called yet */
299: mherm = PETSC_TRUE;
300: for (i=0;i<st->nmat;i++) {
301: MatIsSymmetricKnown(st->A[i],&set,&aherm);
302: if (!set) aherm = PETSC_FALSE;
303: if (!aherm) { mherm = PETSC_FALSE; break; }
304: if (PetscRealPart(st->sigma)==0.0) break;
305: }
306: MatSetOption(M,MAT_SYMMETRIC,mherm);
307: #if defined(PETSC_USE_COMPLEX)
308: if (PetscImaginaryPart(st->sigma)==0.0) {
309: mherm = PETSC_TRUE;
310: for (i=0;i<st->nmat;i++) {
311: MatIsHermitianKnown(st->A[i],&set,&aherm);
312: if (!set) aherm = PETSC_FALSE;
313: if (!aherm) { mherm = PETSC_FALSE; break; }
314: if (PetscRealPart(st->sigma)==0.0) break;
315: }
316: MatSetOption(M,MAT_HERMITIAN,mherm);
317: }
318: #endif
319: return(0);
320: }
322: PetscErrorCode STCheckFactorPackage(ST st)323: {
325: PC pc;
326: PetscMPIInt size;
327: PetscBool flg;
328: MatSolverType stype;
331: MPI_Comm_size(PetscObjectComm((PetscObject)st),&size);
332: if (size==1) return(0);
333: KSPGetPC(st->ksp,&pc);
334: PCFactorGetMatSolverType(pc,&stype);
335: if (stype) { /* currently selected PC is a factorization */
336: PetscStrcmp(stype,MATSOLVERPETSC,&flg);
337: if (flg) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"You chose to solve linear systems with a factorization, but in parallel runs you need to select an external package; see the users guide for details");
338: }
339: return(0);
340: }
342: /*@
343: STSetKSP - Sets the KSP object associated with the spectral
344: transformation.
346: Collective on st
348: Input Parameters:
349: + st - the spectral transformation context
350: - ksp - the linear system context
352: Level: advanced
353: @*/
354: PetscErrorCode STSetKSP(ST st,KSP ksp)355: {
362: PetscObjectReference((PetscObject)ksp);
363: KSPDestroy(&st->ksp);
364: st->ksp = ksp;
365: PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
366: return(0);
367: }
369: /*@
370: STGetKSP - Gets the KSP object associated with the spectral
371: transformation.
373: Not Collective
375: Input Parameter:
376: . st - the spectral transformation context
378: Output Parameter:
379: . ksp - the linear system context
381: Level: intermediate
382: @*/
383: PetscErrorCode STGetKSP(ST st,KSP* ksp)384: {
390: if (!st->ksp) {
391: KSPCreate(PetscObjectComm((PetscObject)st),&st->ksp);
392: PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1);
393: KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
394: KSPAppendOptionsPrefix(st->ksp,"st_");
395: PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
396: PetscObjectSetOptions((PetscObject)st->ksp,((PetscObject)st)->options);
397: KSPSetTolerances(st->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
398: }
399: *ksp = st->ksp;
400: return(0);
401: }
403: PetscErrorCode STCheckNullSpace_Default(ST st,BV V)404: {
406: PetscInt nc,i,c;
407: PetscReal norm;
408: Vec *T,w,vi;
409: Mat A;
410: PC pc;
411: MatNullSpace nullsp;
414: BVGetNumConstraints(V,&nc);
415: PetscMalloc1(nc,&T);
416: if (!st->ksp) { STGetKSP(st,&st->ksp); }
417: KSPGetPC(st->ksp,&pc);
418: PCGetOperators(pc,&A,NULL);
419: MatCreateVecs(A,NULL,&w);
420: c = 0;
421: for (i=0;i<nc;i++) {
422: BVGetColumn(V,-nc+i,&vi);
423: MatMult(A,vi,w);
424: VecNorm(w,NORM_2,&norm);
425: if (norm < 1e-8) {
426: PetscInfo2(st,"Vector %D norm=%g\n",i,(double)norm);
427: BVCreateVec(V,T+c);
428: VecCopy(vi,T[c]);
429: c++;
430: }
431: BVRestoreColumn(V,-nc+i,&vi);
432: }
433: VecDestroy(&w);
434: if (c>0) {
435: MatNullSpaceCreate(PetscObjectComm((PetscObject)st),PETSC_FALSE,c,T,&nullsp);
436: MatSetNullSpace(A,nullsp);
437: MatNullSpaceDestroy(&nullsp);
438: VecDestroyVecs(c,&T);
439: } else {
440: PetscFree(T);
441: }
442: return(0);
443: }
445: /*@
446: STCheckNullSpace - Given a basis vectors object, this function tests each
447: of its constraint vectors to be a nullspace vector of the coefficient
448: matrix of the associated KSP object. All these nullspace vectors are passed
449: to the KSP object.
451: Collective on st
453: Input Parameters:
454: + st - the spectral transformation context
455: - V - basis vectors to be checked
457: Note:
458: This function allows to handle singular pencils and to solve some problems
459: in which the nullspace is important (see the users guide for details).
461: Level: developer
463: .seealso: EPSSetDeflationSpace()
464: @*/
465: PetscErrorCode STCheckNullSpace(ST st,BV V)466: {
468: PetscInt nc;
475: if (!st->state) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Must call STSolve() first");
477: BVGetNumConstraints(V,&nc);
478: if (nc && st->ops->checknullspace) {
479: (*st->ops->checknullspace)(st,V);
480: }
481: return(0);
482: }