Actual source code: stsles.c

  2: /*
  3:     The ST (spectral transformation) interface routines related to the
  4:     KSP object associated to it.
  5: */

 7:  #include src/st/stimpl.h

 11: /*
 12:    STAssociatedKSPSolve - Solves the linear system of equations associated
 13:    to the spectral transformation.

 15:    Input Parameters:
 16: .  st - the spectral transformation context
 17: .  b  - right hand side vector

 19:    Output  Parameter:
 20: .  x - computed solution
 21: */
 22: PetscErrorCode STAssociatedKSPSolve(ST st,Vec b,Vec x)
 23: {
 25:   PetscInt       its;
 26:   KSPConvergedReason reason;

 32:   if (!st->ksp) { SETERRQ(PETSC_ERR_SUP,"ST has no associated KSP"); }
 33:   KSPSolve(st->ksp,b,x);
 34:   KSPGetConvergedReason(st->ksp,&reason);
 35:   if (reason<0) { SETERRQ1(0,"Warning: KSP did not converge (%d)",reason); }
 36:   KSPGetIterationNumber(st->ksp,&its);
 37:   st->lineariterations += its;
 38:   PetscInfo1(st,"Linear solve iterations=%d\n",its);
 39:   return(0);
 40: }

 44: /*
 45:    STAssociatedKSPSolveTranspose - Solves the transpose of the linear 
 46:    system of equations associated to the spectral transformation.

 48:    Input Parameters:
 49: .  st - the spectral transformation context
 50: .  b  - right hand side vector

 52:    Output  Parameter:
 53: .  x - computed solution
 54: */
 55: PetscErrorCode STAssociatedKSPSolveTranspose(ST st,Vec b,Vec x)
 56: {
 58:   PetscInt       its;
 59:   KSPConvergedReason reason;

 65:   if (!st->ksp) { SETERRQ(PETSC_ERR_SUP,"ST has no associated KSP"); }
 66:   KSPSolveTranspose(st->ksp,b,x);
 67:   KSPGetConvergedReason(st->ksp,&reason);
 68:   if (reason<0) { SETERRQ1(0,"Warning: KSP did not converge (%d)",reason); }
 69:   KSPGetIterationNumber(st->ksp,&its);
 70:   st->lineariterations += its;
 71:   PetscInfo1(st,"Linear solve iterations=%d\n",its);
 72:   return(0);
 73: }

 77: /*@
 78:    STSetKSP - Sets the KSP object associated with the spectral 
 79:    transformation.

 81:    Not collective

 83:    Input Parameters:
 84: +  st   - the spectral transformation context
 85: -  ksp  - the linear system context

 87:    Level: advanced

 89: @*/
 90: PetscErrorCode STSetKSP(ST st,KSP ksp)
 91: {

 98:   if (st->ksp) {
 99:     KSPDestroy(st->ksp);
100:   }
101:   st->ksp = ksp;
102:   PetscObjectReference((PetscObject)st->ksp);
103:   return(0);
104: }

108: /*@
109:    STGetKSP - Gets the KSP object associated with the spectral
110:    transformation.

112:    Not collective

114:    Input Parameter:
115: .  st - the spectral transformation context

117:    Output Parameter:
118: .  ksp  - the linear system context

120:    Notes:
121:    On output, the value of ksp can be PETSC_NULL if the combination of 
122:    eigenproblem type and selected transformation does not require to 
123:    solve a linear system of equations.
124:    
125:    Level: intermediate

127: @*/
128: PetscErrorCode STGetKSP(ST st,KSP* ksp)
129: {
132:   if (!st->type_name) { SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call STSetType first"); }
133:   if (ksp)  *ksp = st->ksp;
134:   return(0);
135: }

139: /*@
140:    STGetNumberLinearIterations - Gets the total number of linear iterations
141:    used by the ST object.

143:    Not Collective

145:    Input Parameter:
146: .  st - the spectral transformation context

148:    Output Parameter:
149: .  lits - number of linear iterations

151:    Level: intermediate

153: .seealso: STResetNumberLinearIterations()
154: @*/
155: PetscErrorCode STGetNumberLinearIterations(ST st,int* lits)
156: {
160:   *lits = st->lineariterations;
161:   return(0);
162: }

166: /*@
167:    STResetNumberLinearIterations - Resets the counter for total number of 
168:    linear iterations used by the ST object.

170:    Collective on ST

172:    Input Parameter:
173: .  st - the spectral transformation context

175:    Level: intermediate

177: .seealso: STGetNumberLinearIterations()
178: @*/
179: PetscErrorCode STResetNumberLinearIterations(ST st)
180: {
183:   st->lineariterations = 0;
184:   return(0);
185: }

189: PetscErrorCode STCheckNullSpace_Default(ST st,int n,const Vec V[])
190: {
192:   int            i,c;
193:   PetscReal      norm;
194:   Vec            *T,w;
195:   Mat            A;
196:   PC             pc;
197:   MatNullSpace   nullsp;
198: 
200:   PetscMalloc(n*sizeof(Vec),&T);
201:   KSPGetPC(st->ksp,&pc);
202:   PCGetOperators(pc,&A,PETSC_NULL,PETSC_NULL);
203:   MatGetVecs(A,PETSC_NULL,&w);
204:   c = 0;
205:   for (i=0;i<n;i++) {
206:     MatMult(A,V[i],w);
207:     VecNorm(w,NORM_2,&norm);
208:     if (norm < 1e-8) {
209:       PetscInfo2(st,"Vector %i norm=%g\n",i,norm);
210:       T[c] = V[i];
211:       c++;
212:     }
213:   }
214:   VecDestroy(w);
215:   if (c>0) {
216:     MatNullSpaceCreate(st->comm,PETSC_FALSE,c,T,&nullsp);
217:     KSPSetNullSpace(st->ksp,nullsp);
218:     MatNullSpaceDestroy(nullsp);
219:   }
220:   PetscFree(T);
221:   return(0);
222: }

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

231:    Collective on ST

233:    Input Parameters:
234: +  st - the spectral transformation context
235: .  n  - number of vectors
236: -  V  - vectors to be checked

238:    Note:
239:    This function allows to handle singular pencils and to solve some problems
240:    in which the nullspace is important (see the users guide for details).
241:    
242:    Level: developer

244: .seealso: EPSAttachDeflationSpace()
245: @*/
246: PetscErrorCode STCheckNullSpace(ST st,int n,const Vec V[])
247: {

251:   if (n>0 && st->checknullspace) {
252:     (*st->checknullspace)(st,n,V);
253:   }
254:   return(0);
255: }