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-2007, Universidad Politecnica de Valencia, Spain

  9:       This file is part of SLEPc. See the README file for conditions of use
 10:       and additional information.
 11:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 12: */

 14:  #include src/st/stimpl.h

 18: /*
 19:    STAssociatedKSPSolve - Solves the linear system of equations associated
 20:    to the spectral transformation.

 22:    Input Parameters:
 23: .  st - the spectral transformation context
 24: .  b  - right hand side vector

 26:    Output  Parameter:
 27: .  x - computed solution
 28: */
 29: PetscErrorCode STAssociatedKSPSolve(ST st,Vec b,Vec x)
 30: {
 32:   PetscInt       its;
 33:   KSPConvergedReason reason;

 39:   if (!st->ksp) { SETERRQ(PETSC_ERR_SUP,"ST has no associated KSP"); }
 40:   KSPSolve(st->ksp,b,x);
 41:   KSPGetConvergedReason(st->ksp,&reason);
 42:   if (reason<0) { SETERRQ1(0,"Warning: KSP did not converge (%d)",reason); }
 43:   KSPGetIterationNumber(st->ksp,&its);
 44:   st->lineariterations += its;
 45:   PetscInfo1(st,"Linear solve iterations=%d\n",its);
 46:   return(0);
 47: }

 51: /*
 52:    STAssociatedKSPSolveTranspose - Solves the transpose of the linear 
 53:    system of equations associated to the spectral transformation.

 55:    Input Parameters:
 56: .  st - the spectral transformation context
 57: .  b  - right hand side vector

 59:    Output  Parameter:
 60: .  x - computed solution
 61: */
 62: PetscErrorCode STAssociatedKSPSolveTranspose(ST st,Vec b,Vec x)
 63: {
 65:   PetscInt       its;
 66:   KSPConvergedReason reason;

 72:   if (!st->ksp) { SETERRQ(PETSC_ERR_SUP,"ST has no associated KSP"); }
 73:   KSPSolveTranspose(st->ksp,b,x);
 74:   KSPGetConvergedReason(st->ksp,&reason);
 75:   if (reason<0) { SETERRQ1(0,"Warning: KSP did not converge (%d)",reason); }
 76:   KSPGetIterationNumber(st->ksp,&its);
 77:   st->lineariterations += its;
 78:   PetscInfo1(st,"Linear solve iterations=%d\n",its);
 79:   return(0);
 80: }

 84: /*@
 85:    STSetKSP - Sets the KSP object associated with the spectral 
 86:    transformation.

 88:    Not collective

 90:    Input Parameters:
 91: +  st   - the spectral transformation context
 92: -  ksp  - the linear system context

 94:    Level: advanced

 96: @*/
 97: PetscErrorCode STSetKSP(ST st,KSP ksp)
 98: {

105:   PetscObjectReference((PetscObject)ksp);
106:   if (st->ksp) {
107:     KSPDestroy(st->ksp);
108:   }
109:   st->ksp = ksp;
110:   return(0);
111: }

115: /*@
116:    STGetKSP - Gets the KSP object associated with the spectral
117:    transformation.

119:    Not collective

121:    Input Parameter:
122: .  st - the spectral transformation context

124:    Output Parameter:
125: .  ksp  - the linear system context

127:    Notes:
128:    On output, the value of ksp can be PETSC_NULL if the combination of 
129:    eigenproblem type and selected transformation does not require to 
130:    solve a linear system of equations.
131:    
132:    Level: intermediate

134: @*/
135: PetscErrorCode STGetKSP(ST st,KSP* ksp)
136: {
139:   if (!st->type_name) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call STSetType first"); }
140:   if (ksp)  *ksp = st->ksp;
141:   return(0);
142: }

146: /*@
147:    STGetOperationCounters - Gets the total number of operator applications
148:    and linear solver iterations used by the ST object.

150:    Not Collective

152:    Input Parameter:
153: .  st - the spectral transformation context

155:    Output Parameter:
156: +  ops  - number of operator applications
157: -  lits - number of linear solver iterations

159:    Notes:
160:    Any output parameter may be PETSC_NULL on input if not needed. 
161:    
162:    Level: intermediate

164: .seealso: STResetOperationCounters()
165: @*/
166: PetscErrorCode STGetOperationCounters(ST st,int* ops,int* lits)
167: {
170:   if (ops) *ops = st->applys;
171:   if (lits) *lits = st->lineariterations;
172:   return(0);
173: }

177: /*@
178:    STResetOperationCounters - Resets the counters for operator applications,
179:    inner product operations and total number of linear iterations used by 
180:    the ST object.

182:    Collective on ST

184:    Input Parameter:
185: .  st - the spectral transformation context

187:    Level: intermediate

189: .seealso: STGetOperationCounters()
190: @*/
191: PetscErrorCode STResetOperationCounters(ST st)
192: {
195:   st->lineariterations = 0;
196:   st->applys = 0;
197:   return(0);
198: }

202: PetscErrorCode STCheckNullSpace_Default(ST st,int n,const Vec V[])
203: {
205:   int            i,c;
206:   PetscReal      norm;
207:   Vec            *T,w;
208:   Mat            A;
209:   PC             pc;
210:   MatNullSpace   nullsp;
211: 
213:   PetscMalloc(n*sizeof(Vec),&T);
214:   KSPGetPC(st->ksp,&pc);
215:   PCGetOperators(pc,&A,PETSC_NULL,PETSC_NULL);
216:   MatGetVecs(A,PETSC_NULL,&w);
217:   c = 0;
218:   for (i=0;i<n;i++) {
219:     MatMult(A,V[i],w);
220:     VecNorm(w,NORM_2,&norm);
221:     if (norm < 1e-8) {
222:       PetscInfo2(st,"Vector %i norm=%g\n",i,norm);
223:       T[c] = V[i];
224:       c++;
225:     }
226:   }
227:   VecDestroy(w);
228:   if (c>0) {
229:     MatNullSpaceCreate(st->comm,PETSC_FALSE,c,T,&nullsp);
230:     KSPSetNullSpace(st->ksp,nullsp);
231:     MatNullSpaceDestroy(nullsp);
232:   }
233:   PetscFree(T);
234:   return(0);
235: }

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

244:    Collective on ST

246:    Input Parameters:
247: +  st - the spectral transformation context
248: .  n  - number of vectors
249: -  V  - vectors to be checked

251:    Note:
252:    This function allows to handle singular pencils and to solve some problems
253:    in which the nullspace is important (see the users guide for details).
254:    
255:    Level: developer

257: .seealso: EPSAttachDeflationSpace()
258: @*/
259: PetscErrorCode STCheckNullSpace(ST st,int n,const Vec V[])
260: {

264:   if (n>0 && st->checknullspace) {
265:     (*st->checknullspace)(st,n,V);
266:   }
267:   return(0);
268: }