Actual source code: stsles.c
1: /*
2: The ST (spectral transformation) interface routines related to the
3: KSP object associated to it.
5: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
6: SLEPc - Scalable Library for Eigenvalue Problem Computations
7: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
9: This file is part of SLEPc.
10:
11: SLEPc is free software: you can redistribute it and/or modify it under the
12: terms of version 3 of the GNU Lesser General Public License as published by
13: the Free Software Foundation.
15: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
16: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
17: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
18: more details.
20: You should have received a copy of the GNU Lesser General Public License
21: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
22: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: */
25: #include <slepc-private/stimpl.h> /*I "slepcst.h" I*/
26: #include <slepcsys.h>
30: /*@
31: STAssociatedKSPSolve - Solves the linear system of equations associated
32: to the spectral transformation.
34: Collective on ST
36: Input Parameters:
37: . st - the spectral transformation context
38: . b - right hand side vector
40: Output Parameter:
41: . x - computed solution
43: Level: developer
45: .seealso: STAssociatedKSPSolveTranspose()
46: @*/
47: PetscErrorCode STAssociatedKSPSolve(ST st,Vec b,Vec x)
48: {
49: PetscErrorCode ierr;
50: PetscInt its;
51: KSPConvergedReason reason;
54: if (!st->ksp) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST has no associated KSP");
55: KSPSolve(st->ksp,b,x);
56: KSPGetConvergedReason(st->ksp,&reason);
57: if (reason<0) SETERRQ1(((PetscObject)st)->comm,PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
58: KSPGetIterationNumber(st->ksp,&its);
59: st->lineariterations += its;
60: PetscInfo1(st,"Linear solve iterations=%D\n",its);
61: return(0);
62: }
66: /*@
67: STAssociatedKSPSolveTranspose - Solves the transpose of the linear
68: system of equations associated to the spectral transformation.
70: Input Parameters:
71: . st - the spectral transformation context
72: . b - right hand side vector
74: Output Parameter:
75: . x - computed solution
77: Level: developer
79: .seealso: STAssociatedKSPSolve()
80: @*/
81: PetscErrorCode STAssociatedKSPSolveTranspose(ST st,Vec b,Vec x)
82: {
84: PetscInt its;
85: KSPConvergedReason reason;
88: if (!st->ksp) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST has no associated KSP");
89: KSPSolveTranspose(st->ksp,b,x);
90: KSPGetConvergedReason(st->ksp,&reason);
91: if (reason<0) SETERRQ1(((PetscObject)st)->comm,PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
92: KSPGetIterationNumber(st->ksp,&its);
93: st->lineariterations += its;
94: PetscInfo1(st,"Linear solve iterations=%D\n",its);
95: return(0);
96: }
100: /*
101: STMatSetHermitian - Sets the Hermitian flag to the ST matrix.
103: Input Parameters:
104: . st - the spectral transformation context
105: . M - matrix
106: */
107: PetscErrorCode STMatSetHermitian(ST st,Mat M)
108: {
109: #if defined(PETSC_USE_COMPLEX)
111: PetscBool set,aherm,bherm,mherm;
112: #endif
115: #if defined(PETSC_USE_COMPLEX)
116: MatIsHermitianKnown(st->A,&set,&aherm);
117: if (!set) aherm = PETSC_FALSE;
118: mherm = aherm;
119: if (st->B) {
120: MatIsHermitianKnown(st->B,&set,&bherm);
121: if (!set) bherm = PETSC_FALSE;
122: mherm = (mherm && bherm)? PETSC_TRUE: PETSC_FALSE;
123: }
124: mherm = (mherm && PetscImaginaryPart(st->sigma)==0.0)? PETSC_TRUE: PETSC_FALSE;
125: MatSetOption(M,MAT_HERMITIAN,mherm);
126: #endif
127: return(0);
128: }
133: /*@
134: STSetKSP - Sets the KSP object associated with the spectral
135: transformation.
137: Collective on ST
139: Input Parameters:
140: + st - the spectral transformation context
141: - ksp - the linear system context
143: Level: advanced
145: @*/
146: PetscErrorCode STSetKSP(ST st,KSP ksp)
147: {
154: PetscObjectReference((PetscObject)ksp);
155: KSPDestroy(&st->ksp);
156: st->ksp = ksp;
157: PetscLogObjectParent(st,st->ksp);
158: return(0);
159: }
163: /*@
164: STGetKSP - Gets the KSP object associated with the spectral
165: transformation.
167: Not Collective
169: Input Parameter:
170: . st - the spectral transformation context
172: Output Parameter:
173: . ksp - the linear system context
175: Notes:
176: On output, the value of ksp can be PETSC_NULL if the combination of
177: eigenproblem type and selected transformation does not require to
178: solve a linear system of equations.
179:
180: Level: intermediate
182: @*/
183: PetscErrorCode STGetKSP(ST st,KSP* ksp)
184: {
190: if (!st->ksp) {
191: KSPCreate(((PetscObject)st)->comm,&st->ksp);
192: KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
193: KSPAppendOptionsPrefix(st->ksp,"st_");
194: PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1);
195: PetscLogObjectParent(st,st->ksp);
196: KSPSetTolerances(st->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
197: }
198: *ksp = st->ksp;
199: return(0);
200: }
204: /*@
205: STGetOperationCounters - Gets the total number of operator applications
206: and linear solver iterations used by the ST object.
208: Not Collective
210: Input Parameter:
211: . st - the spectral transformation context
213: Output Parameter:
214: + ops - number of operator applications
215: - lits - number of linear solver iterations
217: Notes:
218: Any output parameter may be PETSC_NULL on input if not needed.
219:
220: Level: intermediate
222: .seealso: STResetOperationCounters()
223: @*/
224: PetscErrorCode STGetOperationCounters(ST st,PetscInt* ops,PetscInt* lits)
225: {
228: if (ops) *ops = st->applys;
229: if (lits) *lits = st->lineariterations;
230: return(0);
231: }
235: /*@
236: STResetOperationCounters - Resets the counters for operator applications,
237: inner product operations and total number of linear iterations used by
238: the ST object.
240: Logically Collective on ST
242: Input Parameter:
243: . st - the spectral transformation context
245: Level: intermediate
247: .seealso: STGetOperationCounters()
248: @*/
249: PetscErrorCode STResetOperationCounters(ST st)
250: {
253: st->lineariterations = 0;
254: st->applys = 0;
255: return(0);
256: }
260: PetscErrorCode STCheckNullSpace_Default(ST st,PetscInt n,const Vec V[])
261: {
263: PetscInt i,c;
264: PetscReal norm;
265: Vec *T,w;
266: Mat A;
267: PC pc;
268: MatNullSpace nullsp;
269:
271: PetscMalloc(n*sizeof(Vec),&T);
272: if (!st->ksp) { STGetKSP(st,&st->ksp); }
273: KSPGetPC(st->ksp,&pc);
274: PCGetOperators(pc,&A,PETSC_NULL,PETSC_NULL);
275: MatGetVecs(A,PETSC_NULL,&w);
276: c = 0;
277: for (i=0;i<n;i++) {
278: MatMult(A,V[i],w);
279: VecNorm(w,NORM_2,&norm);
280: if (norm < 1e-8) {
281: PetscInfo2(st,"Vector %D norm=%g\n",i,(double)norm);
282: T[c] = V[i];
283: c++;
284: }
285: }
286: VecDestroy(&w);
287: if (c>0) {
288: MatNullSpaceCreate(((PetscObject)st)->comm,PETSC_FALSE,c,T,&nullsp);
289: KSPSetNullSpace(st->ksp,nullsp);
290: MatNullSpaceDestroy(&nullsp);
291: }
292: PetscFree(T);
293: return(0);
294: }
298: /*@
299: STCheckNullSpace - Given a set of vectors, this function tests each of
300: them to be a nullspace vector of the coefficient matrix of the associated
301: KSP object. All these nullspace vectors are passed to the KSP object.
303: Collective on ST
305: Input Parameters:
306: + st - the spectral transformation context
307: . n - number of vectors
308: - V - vectors to be checked
310: Note:
311: This function allows to handle singular pencils and to solve some problems
312: in which the nullspace is important (see the users guide for details).
313:
314: Level: developer
316: .seealso: EPSSetDeflationSpace()
317: @*/
318: PetscErrorCode STCheckNullSpace(ST st,PetscInt n,const Vec V[])
319: {
325: if (n>0 && st->ops->checknullspace) {
328: (*st->ops->checknullspace)(st,n,V);
329: }
330: return(0);
331: }