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-2011, 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 <private/stimpl.h> /*I "slepcst.h" I*/
29: /*
30: STAssociatedKSPSolve - Solves the linear system of equations associated
31: to the spectral transformation.
33: Input Parameters:
34: . st - the spectral transformation context
35: . b - right hand side vector
37: Output Parameter:
38: . x - computed solution
39: */
40: PetscErrorCode STAssociatedKSPSolve(ST st,Vec b,Vec x)
41: {
42: PetscErrorCode ierr;
43: PetscInt its;
44: KSPConvergedReason reason;
47: if (!st->ksp) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST has no associated KSP");
48: KSPSolve(st->ksp,b,x);
49: KSPGetConvergedReason(st->ksp,&reason);
50: if (reason<0) SETERRQ1(((PetscObject)st)->comm,PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
51: KSPGetIterationNumber(st->ksp,&its);
52: st->lineariterations += its;
53: PetscInfo1(st,"Linear solve iterations=%D\n",its);
54: return(0);
55: }
59: /*
60: STAssociatedKSPSolveTranspose - Solves the transpose of the linear
61: system of equations associated to the spectral transformation.
63: Input Parameters:
64: . st - the spectral transformation context
65: . b - right hand side vector
67: Output Parameter:
68: . x - computed solution
69: */
70: PetscErrorCode STAssociatedKSPSolveTranspose(ST st,Vec b,Vec x)
71: {
73: PetscInt its;
74: KSPConvergedReason reason;
77: if (!st->ksp) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST has no associated KSP");
78: KSPSolveTranspose(st->ksp,b,x);
79: KSPGetConvergedReason(st->ksp,&reason);
80: if (reason<0) SETERRQ1(((PetscObject)st)->comm,PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
81: KSPGetIterationNumber(st->ksp,&its);
82: st->lineariterations += its;
83: PetscInfo1(st,"Linear solve iterations=%D\n",its);
84: return(0);
85: }
89: /*@
90: STSetKSP - Sets the KSP object associated with the spectral
91: transformation.
93: Collective on ST
95: Input Parameters:
96: + st - the spectral transformation context
97: - ksp - the linear system context
99: Level: advanced
101: @*/
102: PetscErrorCode STSetKSP(ST st,KSP ksp)
103: {
110: PetscObjectReference((PetscObject)ksp);
111: KSPDestroy(&st->ksp);
112: st->ksp = ksp;
113: PetscLogObjectParent(st,st->ksp);
114: return(0);
115: }
119: /*@
120: STGetKSP - Gets the KSP object associated with the spectral
121: transformation.
123: Not Collective
125: Input Parameter:
126: . st - the spectral transformation context
128: Output Parameter:
129: . ksp - the linear system context
131: Notes:
132: On output, the value of ksp can be PETSC_NULL if the combination of
133: eigenproblem type and selected transformation does not require to
134: solve a linear system of equations.
135:
136: Level: intermediate
138: @*/
139: PetscErrorCode STGetKSP(ST st,KSP* ksp)
140: {
146: if (!st->ksp) {
147: KSPCreate(((PetscObject)st)->comm,&st->ksp);
148: KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
149: KSPAppendOptionsPrefix(st->ksp,"st_");
150: PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1);
151: PetscLogObjectParent(st,st->ksp);
152: }
153: *ksp = st->ksp;
154: return(0);
155: }
159: /*@
160: STGetOperationCounters - Gets the total number of operator applications
161: and linear solver iterations used by the ST object.
163: Not Collective
165: Input Parameter:
166: . st - the spectral transformation context
168: Output Parameter:
169: + ops - number of operator applications
170: - lits - number of linear solver iterations
172: Notes:
173: Any output parameter may be PETSC_NULL on input if not needed.
174:
175: Level: intermediate
177: .seealso: STResetOperationCounters()
178: @*/
179: PetscErrorCode STGetOperationCounters(ST st,PetscInt* ops,PetscInt* lits)
180: {
183: if (ops) *ops = st->applys;
184: if (lits) *lits = st->lineariterations;
185: return(0);
186: }
190: /*@
191: STResetOperationCounters - Resets the counters for operator applications,
192: inner product operations and total number of linear iterations used by
193: the ST object.
195: Logically Collective on ST
197: Input Parameter:
198: . st - the spectral transformation context
200: Level: intermediate
202: .seealso: STGetOperationCounters()
203: @*/
204: PetscErrorCode STResetOperationCounters(ST st)
205: {
208: st->lineariterations = 0;
209: st->applys = 0;
210: return(0);
211: }
215: PetscErrorCode STCheckNullSpace_Default(ST st,PetscInt n,const Vec V[])
216: {
218: PetscInt i,c;
219: PetscReal norm;
220: Vec *T,w;
221: Mat A;
222: PC pc;
223: MatNullSpace nullsp;
224:
226: PetscMalloc(n*sizeof(Vec),&T);
227: if (!st->ksp) { STGetKSP(st,&st->ksp); }
228: KSPGetPC(st->ksp,&pc);
229: PCGetOperators(pc,&A,PETSC_NULL,PETSC_NULL);
230: MatGetVecs(A,PETSC_NULL,&w);
231: c = 0;
232: for (i=0;i<n;i++) {
233: MatMult(A,V[i],w);
234: VecNorm(w,NORM_2,&norm);
235: if (norm < 1e-8) {
236: PetscInfo2(st,"Vector %D norm=%g\n",i,(double)norm);
237: T[c] = V[i];
238: c++;
239: }
240: }
241: VecDestroy(&w);
242: if (c>0) {
243: MatNullSpaceCreate(((PetscObject)st)->comm,PETSC_FALSE,c,T,&nullsp);
244: KSPSetNullSpace(st->ksp,nullsp);
245: MatNullSpaceDestroy(&nullsp);
246: }
247: PetscFree(T);
248: return(0);
249: }
253: /*@
254: STCheckNullSpace - Given a set of vectors, this function tests each of
255: them to be a nullspace vector of the coefficient matrix of the associated
256: KSP object. All these nullspace vectors are passed to the KSP object.
258: Collective on ST
260: Input Parameters:
261: + st - the spectral transformation context
262: . n - number of vectors
263: - V - vectors to be checked
265: Note:
266: This function allows to handle singular pencils and to solve some problems
267: in which the nullspace is important (see the users guide for details).
268:
269: Level: developer
271: .seealso: EPSSetDeflationSpace()
272: @*/
273: PetscErrorCode STCheckNullSpace(ST st,PetscInt n,const Vec V[])
274: {
280: if (n>0 && st->ops->checknullspace) {
283: (*st->ops->checknullspace)(st,n,V);
284: }
285: return(0);
286: }