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-2010, Universidad 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

 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: {
 43:   PetscInt       its;
 44:   KSPConvergedReason reason;

 50:   if (!st->ksp) { SETERRQ(PETSC_ERR_SUP,"ST has no associated KSP"); }
 51:   KSPSolve(st->ksp,b,x);
 52:   KSPGetConvergedReason(st->ksp,&reason);
 53:   if (reason<0) { SETERRQ1(0,"Warning: KSP did not converge (%d)",reason); }
 54:   KSPGetIterationNumber(st->ksp,&its);
 55:   st->lineariterations += its;
 56:   PetscInfo1(st,"Linear solve iterations=%d\n",its);
 57:   return(0);
 58: }

 62: /*
 63:    STAssociatedKSPSolveTranspose - Solves the transpose of the linear 
 64:    system of equations associated to the spectral transformation.

 66:    Input Parameters:
 67: .  st - the spectral transformation context
 68: .  b  - right hand side vector

 70:    Output  Parameter:
 71: .  x - computed solution
 72: */
 73: PetscErrorCode STAssociatedKSPSolveTranspose(ST st,Vec b,Vec x)
 74: {
 76:   PetscInt       its;
 77:   KSPConvergedReason reason;

 83:   if (!st->ksp) { SETERRQ(PETSC_ERR_SUP,"ST has no associated KSP"); }
 84:   KSPSolveTranspose(st->ksp,b,x);
 85:   KSPGetConvergedReason(st->ksp,&reason);
 86:   if (reason<0) { SETERRQ1(0,"Warning: KSP did not converge (%d)",reason); }
 87:   KSPGetIterationNumber(st->ksp,&its);
 88:   st->lineariterations += its;
 89:   PetscInfo1(st,"Linear solve iterations=%d\n",its);
 90:   return(0);
 91: }

 95: /*@
 96:    STSetKSP - Sets the KSP object associated with the spectral 
 97:    transformation.

 99:    Not collective

101:    Input Parameters:
102: +  st   - the spectral transformation context
103: -  ksp  - the linear system context

105:    Level: advanced

107: @*/
108: PetscErrorCode STSetKSP(ST st,KSP ksp)
109: {

116:   PetscObjectReference((PetscObject)ksp);
117:   if (st->ksp) {
118:     KSPDestroy(st->ksp);
119:   }
120:   st->ksp = ksp;
121:   return(0);
122: }

126: /*@
127:    STGetKSP - Gets the KSP object associated with the spectral
128:    transformation.

130:    Not collective

132:    Input Parameter:
133: .  st - the spectral transformation context

135:    Output Parameter:
136: .  ksp  - the linear system context

138:    Notes:
139:    On output, the value of ksp can be PETSC_NULL if the combination of 
140:    eigenproblem type and selected transformation does not require to 
141:    solve a linear system of equations.
142:    
143:    Level: intermediate

145: @*/
146: PetscErrorCode STGetKSP(ST st,KSP* ksp)
147: {
150:   if (!((PetscObject)st)->type_name) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call STSetType first"); }
151:   if (ksp)  *ksp = st->ksp;
152:   return(0);
153: }

157: /*@
158:    STGetOperationCounters - Gets the total number of operator applications
159:    and linear solver iterations used by the ST object.

161:    Not Collective

163:    Input Parameter:
164: .  st - the spectral transformation context

166:    Output Parameter:
167: +  ops  - number of operator applications
168: -  lits - number of linear solver iterations

170:    Notes:
171:    Any output parameter may be PETSC_NULL on input if not needed. 
172:    
173:    Level: intermediate

175: .seealso: STResetOperationCounters()
176: @*/
177: PetscErrorCode STGetOperationCounters(ST st,PetscInt* ops,PetscInt* lits)
178: {
181:   if (ops) *ops = st->applys;
182:   if (lits) *lits = st->lineariterations;
183:   return(0);
184: }

188: /*@
189:    STResetOperationCounters - Resets the counters for operator applications,
190:    inner product operations and total number of linear iterations used by 
191:    the ST object.

193:    Collective on ST

195:    Input Parameter:
196: .  st - the spectral transformation context

198:    Level: intermediate

200: .seealso: STGetOperationCounters()
201: @*/
202: PetscErrorCode STResetOperationCounters(ST st)
203: {
206:   st->lineariterations = 0;
207:   st->applys = 0;
208:   return(0);
209: }

213: PetscErrorCode STCheckNullSpace_Default(ST st,PetscInt n,const Vec V[])
214: {
216:   PetscInt       i,c;
217:   PetscReal      norm;
218:   Vec            *T,w;
219:   Mat            A;
220:   PC             pc;
221:   MatNullSpace   nullsp;
222: 
224:   PetscMalloc(n*sizeof(Vec),&T);
225:   KSPGetPC(st->ksp,&pc);
226:   PCGetOperators(pc,&A,PETSC_NULL,PETSC_NULL);
227:   MatGetVecs(A,PETSC_NULL,&w);
228:   c = 0;
229:   for (i=0;i<n;i++) {
230:     MatMult(A,V[i],w);
231:     VecNorm(w,NORM_2,&norm);
232:     if (norm < 1e-8) {
233:       PetscInfo2(st,"Vector %i norm=%g\n",i,norm);
234:       T[c] = V[i];
235:       c++;
236:     }
237:   }
238:   VecDestroy(w);
239:   if (c>0) {
240:     MatNullSpaceCreate(((PetscObject)st)->comm,PETSC_FALSE,c,T,&nullsp);
241:     KSPSetNullSpace(st->ksp,nullsp);
242:     MatNullSpaceDestroy(nullsp);
243:   }
244:   PetscFree(T);
245:   return(0);
246: }

250: /*@
251:    STCheckNullSpace - Given a set of vectors, this function tests each of
252:    them to be a nullspace vector of the coefficient matrix of the associated
253:    KSP object. All these nullspace vectors are passed to the KSP object.

255:    Collective on ST

257:    Input Parameters:
258: +  st - the spectral transformation context
259: .  n  - number of vectors
260: -  V  - vectors to be checked

262:    Note:
263:    This function allows to handle singular pencils and to solve some problems
264:    in which the nullspace is important (see the users guide for details).
265:    
266:    Level: developer

268: .seealso: EPSSetDeflationSpace()
269: @*/
270: PetscErrorCode STCheckNullSpace(ST st,PetscInt n,const Vec V[])
271: {

275:   if (n>0 && st->checknullspace) {
276:     (*st->checknullspace)(st,n,V);
277:   }
278:   return(0);
279: }