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*/
 26: #include <slepcsys.h>

 30: /*
 31:    STAssociatedKSPSolve - Solves the linear system of equations associated
 32:    to the spectral transformation.

 34:    Input Parameters:
 35: .  st - the spectral transformation context
 36: .  b  - right hand side vector

 38:    Output  Parameter:
 39: .  x - computed solution
 40: */
 41: PetscErrorCode STAssociatedKSPSolve(ST st,Vec b,Vec x)
 42: {
 43:   PetscErrorCode     ierr;
 44:   PetscInt           its;
 45:   KSPConvergedReason reason;

 48:   if (!st->ksp) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST has no associated KSP");
 49:   KSPSolve(st->ksp,b,x);
 50:   KSPGetConvergedReason(st->ksp,&reason);
 51:   if (reason<0) SETERRQ1(((PetscObject)st)->comm,PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
 52:   KSPGetIterationNumber(st->ksp,&its);
 53:   st->lineariterations += its;
 54:   PetscInfo1(st,"Linear solve iterations=%D\n",its);
 55:   return(0);
 56: }

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

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

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

 78:   if (!st->ksp) SETERRQ(((PetscObject)st)->comm,PETSC_ERR_SUP,"ST has no associated KSP");
 79:   KSPSolveTranspose(st->ksp,b,x);
 80:   KSPGetConvergedReason(st->ksp,&reason);
 81:   if (reason<0) SETERRQ1(((PetscObject)st)->comm,PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
 82:   KSPGetIterationNumber(st->ksp,&its);
 83:   st->lineariterations += its;
 84:   PetscInfo1(st,"Linear solve iterations=%D\n",its);
 85:   return(0);
 86: }

 90: /*@
 91:    STSetKSP - Sets the KSP object associated with the spectral 
 92:    transformation.

 94:    Collective on ST

 96:    Input Parameters:
 97: +  st   - the spectral transformation context
 98: -  ksp  - the linear system context

100:    Level: advanced

102: @*/
103: PetscErrorCode STSetKSP(ST st,KSP ksp)
104: {

111:   PetscObjectReference((PetscObject)ksp);
112:   KSPDestroy(&st->ksp);
113:   st->ksp = ksp;
114:   PetscLogObjectParent(st,st->ksp);
115:   return(0);
116: }

120: /*@
121:    STGetKSP - Gets the KSP object associated with the spectral
122:    transformation.

124:    Not Collective

126:    Input Parameter:
127: .  st - the spectral transformation context

129:    Output Parameter:
130: .  ksp  - the linear system context

132:    Notes:
133:    On output, the value of ksp can be PETSC_NULL if the combination of 
134:    eigenproblem type and selected transformation does not require to 
135:    solve a linear system of equations.
136:    
137:    Level: intermediate

139: @*/
140: PetscErrorCode STGetKSP(ST st,KSP* ksp)
141: {

147:   if (!st->ksp) {
148:     KSPCreate(((PetscObject)st)->comm,&st->ksp);
149:     KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
150:     KSPAppendOptionsPrefix(st->ksp,"st_");
151:     PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1);
152:     PetscLogObjectParent(st,st->ksp);
153:     KSPSetTolerances(st->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
154:   }
155:   *ksp = st->ksp;
156:   return(0);
157: }

161: /*@
162:    STGetOperationCounters - Gets the total number of operator applications
163:    and linear solver iterations used by the ST object.

165:    Not Collective

167:    Input Parameter:
168: .  st - the spectral transformation context

170:    Output Parameter:
171: +  ops  - number of operator applications
172: -  lits - number of linear solver iterations

174:    Notes:
175:    Any output parameter may be PETSC_NULL on input if not needed. 
176:    
177:    Level: intermediate

179: .seealso: STResetOperationCounters()
180: @*/
181: PetscErrorCode STGetOperationCounters(ST st,PetscInt* ops,PetscInt* lits)
182: {
185:   if (ops) *ops = st->applys;
186:   if (lits) *lits = st->lineariterations;
187:   return(0);
188: }

192: /*@
193:    STResetOperationCounters - Resets the counters for operator applications,
194:    inner product operations and total number of linear iterations used by 
195:    the ST object.

197:    Logically Collective on ST

199:    Input Parameter:
200: .  st - the spectral transformation context

202:    Level: intermediate

204: .seealso: STGetOperationCounters()
205: @*/
206: PetscErrorCode STResetOperationCounters(ST st)
207: {
210:   st->lineariterations = 0;
211:   st->applys = 0;
212:   return(0);
213: }

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

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

260:    Collective on ST

262:    Input Parameters:
263: +  st - the spectral transformation context
264: .  n  - number of vectors
265: -  V  - vectors to be checked

267:    Note:
268:    This function allows to handle singular pencils and to solve some problems
269:    in which the nullspace is important (see the users guide for details).
270:    
271:    Level: developer

273: .seealso: EPSSetDeflationSpace()
274: @*/
275: PetscErrorCode STCheckNullSpace(ST st,PetscInt n,const Vec V[])
276: {

282:   if (n>0 && st->ops->checknullspace) {
285:     (*st->ops->checknullspace)(st,n,V);
286:   }
287:   return(0);
288: }