Actual source code: iterativ.c

petsc-3.12.5 2020-03-29
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 truly 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:   KSP            kspinner;
405:   PetscBool      flg;

408:   KSPGetPC(ksp, &pc);

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

420:   /* if pc is ksp */
421:   PetscObjectTypeCompare((PetscObject)pc,PCKSP,&flg);
422:   if (flg) {
423:     PCKSPGetKSP(pc, &kspinner);
424:     KSPSetTolerances(kspinner, inner_rtol, outer_abstol, outer_dtol, outer_maxits);
425:     return(0);
426:   }

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

440:   /* if pc is deflation*/
441:   PetscObjectTypeCompare((PetscObject)pc,PCDEFLATION,&flg);
442:   if (flg) {
443:     PCDeflationGetCoarseKSP(pc,&kspinner);
444:     KSPSetTolerances(kspinner,inner_rtol,outer_abstol,outer_dtol,PETSC_DEFAULT);
445:     return(0);
446:   }

448:   /* todo: dynamic tolerance may apply to other types of pc too */
449:   return(0);
450: }

452: /*
453:   Destroy the dummy context used in KSPMonitorDynamicTolerance()
454: */
455: PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy)
456: {

460:   PetscFree(*dummy);
461:   return(0);
462: }

464: /*
465:   Default (short) KSP Monitor, same as KSPMonitorDefault() except
466:   it prints fewer digits of the residual as the residual gets smaller.
467:   This is because the later digits are meaningless and are often
468:   different on different machines; by using this routine different
469:   machines will usually generate the same output.

471:   Deprecated: Intentionally has no manual page
472: */
473: PetscErrorCode  KSPMonitorDefaultShort(KSP ksp,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *dummy)
474: {
476:   PetscViewer    viewer = dummy->viewer;

480:   PetscViewerPushFormat(viewer,dummy->format);
481:   PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
482:   if (its == 0 && ((PetscObject)ksp)->prefix) {
483:     PetscViewerASCIIPrintf(viewer,"  Residual norms for %s solve.\n",((PetscObject)ksp)->prefix);
484:   }

486:   if (fnorm > 1.e-9) {
487:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %g \n",its,(double)fnorm);
488:   } else if (fnorm > 1.e-11) {
489:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm %5.3e \n",its,(double)fnorm);
490:   } else {
491:     PetscViewerASCIIPrintf(viewer,"%3D KSP Residual norm < 1.e-11\n",its);
492:   }
493:   PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
494:   PetscViewerPopFormat(viewer);
495:   return(0);
496: }

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

502:    Collective on ksp

504:    Input Parameters:
505: +  ksp   - iterative context
506: .  n     - iteration number
507: .  rnorm - 2-norm residual value (may be estimated)
508: -  dummy - unused convergence context

510:    Returns:
511: .  reason - KSP_CONVERGED_ITERATING, KSP_CONVERGED_ITS

513:    Notes:
514:    This should be used as the convergence test with the option
515:    KSPSetNormType(ksp,KSP_NORM_NONE), since norms of the residual are
516:    not computed. Convergence is then declared after the maximum number
517:    of iterations have been reached. Useful when one is using CG or
518:    BiCGStab as a smoother.

520:    Level: advanced

522: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPSetNormType()
523: @*/
524: PetscErrorCode  KSPConvergedSkip(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *dummy)
525: {
529:   *reason = KSP_CONVERGED_ITERATING;
530:   if (n >= ksp->max_it) *reason = KSP_CONVERGED_ITS;
531:   return(0);
532: }


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

538:    Note Collective

540:    Output Parameter:
541: .  ctx - convergence context

543:    Level: intermediate

545: .seealso: KSPConvergedDefault(), KSPConvergedDefaultDestroy(), KSPSetConvergenceTest(), KSPSetTolerances(),
546:           KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm()
547: @*/
548: PetscErrorCode  KSPConvergedDefaultCreate(void **ctx)
549: {
550:   PetscErrorCode         ierr;
551:   KSPConvergedDefaultCtx *cctx;

554:   PetscNew(&cctx);
555:   *ctx = cctx;
556:   return(0);
557: }

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

564:    Collective on ksp

566:    Input Parameters:
567: .  ksp   - iterative context

569:    Options Database:
570: .   -ksp_converged_use_initial_residual_norm

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

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

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

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


583:    Level: intermediate

585: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUMIRNorm()
586: @*/
587: PetscErrorCode  KSPConvergedDefaultSetUIRNorm(KSP ksp)
588: {
589:   KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx*) ksp->cnvP;

593:   if (ksp->converged != KSPConvergedDefault) return(0);
594:   if (ctx->mininitialrtol) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
595:   ctx->initialrtol = PETSC_TRUE;
596:   return(0);
597: }

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

604:    Collective on ksp

606:    Input Parameters:
607: .  ksp   - iterative context

609:    Options Database:
610: .   -ksp_converged_use_min_initial_residual_norm

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

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

617:    Level: intermediate

619: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUIRNorm()
620: @*/
621: PetscErrorCode  KSPConvergedDefaultSetUMIRNorm(KSP ksp)
622: {
623:   KSPConvergedDefaultCtx *ctx = (KSPConvergedDefaultCtx*) ksp->cnvP;

627:   if (ksp->converged != KSPConvergedDefault) return(0);
628:   if (ctx->initialrtol) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"Cannot use KSPConvergedDefaultSetUIRNorm() and KSPConvergedDefaultSetUMIRNorm() together");
629:   ctx->mininitialrtol = PETSC_TRUE;
630:   return(0);
631: }

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

636:    Collective on ksp

638:    Input Parameters:
639: +  ksp   - iterative context
640: .  n     - iteration number
641: .  rnorm - residual norm (may be estimated, depending on the method may be the preconditioned residual norm)
642: -  ctx - convergence context which must be created by KSPConvergedDefaultCreate()

644:    Output Parameter:
645: +   positive - if the iteration has converged;
646: .   negative - if residual norm exceeds divergence threshold;
647: -   0 - otherwise.

649:    Notes:
650:    KSPConvergedDefault() reaches convergence when   rnorm < MAX (rtol * rnorm_0, abstol);
651:    Divergence is detected if  rnorm > dtol * rnorm_0,

653:    where:
654: +     rtol = relative tolerance,
655: .     abstol = absolute tolerance.
656: .     dtol = divergence tolerance,
657: -     rnorm_0 is the two norm of the right hand side (or the preconditioned norm, depending on what was set with
658:           KSPSetNormType(). When initial guess is non-zero you
659:           can call KSPConvergedDefaultSetUIRNorm() to use the norm of (b - A*(initial guess))
660:           as the starting point for relative norm convergence testing, that is as rnorm_0

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

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

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

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

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

672:    Level: intermediate

674: .seealso: KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(), KSPConvergedReason, KSPGetConvergedReason(),
675:           KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm(), KSPConvergedDefaultCreate(), KSPConvergedDefaultDestroy()
676: @*/
677: PetscErrorCode  KSPConvergedDefault(KSP ksp,PetscInt n,PetscReal rnorm,KSPConvergedReason *reason,void *ctx)
678: {
679:   PetscErrorCode         ierr;
680:   KSPConvergedDefaultCtx *cctx = (KSPConvergedDefaultCtx*) ctx;
681:   KSPNormType            normtype;

686:   *reason = KSP_CONVERGED_ITERATING;

688:   KSPGetNormType(ksp,&normtype);
689:   if (normtype == KSP_NORM_NONE) return(0);

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

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

730:   if (PetscIsInfOrNanReal(rnorm)) {
731:     PCFailedReason pcreason;
732:     PetscInt       sendbuf,pcreason_max;
733:     PCGetFailedReason(ksp->pc,&pcreason);
734:     sendbuf = (PetscInt)pcreason;
735:     MPI_Allreduce(&sendbuf,&pcreason_max,1,MPIU_INT,MPIU_MAX,PetscObjectComm((PetscObject)ksp));
736:     if (pcreason_max) {
737:       *reason = KSP_DIVERGED_PC_FAILED;
738:       VecSetInf(ksp->vec_sol);
739:       PetscInfo(ksp,"Linear solver pcsetup fails, declaring divergence \n");
740:     } else {
741:       *reason = KSP_DIVERGED_NANORINF;
742:       PetscInfo(ksp,"Linear solver has created a not a number (NaN) as the residual norm, declaring divergence \n");
743:     }
744:   } else if (rnorm <= ksp->ttol) {
745:     if (rnorm < ksp->abstol) {
746:       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);
747:       *reason = KSP_CONVERGED_ATOL;
748:     } else {
749:       if (cctx->initialrtol) {
750:         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);
751:       } else {
752:         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);
753:       }
754:       *reason = KSP_CONVERGED_RTOL;
755:     }
756:   } else if (rnorm >= ksp->divtol*ksp->rnorm0) {
757:     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);
758:     *reason = KSP_DIVERGED_DTOL;
759:   }
760:   return(0);
761: }

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

766:    Not Collective

768:    Input Parameters:
769: .  ctx - convergence context

771:    Level: intermediate

773: .seealso: KSPConvergedDefault(), KSPConvergedDefaultCreate(), KSPSetConvergenceTest(), KSPSetTolerances(), KSPConvergedSkip(),
774:           KSPConvergedReason, KSPGetConvergedReason(), KSPConvergedDefaultSetUIRNorm(), KSPConvergedDefaultSetUMIRNorm()
775: @*/
776: PetscErrorCode  KSPConvergedDefaultDestroy(void *ctx)
777: {
778:   PetscErrorCode         ierr;
779:   KSPConvergedDefaultCtx *cctx = (KSPConvergedDefaultCtx*) ctx;

782:   VecDestroy(&cctx->work);
783:   PetscFree(ctx);
784:   return(0);
785: }

787: /*
788:    KSPBuildSolutionDefault - Default code to create/move the solution.

790:    Collective on ksp

792:    Input Parameters:
793: +  ksp - iterative context
794: -  v   - pointer to the user's vector

796:    Output Parameter:
797: .  V - pointer to a vector containing the solution

799:    Level: advanced

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

803: .seealso: KSPGetSolution(), KSPBuildResidualDefault()
804: */
805: PetscErrorCode KSPBuildSolutionDefault(KSP ksp,Vec v,Vec *V)
806: {

810:   if (ksp->pc_side == PC_RIGHT) {
811:     if (ksp->pc) {
812:       if (v) {
813:         KSP_PCApply(ksp,ksp->vec_sol,v); *V = v;
814:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with right preconditioner");
815:     } else {
816:       if (v) {
817:         VecCopy(ksp->vec_sol,v); *V = v;
818:       } else *V = ksp->vec_sol;
819:     }
820:   } else if (ksp->pc_side == PC_SYMMETRIC) {
821:     if (ksp->pc) {
822:       if (ksp->transpose_solve) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with symmetric preconditioner and transpose solve");
823:       if (v) {
824:         PCApplySymmetricRight(ksp->pc,ksp->vec_sol,v);
825:         *V = v;
826:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Not working with symmetric preconditioner");
827:     } else {
828:       if (v) {
829:         VecCopy(ksp->vec_sol,v); *V = v;
830:       } else *V = ksp->vec_sol;
831:     }
832:   } else {
833:     if (v) {
834:       VecCopy(ksp->vec_sol,v); *V = v;
835:     } else *V = ksp->vec_sol;
836:   }
837:   return(0);
838: }

840: /*
841:    KSPBuildResidualDefault - Default code to compute the residual.

843:    Collecive on ksp

845:    Input Parameters:
846: .  ksp - iterative context
847: .  t   - pointer to temporary vector
848: .  v   - pointer to user vector

850:    Output Parameter:
851: .  V - pointer to a vector containing the residual

853:    Level: advanced

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

857: .seealso: KSPBuildSolutionDefault()
858: */
859: PetscErrorCode KSPBuildResidualDefault(KSP ksp,Vec t,Vec v,Vec *V)
860: {
862:   Mat            Amat,Pmat;

865:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
866:   PCGetOperators(ksp->pc,&Amat,&Pmat);
867:   KSPBuildSolution(ksp,t,NULL);
868:   KSP_MatMult(ksp,Amat,t,v);
869:   VecAYPX(v,-1.0,ksp->vec_rhs);
870:   *V   = v;
871:   return(0);
872: }

874: /*@C
875:   KSPCreateVecs - Gets a number of work vectors.

877:   Collective on ksp

879:   Input Parameters:
880: + ksp  - iterative context
881: . rightn  - number of right work vectors
882: - leftn   - number of left work vectors to allocate

884:   Output Parameter:
885: +  right - the array of vectors created
886: -  left - the array of left vectors

888:    Note: The right vector has as many elements as the matrix has columns. The left
889:      vector has as many elements as the matrix has rows.

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

893:    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
894:                     that does not exist tries to get them from the DM (if it is provided).

896:    Level: advanced

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

900: @*/
901: PetscErrorCode KSPCreateVecs(KSP ksp,PetscInt rightn, Vec **right,PetscInt leftn,Vec **left)
902: {
904:   Vec            vecr = NULL,vecl = NULL;
905:   PetscBool      matset,pmatset;
906:   Mat            mat = NULL;

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

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

975:   Collective on ksp

977:   Input Parameters:
978: + ksp  - iterative context
979: - nw   - number of work vectors to allocate

981:   Level: developer

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

985: @*/
986: PetscErrorCode KSPSetWorkVecs(KSP ksp,PetscInt nw)
987: {

991:   VecDestroyVecs(ksp->nwork,&ksp->work);
992:   ksp->nwork = nw;
993:   KSPCreateVecs(ksp,nw,&ksp->work,0,NULL);
994:   PetscLogObjectParents(ksp,nw,ksp->work);
995:   return(0);
996: }

998: /*
999:   KSPDestroyDefault - Destroys a iterative context variable for methods with
1000:   no separate context.  Preferred calling sequence KSPDestroy().

1002:   Input Parameter:
1003: . ksp - the iterative context

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

1007: */
1008: PetscErrorCode KSPDestroyDefault(KSP ksp)
1009: {

1014:   PetscFree(ksp->data);
1015:   return(0);
1016: }

1018: /*@
1019:    KSPGetConvergedReason - Gets the reason the KSP iteration was stopped.

1021:    Not Collective

1023:    Input Parameter:
1024: .  ksp - the KSP context

1026:    Output Parameter:
1027: .  reason - negative value indicates diverged, positive value converged, see KSPConvergedReason

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

1043:    Options Database:
1044: .   -ksp_converged_reason - prints the reason to standard out

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

1049:    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
1050:    solvers which are used by the SNESNEWTONTR (trust region) solver.

1052:    Level: intermediate

1054: .seealso: KSPSetConvergenceTest(), KSPConvergedDefault(), KSPSetTolerances(), KSPConvergedReason
1055: @*/
1056: PetscErrorCode  KSPGetConvergedReason(KSP ksp,KSPConvergedReason *reason)
1057: {
1061:   *reason = ksp->reason;
1062:   return(0);
1063: }

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

1069:    Logically Collective on ksp

1071:    Input Parameters:
1072: +  ksp - the preconditioner context
1073: -  dm - the dm, cannot be NULL

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

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

1084:    Level: intermediate

1086: .seealso: KSPGetDM(), KSPSetDMActive(), KSPSetComputeOperators(), KSPSetComputeRHS(), KSPSetComputeInitialGuess(), DMKSPSetComputeOperators(), DMKSPSetComputeRHS(), DMKSPSetComputeInitialGuess()
1087: @*/
1088: PetscErrorCode  KSPSetDM(KSP ksp,DM dm)
1089: {
1091:   PC             pc;

1096:   PetscObjectReference((PetscObject)dm);
1097:   if (ksp->dm) {                /* Move the DMSNES context over to the new DM unless the new DM already has one */
1098:     if (ksp->dm->dmksp && !dm->dmksp) {
1099:       DMKSP kdm;
1100:       DMCopyDMKSP(ksp->dm,dm);
1101:       DMGetDMKSP(ksp->dm,&kdm);
1102:       if (kdm->originaldm == ksp->dm) kdm->originaldm = dm; /* Grant write privileges to the replacement DM */
1103:     }
1104:     DMDestroy(&ksp->dm);
1105:   }
1106:   ksp->dm       = dm;
1107:   ksp->dmAuto   = PETSC_FALSE;
1108:   KSPGetPC(ksp,&pc);
1109:   PCSetDM(pc,dm);
1110:   ksp->dmActive = PETSC_TRUE;
1111:   return(0);
1112: }

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

1117:    Logically Collective on ksp

1119:    Input Parameters:
1120: +  ksp - the preconditioner context
1121: -  flg - use the DM

1123:    Notes:
1124:    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.

1126:    Level: intermediate

1128: .seealso: KSPGetDM(), KSPSetDM(), SNESSetDM(), KSPSetComputeOperators(), KSPSetComputeRHS(), KSPSetComputeInitialGuess()
1129: @*/
1130: PetscErrorCode  KSPSetDMActive(KSP ksp,PetscBool flg)
1131: {
1135:   ksp->dmActive = flg;
1136:   return(0);
1137: }

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

1142:    Not Collective

1144:    Input Parameter:
1145: . ksp - the preconditioner context

1147:    Output Parameter:
1148: .  dm - the dm

1150:    Level: intermediate


1153: .seealso: KSPSetDM(), KSPSetDMActive()
1154: @*/
1155: PetscErrorCode  KSPGetDM(KSP ksp,DM *dm)
1156: {

1161:   if (!ksp->dm) {
1162:     DMShellCreate(PetscObjectComm((PetscObject)ksp),&ksp->dm);
1163:     ksp->dmAuto = PETSC_TRUE;
1164:   }
1165:   *dm = ksp->dm;
1166:   return(0);
1167: }

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

1172:    Logically Collective on ksp

1174:    Input Parameters:
1175: +  ksp - the KSP context
1176: -  usrP - optional user context

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

1182:    Level: intermediate

1184: .seealso: KSPGetApplicationContext()
1185: @*/
1186: PetscErrorCode  KSPSetApplicationContext(KSP ksp,void *usrP)
1187: {
1189:   PC             pc;

1193:   ksp->user = usrP;
1194:   KSPGetPC(ksp,&pc);
1195:   PCSetApplicationContext(pc,usrP);
1196:   return(0);
1197: }

1199: /*@
1200:    KSPGetApplicationContext - Gets the user-defined context for the linear solver.

1202:    Not Collective

1204:    Input Parameter:
1205: .  ksp - KSP context

1207:    Output Parameter:
1208: .  usrP - user context

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

1214:    Level: intermediate

1216: .seealso: KSPSetApplicationContext()
1217: @*/
1218: PetscErrorCode  KSPGetApplicationContext(KSP ksp,void *usrP)
1219: {
1222:   *(void**)usrP = ksp->user;
1223:   return(0);
1224: }

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

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

1232:    Collective on ksp

1234:    Input Parameter:
1235: +  ksp - the linear solver (KSP) context.
1236: .  pc - the preconditioner context
1237: -  vec - a vector that will be initialized with Inf to indicate lack of convergence

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

1241:    Level: developer

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

1245: .seealso: KSPCreate(), KSPSetType(), KSP, KSPCheckNorm(), KSPCheckDot()
1246: @*/
1247: PetscErrorCode KSPCheckSolve(KSP ksp,PC pc,Vec vec)
1248: {
1249:   PetscErrorCode     ierr;
1250:   PCFailedReason     pcreason;
1251:   PC                 subpc;

1254:   KSPGetPC(ksp,&subpc);
1255:   PCGetFailedReason(subpc,&pcreason);
1256:   if (pcreason || (ksp->reason < 0 && ksp->reason != KSP_DIVERGED_ITS)) {
1257:     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]);
1258:     else {
1259:       PetscInfo2(ksp,"Detected not converged in KSP inner solve: KSP reason %s PC reason %s\n",KSPConvergedReasons[ksp->reason],PCFailedReasons[pcreason]);
1260:       pc->failedreason = PC_SUBPC_ERROR;
1261:       VecSetInf(vec);
1262:     }
1263:   }
1264:   return(0);
1265: }
1266: