Actual source code: iterativ.c

petsc-master 2019-06-19
Report Typos and Errors
  1: /*
  2:    This file contains some simple default routines.
  3:    These routines should be SHORT, since they will be included in every
  4:    executable image that uses the iterative routines (note that, through
  5:    the registry system, we provide a way to load only the truely necessary
  6:    files)
  7:  */
  8:  #include <petsc/private/kspimpl.h>
  9:  #include <petscdmshell.h>

 11: /*@
 12:    KSPGetResidualNorm - Gets the last (approximate preconditioned)
 13:    residual norm that has been computed.

 15:    Not Collective

 17:    Input Parameters:
 18: .  ksp - the iterative context

 20:    Output Parameters:
 21: .  rnorm - residual norm

 23:    Level: intermediate

 25: .seealso: KSPBuildResidual()
 26: @*/
 27: PetscErrorCode  KSPGetResidualNorm(KSP ksp,PetscReal *rnorm)
 28: {
 32:   *rnorm = ksp->rnorm;
 33:   return(0);
 34: }

 36: /*@
 37:    KSPGetIterationNumber - Gets the current iteration number; if the
 38:          KSPSolve() is complete, returns the number of iterations
 39:          used.

 41:    Not Collective

 43:    Input Parameters:
 44: .  ksp - the iterative context

 46:    Output Parameters:
 47: .  its - number of iterations

 49:    Level: intermediate

 51:    Notes:
 52:       During the ith iteration this returns i-1
 53: .seealso: KSPBuildResidual(), KSPGetResidualNorm(), KSPGetTotalIterations()
 54: @*/
 55: PetscErrorCode  KSPGetIterationNumber(KSP ksp,PetscInt *its)
 56: {
 60:   *its = ksp->its;
 61:   return(0);
 62: }

 64: /*@
 65:    KSPGetTotalIterations - Gets the total number of iterations this KSP object has performed since was created, counted over all linear solves

 67:    Not Collective

 69:    Input Parameters:
 70: .  ksp - the iterative context

 72:    Output Parameters:
 73: .  its - total number of iterations

 75:    Level: intermediate

 77:    Notes:
 78:     Use KSPGetIterationNumber() to get the count for the most recent solve only
 79:    If this is called within a linear solve (such as in a KSPMonitor routine) then it does not include iterations within that current solve

 81: .seealso: KSPBuildResidual(), KSPGetResidualNorm(), KSPGetIterationNumber()
 82: @*/
 83: PetscErrorCode  KSPGetTotalIterations(KSP ksp,PetscInt *its)
 84: {
 88:   *its = ksp->totalits;
 89:   return(0);
 90: }

 92: /*@C
 93:     KSPMonitorSingularValue - Prints the two norm of the true residual and
 94:     estimation of the extreme singular values of the preconditioned problem
 95:     at each iteration.

 97:     Logically Collective on ksp

 99:     Input Parameters:
100: +   ksp - the iterative context
101: .   n  - the iteration
102: -   rnorm - the two norm of the residual

104:     Options Database Key:
105: .   -ksp_monitor_singular_value - Activates KSPMonitorSingularValue()

107:     Notes:
108:     The CG solver uses the Lanczos technique for eigenvalue computation,
109:     while GMRES uses the Arnoldi technique; other iterative methods do
110:     not currently compute singular values.

112:     Level: intermediate

114: .seealso: KSPComputeExtremeSingularValues()
115: @*/
116: PetscErrorCode  KSPMonitorSingularValue(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat *dummy)
117: {
118:   PetscReal      emin,emax,c;
120:   PetscViewer    viewer = dummy->viewer;

125:   PetscViewerPushFormat(viewer,dummy->format);
126:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
127:   if (!ksp->calc_sings) {
128:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,(double)rnorm);
129:   } else {
130:     KSPComputeExtremeSingularValues(ksp,&emax,&emin);
131:     c    = emax/emin;
132:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e %% max %14.12e min %14.12e max/min %14.12e\n",n,(double)rnorm,(double)emax,(double)emin,(double)c);
133:   }
134:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
135:   PetscViewerPopFormat(viewer);
136:   return(0);
137: }

139: /*@C
140:    KSPMonitorSolution - Monitors progress of the KSP solvers by calling
141:    VecView() for the approximate solution at each iteration.

143:    Collective on ksp

145:    Input Parameters:
146: +  ksp - the KSP context
147: .  its - iteration number
148: .  fgnorm - 2-norm of residual (or gradient)
149: -  dummy - a viewer

151:    Level: intermediate

153:    Notes:
154:     For some Krylov methods such as GMRES constructing the solution at
155:   each iteration is expensive, hence using this will slow the code.

157: .seealso: KSPMonitorSet(), KSPMonitorDefault(), VecView()
158: @*/
159: PetscErrorCode  KSPMonitorSolution(KSP ksp,PetscInt its,PetscReal fgnorm,PetscViewerAndFormat *dummy)
160: {
162:   Vec            x;
163:   PetscViewer    viewer = dummy->viewer;

167:   KSPBuildSolution(ksp,NULL,&x);
168:   PetscViewerPushFormat(viewer,dummy->format);
169:   VecView(x,viewer);
170:   PetscViewerPopFormat(viewer);
171:   return(0);
172: }

174: /*@C
175:    KSPMonitorDefault - Print the residual norm at each iteration of an
176:    iterative solver.

178:    Collective on ksp

180:    Input Parameters:
181: +  ksp   - iterative context
182: .  n     - iteration number
183: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
184: -  dummy - an ASCII PetscViewer

186:    Level: intermediate

188: .seealso: KSPMonitorSet(), KSPMonitorTrueResidualNorm(), KSPMonitorLGResidualNormCreate()
189: @*/
190: PetscErrorCode  KSPMonitorDefault(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat *dummy)
191: {
193:   PetscViewer    viewer =  dummy->viewer;

197:   PetscViewerPushFormat(viewer,dummy->format);
198:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
199:   if (n == 0 && ((PetscObject)ksp)->prefix) {
200:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
201:   }
202:   PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %14.12e \n",n,(double)rnorm);
203:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
204:   PetscViewerPopFormat(viewer);
205:   return(0);
206: }

208: /*@C
209:    KSPMonitorTrueResidualNorm - Prints the true residual norm as well as the preconditioned
210:    residual norm at each iteration of an iterative solver.

212:    Collective on ksp

214:    Input Parameters:
215: +  ksp   - iterative context
216: .  n     - iteration number
217: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
218: -  dummy - an ASCII PetscViewer

220:    Options Database Key:
221: .  -ksp_monitor_true_residual - Activates KSPMonitorTrueResidualNorm()

223:    Notes:
224:    When using right preconditioning, these values are equivalent.

226:    Level: intermediate

228: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(),KSPMonitorTrueResidualMaxNorm()
229: @*/
230: PetscErrorCode  KSPMonitorTrueResidualNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat *dummy)
231: {
233:   Vec            resid;
234:   PetscReal      truenorm,bnorm;
235:   PetscViewer    viewer = dummy->viewer;
236:   char           normtype[256];

240:   PetscViewerPushFormat(viewer,dummy->format);
241:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
242:   if (n == 0 && ((PetscObject)ksp)->prefix) {
243:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
244:   }
245:   KSPBuildResidual(ksp,NULL,NULL,&resid);
246:   VecNorm(resid,NORM_2,&truenorm);
247:   VecDestroy(&resid);
248:   VecNorm(ksp->vec_rhs,NORM_2,&bnorm);
249:   PetscStrncpy(normtype,KSPNormTypes[ksp->normtype],sizeof(normtype));
250:   PetscStrtolower(normtype);
251:   PetscViewerASCIIPrintf(viewer,"%3D KSP %s resid norm %14.12e true resid norm %14.12e ||r(i)||/||b|| %14.12e\n",n,normtype,(double)rnorm,(double)truenorm,(double)(truenorm/bnorm));
252:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
253:   PetscViewerPopFormat(viewer);
254:   return(0);
255: }

257: /*@C
258:    KSPMonitorTrueResidualMaxNorm - Prints the true residual max norm each iteration of an iterative solver.

260:    Collective on ksp

262:    Input Parameters:
263: +  ksp   - iterative context
264: .  n     - iteration number
265: .  rnorm - norm (preconditioned) residual value (may be estimated).
266: -  dummy - an ASCII viewer

268:    Options Database Key:
269: .  -ksp_monitor_max - Activates KSPMonitorTrueResidualMaxNorm()

271:    Notes:
272:    This could be implemented (better) with a flag in ksp.

274:    Level: intermediate

276: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(),KSPMonitorTrueResidualNorm()
277: @*/
278: PetscErrorCode  KSPMonitorTrueResidualMaxNorm(KSP ksp,PetscInt n,PetscReal rnorm,PetscViewerAndFormat *dummy)
279: {
281:   Vec            resid;
282:   PetscReal      truenorm,bnorm;
283:   PetscViewer    viewer = dummy->viewer;
284:   char           normtype[256];

288:   PetscViewerPushFormat(viewer,dummy->format);
289:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
290:   if (n == 0 && ((PetscObject)ksp)->prefix) {
291:     PetscViewerASCIIPrintf(viewer,"  Residual norms (max) for %s solve.\n",((PetscObject)ksp)->prefix);
292:   }
293:   KSPBuildResidual(ksp,NULL,NULL,&resid);
294:   VecNorm(resid,NORM_INFINITY,&truenorm);
295:   VecDestroy(&resid);
296:   VecNorm(ksp->vec_rhs,NORM_INFINITY,&bnorm);
297:   PetscStrncpy(normtype,KSPNormTypes[ksp->normtype],sizeof(normtype));
298:   PetscStrtolower(normtype);
299:   PetscViewerASCIIPrintf(viewer,"%3D KSP true resid max norm %14.12e ||r(i)||/||b|| %14.12e\n",n,(double)truenorm,(double)(truenorm/bnorm));
300:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
301:   PetscViewerPopFormat(viewer);
302:   return(0);
303: }

305: PetscErrorCode  KSPMonitorRange_Private(KSP ksp,PetscInt it,PetscReal *per)
306: {
308:   Vec            resid;
309:   PetscReal      rmax,pwork;
310:   PetscInt       i,n,N;
311:   const PetscScalar *r;

314:   KSPBuildResidual(ksp,NULL,NULL,&resid);
315:   VecNorm(resid,NORM_INFINITY,&rmax);
316:   VecGetLocalSize(resid,&n);
317:   VecGetSize(resid,&N);
318:   VecGetArrayRead(resid,&r);
319:   pwork = 0.0;
320:   for (i=0; i<n; i++) pwork += (PetscAbsScalar(r[i]) > .20*rmax);
321:   VecRestoreArrayRead(resid,&r);
322:   VecDestroy(&resid);
323:   MPIU_Allreduce(&pwork,per,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)ksp));
324:   *per = *per/N;
325:   return(0);
326: }

328: /*@C
329:    KSPMonitorRange - Prints the percentage of residual elements that are more then 10 percent of the maximum value.

331:    Collective on ksp

333:    Input Parameters:
334: +  ksp   - iterative context
335: .  it    - iteration number
336: .  rnorm - 2-norm (preconditioned) residual value (may be estimated).
337: -  dummy - an ASCII viewer

339:    Options Database Key:
340: .  -ksp_monitor_range - Activates KSPMonitorRange()

342:    Level: intermediate

344: .seealso: KSPMonitorSet(), KSPMonitorDefault(), KSPMonitorLGResidualNormCreate()
345: @*/
346: PetscErrorCode  KSPMonitorRange(KSP ksp,PetscInt it,PetscReal rnorm,PetscViewerAndFormat *dummy)
347: {
348:   PetscErrorCode   ierr;
349:   PetscReal        perc,rel;
350:   PetscViewer      viewer = dummy->viewer;
351:   /* should be in a MonitorRangeContext */
352:   static PetscReal prev;

356:   PetscViewerPushFormat(viewer,dummy->format);
357:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
358:   if (!it) prev = rnorm;
359:   if (it == 0 && ((PetscObject)ksp)->prefix) {
360:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
361:   }
362:   KSPMonitorRange_Private(ksp,it,&perc);

364:   rel  = (prev - rnorm)/prev;
365:   prev = rnorm;
366:   PetscViewerASCIIPrintf(viewer,"%3D KSP preconditioned resid norm %14.12e Percent values above 20 percent of maximum %5.2f relative decrease %5.2e ratio %5.2e \n",it,(double)rnorm,(double)(100.0*perc),(double)rel,(double)(rel/perc));
367:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
368:   PetscViewerPopFormat(viewer);
369:   return(0);
370: }

372: /*@C
373:    KSPMonitorDynamicTolerance - Recompute the inner tolerance in every
374:    outer iteration in an adaptive way.

376:    Collective on ksp

378:    Input Parameters:
379: +  ksp   - iterative context
380: .  n     - iteration number (not used)
381: .  fnorm - the current residual norm
382: .  dummy - some context as a C struct. fields:
383:              coef: a scaling coefficient. default 1.0. can be passed through
384:                    -sub_ksp_dynamic_tolerance_param
385:              bnrm: norm of the right-hand side. store it to avoid repeated calculation

387:    Notes:
388:    This may be useful for a flexibly preconditioner Krylov method to
389:    control the accuracy of the inner solves needed to gaurantee the
390:    convergence of the outer iterations.

392:    Level: advanced

394: .seealso: KSPMonitorDynamicToleranceDestroy()
395: @*/
396: PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy)
397: {
399:   PC             pc;
400:   PetscReal      outer_rtol, outer_abstol, outer_dtol, inner_rtol;
401:   PetscInt       outer_maxits,nksp,first,i;
402:   KSPDynTolCtx   *scale   = (KSPDynTolCtx*)dummy;
403:   KSP            *subksp = NULL;
404:   PetscBool      isksp;

407:   KSPGetPC(ksp, &pc);

409:   /* compute inner_rtol */
410:   if (scale->bnrm < 0.0) {
411:     Vec b;
412:     KSPGetRhs(ksp, &b);
413:     VecNorm(b, NORM_2, &(scale->bnrm));
414:   }
415:   KSPGetTolerances(ksp, &outer_rtol, &outer_abstol, &outer_dtol, &outer_maxits);
416:   inner_rtol = PetscMin(scale->coef * scale->bnrm * outer_rtol / fnorm, 0.999);
417:   /*PetscPrintf(PETSC_COMM_WORLD, "        Inner rtol = %g\n", (double)inner_rtol);*/

419:   /* if pc is ksp */
420:   PetscObjectTypeCompare((PetscObject)pc,PCKSP,&isksp);
421:   if (isksp) {
422:     KSP kspinner;

424:     PCKSPGetKSP(pc, &kspinner);
425:     KSPSetTolerances(kspinner, inner_rtol, outer_abstol, outer_dtol, outer_maxits);
426:     return(0);
427:   }

429:   /* if pc is bjacobi */
430:   PCBJacobiGetSubKSP(pc, &nksp, &first, &subksp);
431:   if (subksp) {
432:     for (i=0; i<nksp; i++) {
433:       KSPSetTolerances(subksp[i], inner_rtol, outer_abstol, outer_dtol, outer_maxits);
434:     }
435:     return(0);
436:   }

438:   /* todo: dynamic tolerance may apply to other types of pc too */
439:   return(0);
440: }

442: /*
443:   Destroy the dummy context used in KSPMonitorDynamicTolerance()
444: */
445: PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy)
446: {

450:   PetscFree(*dummy);
451:   return(0);
452: }

454: /*
455:   Default (short) KSP Monitor, same as KSPMonitorDefault() except
456:   it prints fewer digits of the residual as the residual gets smaller.
457:   This is because the later digits are meaningless and are often
458:   different on different machines; by using this routine different
459:   machines will usually generate the same output.

461:   Deprecated: Intentionally has no manual page
462: */
463: PetscErrorCode  KSPMonitorDefaultShort(KSP ksp,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *dummy)
464: {
466:   PetscViewer    viewer = dummy->viewer;

470:   PetscViewerPushFormat(viewer,dummy->format);
471:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
472:   if (its == 0 && ((PetscObject)ksp)->prefix) {
473:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
474:   }

476:   if (fnorm > 1.e-9) {
477:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %g \n",its,(double)fnorm);
478:   } else if (fnorm > 1.e-11) {
479:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %5.3e \n",its,(double)fnorm);
480:   } else {
481:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm < 1.e-11\n",its);
482:   }
483:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
484:   PetscViewerPopFormat(viewer);
485:   return(0);
486: }

488: /*@C
489:    KSPConvergedSkip - Convergence test that do not return as converged
490:    until the maximum number of iterations is reached.

492:    Collective on ksp

494:    Input Parameters:
495: +  ksp   - iterative context
496: .  n     - iteration number
497: .  rnorm - 2-norm residual value (may be estimated)
498: -  dummy - unused convergence context

500:    Returns:
501: .  reason - KSP_CONVERGED_ITERATING, KSP_CONVERGED_ITS

503:    Notes:
504:    This should be used as the convergence test with the option
505:    KSPSetNormType(ksp,KSP_NORM_NONE), since norms of the residual are
506:    not computed. Convergence is then declared after the maximum number
507:    of iterations have been reached. Useful when one is using CG or
508:    BiCGStab as a smoother.

510:    Level: advanced

512: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSetNormType()
513: @*/
514: PetscErrorCode  KSPConvergedSkip(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
515: {
519:   *reason = KSP_CONVERGED_ITERATING;
520:   if (n >= ksp->max_it) *reason = KSP_CONVERGED_ITS;
521:   return(0);
522: }


525: /*@C
526:    KSPConvergedDefaultCreate - Creates and initializes the space used by the KSPConvergedDefault() function context

528:    Note Collective

530:    Output Parameter:
531: .  ctx - convergence context

533:    Level: intermediate

535: .seealso: KSPConvergedDefault(), KSPConvergedDefaultDestroy(), KSPSetConvergenceTest(), KSPSetTolerances(),
536:           KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm()
537: @*/
538: PetscErrorCode  KSPConvergedDefaultCreate(void **ctx)
539: {
540:   PetscErrorCode         ierr;
541:   KSPConvergedDefaultCtx *cctx;

544:   PetscNew(&cctx);
545:   *ctx = cctx;
546:   return(0);
547: }

549: /*@
550:    KSPConvergedDefaultSetUIRNorm - makes the default convergence test use || B*(b - A*(initial guess))||
551:       instead of || B*b ||. In the case of right preconditioner or if KSPSetNormType(ksp,KSP_NORM_UNPRECONDIITONED)
552:       is used there is no B in the above formula. UIRNorm is short for Use Initial Residual Norm.

554:    Collective on ksp

556:    Input Parameters:
557: .  ksp   - iterative context

559:    Options Database:
560: .   -ksp_converged_use_initial_residual_norm

562:    Notes:
563:    Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.

565:    The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
566:    are defined in petscksp.h.

568:    If the convergence test is not KSPConvergedDefault() then this is ignored.

570:    If right preconditioning is being used then B does not appear in the above formula.


573:    Level: intermediate

575: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUMIRNorm()
576: @*/
577: PetscErrorCode  KSPConvergedDefaultSetUIRNorm(KSP ksp)
578: {
579:   KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx*) ksp->cnvP;

583:   if (ksp->converged != KSPConvergedDefault) return(0);
584:   if (ctx->mininitialrtol) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
585:   ctx->initialrtol = PETSC_TRUE;
586:   return(0);
587: }

589: /*@
590:    KSPConvergedDefaultSetUMIRNorm - makes the default convergence test use min(|| B*(b - A*(initial guess))||,|| B*b ||)
591:       In the case of right preconditioner or if KSPSetNormType(ksp,KSP_NORM_UNPRECONDIITONED)
592:       is used there is no B in the above formula. UMIRNorm is short for Use Minimum Initial Residual Norm.

594:    Collective on ksp

596:    Input Parameters:
597: .  ksp   - iterative context

599:    Options Database:
600: .   -ksp_converged_use_min_initial_residual_norm

602:    Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.

604:    The precise values of reason are macros such as KSP_CONVERGED_RTOL, which
605:    are defined in petscksp.h.

607:    Level: intermediate

609: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUIRNorm()
610: @*/
611: PetscErrorCode  KSPConvergedDefaultSetUMIRNorm(KSP ksp)
612: {
613:   KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx*) ksp->cnvP;

617:   if (ksp->converged != KSPConvergedDefault) return(0);
618:   if (ctx->initialrtol) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
619:   ctx->mininitialrtol = PETSC_TRUE;
620:   return(0);
621: }

623: /*@C
624:    KSPConvergedDefault - Determines convergence of the linear iterative solvers by default

626:    Collective on ksp

628:    Input Parameters:
629: +  ksp   - iterative context
630: .  n     - iteration number
631: .  rnorm - residual norm (may be estimated, depending on the method may be the preconditioned residual norm)
632: -  ctx - convergence context which must be created by KSPConvergedDefaultCreate()

634:    Output Parameter:
635: +   positive - if the iteration has converged;
636: .   negative - if residual norm exceeds divergence threshold;
637: -   0 - otherwise.

639:    Notes:
640:    KSPConvergedDefault() reaches convergence when   rnorm < MAX (rtol * rnorm_0, abstol);
641:    Divergence is detected if  rnorm > dtol * rnorm_0,

643:    where:
644: +     rtol = relative tolerance,
645: .     abstol = absolute tolerance.
646: .     dtol = divergence tolerance,
647: -     rnorm_0 is the two norm of the right hand side (or the preconditioned norm, depending on what was set with
648:           KSPSetNormType(). When initial guess is non-zero you
649:           can call KSPConvergedDefaultSetUIRNorm() to use the norm of (b - A*(initial guess))
650:           as the starting point for relative norm convergence testing, that is as rnorm_0

652:    Use KSPSetTolerances() to alter the defaults for rtol, abstol, dtol.

654:    Use KSPSetNormType() (or -ksp_norm_type <none,preconditioned,unpreconditioned,natural>) to change the norm used for computing rnorm

656:    The precise values of reason are macros such as KSP_CONVERGED_RTOL, which are defined in petscksp.h.

658:    This routine is used by KSP by default so the user generally never needs call it directly.

660:    Use KSPSetConvergenceTest() to provide your own test instead of using this one.

662:    Level: intermediate

664: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(),
665:           KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm(), KSPConvergedDefaultCreate(), KSPConvergedDefaultDestroy()
666: @*/
667: PetscErrorCode  KSPConvergedDefault(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
668: {
669:   PetscErrorCode         ierr;
670:   KSPConvergedDefaultCtx *cctx = (KSPConvergedDefaultCtx*) ctx;
671:   KSPNormType            normtype;

676:   *reason = KSP_CONVERGED_ITERATING;

678:   KSPGetNormType(ksp,&normtype);
679:   if (normtype == KSP_NORM_NONE) return(0);

681:   if (!cctx) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_NULL,"Convergence context must have been created with KSPConvergedDefaultCreate()");
682:   if (!n) {
683:     /* if user gives initial guess need to compute norm of b */
684:     if (!ksp->guess_zero && !cctx->initialrtol) {
685:       PetscReal snorm = 0.0;
686:       if (ksp->normtype == KSP_NORM_UNPRECONDITIONED || ksp->pc_side == PC_RIGHT) {
687:         PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of RHS\n");
688:         VecNorm(ksp->vec_rhs,NORM_2,&snorm);        /*     <- b'*b */
689:       } else {
690:         Vec z;
691:         /* Should avoid allocating the z vector each time but cannot stash it in cctx because if KSPReset() is called the vector size might change */
692:         VecDuplicate(ksp->vec_rhs,&z);
693:         KSP_PCApply(ksp,ksp->vec_rhs,z);
694:         if (ksp->normtype == KSP_NORM_PRECONDITIONED) {
695:           PetscInfo(ksp,"user has provided nonzero initial guess, computing 2-norm of preconditioned RHS\n");
696:           VecNorm(z,NORM_2,&snorm);                 /*    dp <- b'*B'*B*b */
697:         } else if (ksp->normtype == KSP_NORM_NATURAL) {
698:           PetscScalar norm;
699:           PetscInfo(ksp,"user has provided nonzero initial guess, computing natural norm of RHS\n");
700:           VecDot(ksp->vec_rhs,z,&norm);
701:           snorm = PetscSqrtReal(PetscAbsScalar(norm));                            /*    dp <- b'*B*b */
702:         }
703:         VecDestroy(&z);
704:       }
705:       /* handle special case of zero RHS and nonzero guess */
706:       if (!snorm) {
707:         PetscInfo(ksp,"Special case, user has provided nonzero initial guess and zero RHS\n");
708:         snorm = rnorm;
709:       }
710:       if (cctx->mininitialrtol) ksp->rnorm0 = PetscMin(snorm,rnorm);
711:       else ksp->rnorm0 = snorm;
712:     } else {
713:       ksp->rnorm0 = rnorm;
714:     }
715:     ksp->ttol = PetscMax(ksp->rtol*ksp->rnorm0,ksp->abstol);
716:   }

718:   if (n <= ksp->chknorm) return(0);

720:   if (PetscIsInfOrNanReal(rnorm)) {
721:     PCFailedReason pcreason;
722:     PetscInt       sendbuf,pcreason_max;
723:     PCGetFailedReason(ksp->pc,&pcreason);
724:     sendbuf = (PetscInt)pcreason;
725:     MPI_Allreduce(&sendbuf,&pcreason_max,1,MPIU_INT,MPIU_MAX,PetscObjectComm((PetscObject)ksp));
726:     if (pcreason_max) {
727:       *reason = KSP_DIVERGED_PC_FAILED;
728:       VecSetInf(ksp->vec_sol);
729:       PetscInfo(ksp,"Linear solver pcsetup fails, declaring divergence \n");
730:     } else {
731:       *reason = KSP_DIVERGED_NANORINF;
732:       PetscInfo(ksp,"Linear solver has created a not a number (NaN) as the residual norm, declaring divergence \n");
733:     }
734:   } else if (rnorm <= ksp->ttol) {
735:     if (rnorm < ksp->abstol) {
736:       PetscInfo3(ksp,"Linear solver has converged. Residual norm %14.12e is less than absolute tolerance %14.12e at iteration %D\n",(double)rnorm,(double)ksp->abstol,n);
737:       *reason = KSP_CONVERGED_ATOL;
738:     } else {
739:       if (cctx->initialrtol) {
740:         PetscInfo4(ksp,"Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial residual norm %14.12e at iteration %D\n",(double)rnorm,(double)ksp->rtol,(double)ksp->rnorm0,n);
741:       } else {
742:         PetscInfo4(ksp,"Linear solver has converged. Residual norm %14.12e is less than relative tolerance %14.12e times initial right hand side norm %14.12e at iteration %D\n",(double)rnorm,(double)ksp->rtol,(double)ksp->rnorm0,n);
743:       }
744:       *reason = KSP_CONVERGED_RTOL;
745:     }
746:   } else if (rnorm >= ksp->divtol*ksp->rnorm0) {
747:     PetscInfo3(ksp,"Linear solver is diverging. Initial right hand size norm %14.12e, current residual norm %14.12e at iteration %D\n",(double)ksp->rnorm0,(double)rnorm,n);
748:     *reason = KSP_DIVERGED_DTOL;
749:   }
750:   return(0);
751: }

753: /*@C
754:    KSPConvergedDefaultDestroy - Frees the space used by the KSPConvergedDefault() function context

756:    Not Collective

758:    Input Parameters:
759: .  ctx - convergence context

761:    Level: intermediate

763: .seealso: KSPConvergedDefault(), KSPConvergedDefaultCreate(), KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(),
764:           KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm()
765: @*/
766: PetscErrorCode  KSPConvergedDefaultDestroy(void *ctx)
767: {
768:   PetscErrorCode         ierr;
769:   KSPConvergedDefaultCtx *cctx = (KSPConvergedDefaultCtx*) ctx;

772:   VecDestroy(&cctx->work);
773:   PetscFree(ctx);
774:   return(0);
775: }

777: /*
778:    KSPBuildSolutionDefault - Default code to create/move the solution.

780:    Collective on ksp

782:    Input Parameters:
783: +  ksp - iterative context
784: -  v   - pointer to the user's vector

786:    Output Parameter:
787: .  V - pointer to a vector containing the solution

789:    Level: advanced

791:    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations

793: .seealso: KSPGetSolution(), KSPBuildResidualDefault()
794: */
795: PetscErrorCode KSPBuildSolutionDefault(KSP ksp,Vec v,Vec *V)
796: {

800:   if (ksp->pc_side == PC_RIGHT) {
801:     if (ksp->pc) {
802:       if (v) {
803:         KSP_PCApply(ksp,ksp->vec_sol,v); *V = v;
804:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with right preconditioner");
805:     } else {
806:       if (v) {
807:         VecCopy(ksp->vec_sol,v); *V = v;
808:       } else *V = ksp->vec_sol;
809:     }
810:   } else if (ksp->pc_side == PC_SYMMETRIC) {
811:     if (ksp->pc) {
812:       if (ksp->transpose_solve) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with symmetric preconditioner and transpose solve");
813:       if (v) {
814:         PCApplySymmetricRight(ksp->pc,ksp->vec_sol,v);
815:         *V = v;
816:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with symmetric preconditioner");
817:     } else {
818:       if (v) {
819:         VecCopy(ksp->vec_sol,v); *V = v;
820:       } else *V = ksp->vec_sol;
821:     }
822:   } else {
823:     if (v) {
824:       VecCopy(ksp->vec_sol,v); *V = v;
825:     } else *V = ksp->vec_sol;
826:   }
827:   return(0);
828: }

830: /*
831:    KSPBuildResidualDefault - Default code to compute the residual.

833:    Collecive on ksp

835:    Input Parameters:
836: .  ksp - iterative context
837: .  t   - pointer to temporary vector
838: .  v   - pointer to user vector

840:    Output Parameter:
841: .  V - pointer to a vector containing the residual

843:    Level: advanced

845:    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations

847: .seealso: KSPBuildSolutionDefault()
848: */
849: PetscErrorCode KSPBuildResidualDefault(KSP ksp,Vec t,Vec v,Vec *V)
850: {
852:   Mat            Amat,Pmat;

855:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
856:   PCGetOperators(ksp->pc,&Amat,&Pmat);
857:   KSPBuildSolution(ksp,t,NULL);
858:   KSP_MatMult(ksp,Amat,t,v);
859:   VecAYPX(v,-1.0,ksp->vec_rhs);
860:   *V   = v;
861:   return(0);
862: }

864: /*@C
865:   KSPCreateVecs - Gets a number of work vectors.

867:   Collective on ksp

869:   Input Parameters:
870: + ksp  - iterative context
871: . rightn  - number of right work vectors
872: - leftn   - number of left work vectors to allocate

874:   Output Parameter:
875: +  right - the array of vectors created
876: -  left - the array of left vectors

878:    Note: The right vector has as many elements as the matrix has columns. The left
879:      vector has as many elements as the matrix has rows.

881:    The vectors are new vectors that are not owned by the KSP, they should be destroyed with calls to VecDestroyVecs() when no longer needed.

883:    Developers Note: First tries to duplicate the rhs and solution vectors of the KSP, if they do not exist tries to get them from the matrix, if
884:                     that does not exist tries to get them from the DM (if it is provided).

886:    Level: advanced

888: .seealso:   MatCreateVecs(), VecDestroyVecs()

890: @*/
891: PetscErrorCode KSPCreateVecs(KSP ksp,PetscInt rightn, Vec **right,PetscInt leftn,Vec **left)
892: {
894:   Vec            vecr = NULL,vecl = NULL;
895:   PetscBool      matset,pmatset;
896:   Mat            mat = NULL;

899:   if (rightn) {
900:     if (!right) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"You asked for right vectors but did not pass a pointer to hold them");
901:     if (ksp->vec_sol) vecr = ksp->vec_sol;
902:     else {
903:       if (ksp->pc) {
904:         PCGetOperatorsSet(ksp->pc,&matset,&pmatset);
905:         /* check for mat before pmat because for KSPLSQR pmat may be a different size than mat since pmat maybe mat'*mat */
906:         if (matset) {
907:           PCGetOperators(ksp->pc,&mat,NULL);
908:           MatCreateVecs(mat,&vecr,NULL);
909:         } else if (pmatset) {
910:           PCGetOperators(ksp->pc,NULL,&mat);
911:           MatCreateVecs(mat,&vecr,NULL);
912:         }
913:       }
914:       if (!vecr) {
915:         if (ksp->dm) {
916:           DMGetGlobalVector(ksp->dm,&vecr);
917:         } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You requested a vector from a KSP that cannot provide one");
918:       }
919:     }
920:     VecDuplicateVecs(vecr,rightn,right);
921:     if (!ksp->vec_sol) {
922:       if (mat) {
923:         VecDestroy(&vecr);
924:       } else if (ksp->dm) {
925:         DMRestoreGlobalVector(ksp->dm,&vecr);
926:       }
927:     }
928:   }
929:   if (leftn) {
930:     if (!left) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_INCOMP,"You asked for left vectors but did not pass a pointer to hold them");
931:     if (ksp->vec_rhs) vecl = ksp->vec_rhs;
932:     else {
933:       if (ksp->pc) {
934:         PCGetOperatorsSet(ksp->pc,&matset,&pmatset);
935:         /* check for mat before pmat because for KSPLSQR pmat may be a different size than mat since pmat maybe mat'*mat */
936:         if (matset) {
937:           PCGetOperators(ksp->pc,&mat,NULL);
938:           MatCreateVecs(mat,NULL,&vecl);
939:         } else if (pmatset) {
940:           PCGetOperators(ksp->pc,NULL,&mat);
941:           MatCreateVecs(mat,NULL,&vecl);
942:         }
943:       }
944:       if (!vecl) {
945:         if (ksp->dm) {
946:           DMGetGlobalVector(ksp->dm,&vecl);
947:         } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You requested a vector from a KSP that cannot provide one");
948:       }
949:     }
950:     VecDuplicateVecs(vecl,leftn,left);
951:     if (!ksp->vec_rhs) {
952:       if (mat) {
953:         VecDestroy(&vecl);
954:       } else if (ksp->dm) {
955:         DMRestoreGlobalVector(ksp->dm,&vecl);
956:       }
957:     }
958:   }
959:   return(0);
960: }

962: /*@C
963:   KSPSetWorkVecs - Sets a number of work vectors into a KSP object

965:   Collective on ksp

967:   Input Parameters:
968: + ksp  - iterative context
969: - nw   - number of work vectors to allocate

971:   Level: developer

973:   Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations

975: @*/
976: PetscErrorCode KSPSetWorkVecs(KSP ksp,PetscInt nw)
977: {

981:   VecDestroyVecs(ksp->nwork,&ksp->work);
982:   ksp->nwork = nw;
983:   KSPCreateVecs(ksp,nw,&ksp->work,0,NULL);
984:   PetscLogObjectParents(ksp,nw,ksp->work);
985:   return(0);
986: }

988: /*
989:   KSPDestroyDefault - Destroys a iterative context variable for methods with
990:   no separate context.  Preferred calling sequence KSPDestroy().

992:   Input Parameter:
993: . ksp - the iterative context

995:    Developers Note: This is PETSC_EXTERN because it may be used by user written plugin KSP implementations

997: */
998: PetscErrorCode KSPDestroyDefault(KSP ksp)
999: {

1004:   PetscFree(ksp->data);
1005:   return(0);
1006: }

1008: /*@
1009:    KSPGetConvergedReason - Gets the reason the KSP iteration was stopped.

1011:    Not Collective

1013:    Input Parameter:
1014: .  ksp - the KSP context

1016:    Output Parameter:
1017: .  reason - negative value indicates diverged, positive value converged, see KSPConvergedReason

1019:    Possible values for reason: See also manual page for each reason
1020: $  KSP_CONVERGED_RTOL (residual 2-norm decreased by a factor of rtol, from 2-norm of right hand side)
1021: $  KSP_CONVERGED_ATOL (residual 2-norm less than abstol)
1022: $  KSP_CONVERGED_ITS (used by the preonly preconditioner that always uses ONE iteration, or when the KSPConvergedSkip() convergence test routine is set.
1023: $  KSP_CONVERGED_CG_NEG_CURVE (see note below)
1024: $  KSP_CONVERGED_CG_CONSTRAINED (see note below)
1025: $  KSP_CONVERGED_STEP_LENGTH (see note below)
1026: $  KSP_CONVERGED_ITERATING (returned if the solver is not yet finished)
1027: $  KSP_DIVERGED_ITS  (required more than its to reach convergence)
1028: $  KSP_DIVERGED_DTOL (residual norm increased by a factor of divtol)
1029: $  KSP_DIVERGED_NANORINF (residual norm became Not-a-number or Inf likely due to 0/0)
1030: $  KSP_DIVERGED_BREAKDOWN (generic breakdown in method)
1031: $  KSP_DIVERGED_BREAKDOWN_BICG (Initial residual is orthogonal to preconditioned initial residual. Try a different preconditioner, or a different initial Level.)

1033:    Options Database:
1034: .   -ksp_converged_reason - prints the reason to standard out

1036:    Notes:
1037:     If this routine is called before or doing the KSPSolve() the value of KSP_CONVERGED_ITERATING is returned

1039:    The values  KSP_CONVERGED_CG_NEG_CURVE, KSP_CONVERGED_CG_CONSTRAINED, and KSP_CONVERGED_STEP_LENGTH are returned only by the special KSPNASH, KSPSTCG, and KSPGLTR
1040:    solvers which are used by the SNESNEWTONTR (trust region) solver.

1042:    Level: intermediate

1044: .seealso: KSPSetConvergenceTest(), KSPConvergedDefault(), KSPSetTolerances(), KSPConvergedReason
1045: @*/
1046: PetscErrorCode  KSPGetConvergedReason(KSP ksp,KSPConvergedReason *reason)
1047: {
1051:   *reason = ksp->reason;
1052:   return(0);
1053: }

1055:  #include <petsc/private/dmimpl.h>
1056: /*@
1057:    KSPSetDM - Sets the DM that may be used by some preconditioners

1059:    Logically Collective on ksp

1061:    Input Parameters:
1062: +  ksp - the preconditioner context
1063: -  dm - the dm, cannot be NULL

1065:    Notes:
1066:    If this is used then the KSP will attempt to use the DM to create the matrix and use the routine set with
1067:    DMKSPSetComputeOperators(). Use KSPSetDMActive(ksp,PETSC_FALSE) to instead use the matrix you've provided with
1068:    KSPSetOperators().

1070:    A DM can only be used for solving one problem at a time because information about the problem is stored on the DM,
1071:    even when not using interfaces like DMKSPSetComputeOperators().  Use DMClone() to get a distinct DM when solving
1072:    different problems using the same function space.

1074:    Level: intermediate

1076: .seealso: KSPGetDM(), KSPSetDMActive(), KSPSetComputeOperators(), KSPSetComputeRHS(), KSPSetComputeInitialGuess(), DMKSPSetComputeOperators(), DMKSPSetComputeRHS(), DMKSPSetComputeInitialGuess()
1077: @*/
1078: PetscErrorCode  KSPSetDM(KSP ksp,DM dm)
1079: {
1081:   PC             pc;

1086:   PetscObjectReference((PetscObject)dm);
1087:   if (ksp->dm) {                /* Move the DMSNES context over to the new DM unless the new DM already has one */
1088:     if (ksp->dm->dmksp && !dm->dmksp) {
1089:       DMKSP kdm;
1090:       DMCopyDMKSP(ksp->dm,dm);
1091:       DMGetDMKSP(ksp->dm,&kdm);
1092:       if (kdm->originaldm == ksp->dm) kdm->originaldm = dm; /* Grant write privileges to the replacement DM */
1093:     }
1094:     DMDestroy(&ksp->dm);
1095:   }
1096:   ksp->dm       = dm;
1097:   ksp->dmAuto   = PETSC_FALSE;
1098:   KSPGetPC(ksp,&pc);
1099:   PCSetDM(pc,dm);
1100:   ksp->dmActive = PETSC_TRUE;
1101:   return(0);
1102: }

1104: /*@
1105:    KSPSetDMActive - Indicates the DM should be used to generate the linear system matrix and right hand side

1107:    Logically Collective on ksp

1109:    Input Parameters:
1110: +  ksp - the preconditioner context
1111: -  flg - use the DM

1113:    Notes:
1114:    By default KSPSetDM() sets the DM as active, call KSPSetDMActive(ksp,PETSC_FALSE); after KSPSetDM(ksp,dm) to not have the KSP object use the DM to generate the matrices.

1116:    Level: intermediate

1118: .seealso: KSPGetDM(), KSPSetDM(), SNESSetDM(), KSPSetComputeOperators(), KSPSetComputeRHS(), KSPSetComputeInitialGuess()
1119: @*/
1120: PetscErrorCode  KSPSetDMActive(KSP ksp,PetscBool flg)
1121: {
1125:   ksp->dmActive = flg;
1126:   return(0);
1127: }

1129: /*@
1130:    KSPGetDM - Gets the DM that may be used by some preconditioners

1132:    Not Collective

1134:    Input Parameter:
1135: . ksp - the preconditioner context

1137:    Output Parameter:
1138: .  dm - the dm

1140:    Level: intermediate


1143: .seealso: KSPSetDM(), KSPSetDMActive()
1144: @*/
1145: PetscErrorCode  KSPGetDM(KSP ksp,DM *dm)
1146: {

1151:   if (!ksp->dm) {
1152:     DMShellCreate(PetscObjectComm((PetscObject)ksp),&ksp->dm);
1153:     ksp->dmAuto = PETSC_TRUE;
1154:   }
1155:   *dm = ksp->dm;
1156:   return(0);
1157: }

1159: /*@
1160:    KSPSetApplicationContext - Sets the optional user-defined context for the linear solver.

1162:    Logically Collective on ksp

1164:    Input Parameters:
1165: +  ksp - the KSP context
1166: -  usrP - optional user context

1168:    Fortran Notes:
1169:     To use this from Fortran you must write a Fortran interface definition for this
1170:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

1172:    Level: intermediate

1174: .seealso: KSPGetApplicationContext()
1175: @*/
1176: PetscErrorCode  KSPSetApplicationContext(KSP ksp,void *usrP)
1177: {
1179:   PC             pc;

1183:   ksp->user = usrP;
1184:   KSPGetPC(ksp,&pc);
1185:   PCSetApplicationContext(pc,usrP);
1186:   return(0);
1187: }

1189: /*@
1190:    KSPGetApplicationContext - Gets the user-defined context for the linear solver.

1192:    Not Collective

1194:    Input Parameter:
1195: .  ksp - KSP context

1197:    Output Parameter:
1198: .  usrP - user context

1200:    Fortran Notes:
1201:     To use this from Fortran you must write a Fortran interface definition for this
1202:     function that tells Fortran the Fortran derived data type that you are passing in as the ctx argument.

1204:    Level: intermediate

1206: .seealso: KSPSetApplicationContext()
1207: @*/
1208: PetscErrorCode  KSPGetApplicationContext(KSP ksp,void *usrP)
1209: {
1212:   *(void**)usrP = ksp->user;
1213:   return(0);
1214: }

1216:  #include <petsc/private/pcimpl.h>

1218: /*@
1219:    KSPCheckSolve - Checks if the PCSetUp() or KSPSolve() failed and set the error flag for the outer PC. A KSP_DIVERGED_ITS is
1220:          not considered a failure in this context

1222:    Collective on ksp

1224:    Input Parameter:
1225: +  ksp - the linear solver (KSP) context.
1226: .  pc - the preconditioner context
1227: -  vec - a vector that will be initialized with Inf to indicate lack of convergence

1229:    Notes: this may be called by a subset of the processes in the PC

1231:    Level: developer

1233:    Developer Note: this is used to manage returning from preconditioners whose inner KSP solvers have failed in some way

1235: .seealso: KSPCreate(), KSPSetType(), KSP, KSPCheckNorm(), KSPCheckDot()
1236: @*/
1237: PetscErrorCode KSPCheckSolve(KSP ksp,PC pc,Vec vec)
1238: {
1239:   PetscErrorCode     ierr;
1240:   PCFailedReason     pcreason;
1241:   PC                 subpc;

1244:   KSPGetPC(ksp,&subpc);
1245:   PCGetFailedReason(subpc,&pcreason);
1246:   if (pcreason || (ksp->reason < 0 && ksp->reason != KSP_DIVERGED_ITS)) {
1247:     if (pc->erroriffailure) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_NOT_CONVERGED,"Detected not converged in KSP inner solve: KSP reason %s PC reason %s",KSPConvergedReasons[ksp->reason],PCFailedReasons[pcreason]);
1248:     else {
1249:       PetscInfo2(ksp,"Detected not converged in KSP inner solve: KSP reason %s PC reason %s\n",KSPConvergedReasons[ksp->reason],PCFailedReasons[pcreason]);
1250:       pc->failedreason = PC_SUBPC_ERROR;
1251:       VecSetInf(vec);
1252:     }
1253:   }
1254:   return(0);
1255: }
1256: