Actual source code: itfunc.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2: /*
  3:       Interface KSP routines that the user calls.
  4: */

  6:  #include <petsc/private/kspimpl.h>
  7:  #include <petscdm.h>

  9: PETSC_STATIC_INLINE PetscErrorCode ObjectView(PetscObject obj, PetscViewer viewer, PetscViewerFormat format)
 10: {

 13:   PetscViewerPushFormat(viewer, format);
 14:   PetscObjectView(obj, viewer);
 15:   PetscViewerPopFormat(viewer);
 16:   return(0);
 17: }

 19: /*@
 20:    KSPComputeExtremeSingularValues - Computes the extreme singular values
 21:    for the preconditioned operator. Called after or during KSPSolve().

 23:    Not Collective

 25:    Input Parameter:
 26: .  ksp - iterative context obtained from KSPCreate()

 28:    Output Parameters:
 29: .  emin, emax - extreme singular values

 31:    Options Database Keys:
 32: .  -ksp_compute_singularvalues - compute extreme singular values and print when KSPSolve completes.

 34:    Notes:
 35:    One must call KSPSetComputeSingularValues() before calling KSPSetUp()
 36:    (or use the option -ksp_compute_eigenvalues) in order for this routine to work correctly.

 38:    Many users may just want to use the monitoring routine
 39:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
 40:    to print the extreme singular values at each iteration of the linear solve.

 42:    Estimates of the smallest singular value may be very inaccurate, especially if the Krylov method has not converged.
 43:    The largest singular value is usually accurate to within a few percent if the method has converged, but is still not
 44:    intended for eigenanalysis.

 46:    Disable restarts if using KSPGMRES, otherwise this estimate will only be using those iterations after the last
 47:    restart. See KSPGMRESSetRestart() for more details.

 49:    Level: advanced

 51: .keywords: compute, extreme, singular, values

 53: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeEigenvalues(), KSP
 54: @*/
 55: PetscErrorCode  KSPComputeExtremeSingularValues(KSP ksp,PetscReal *emax,PetscReal *emin)
 56: {

 63:   if (!ksp->calc_sings) SETERRQ(PetscObjectComm((PetscObject)ksp),4,"Singular values not requested before KSPSetUp()");

 65:   if (ksp->ops->computeextremesingularvalues) {
 66:     (*ksp->ops->computeextremesingularvalues)(ksp,emax,emin);
 67:   } else {
 68:     *emin = -1.0;
 69:     *emax = -1.0;
 70:   }
 71:   return(0);
 72: }

 74: /*@
 75:    KSPComputeEigenvalues - Computes the extreme eigenvalues for the
 76:    preconditioned operator. Called after or during KSPSolve().

 78:    Not Collective

 80:    Input Parameter:
 81: +  ksp - iterative context obtained from KSPCreate()
 82: -  n - size of arrays r and c. The number of eigenvalues computed (neig) will, in
 83:        general, be less than this.

 85:    Output Parameters:
 86: +  r - real part of computed eigenvalues, provided by user with a dimension of at least n
 87: .  c - complex part of computed eigenvalues, provided by user with a dimension of at least n
 88: -  neig - actual number of eigenvalues computed (will be less than or equal to n)

 90:    Options Database Keys:
 91: +  -ksp_compute_eigenvalues - Prints eigenvalues to stdout
 92: -  -ksp_plot_eigenvalues - Plots eigenvalues in an x-window display

 94:    Notes:
 95:    The number of eigenvalues estimated depends on the size of the Krylov space
 96:    generated during the KSPSolve() ; for example, with
 97:    CG it corresponds to the number of CG iterations, for GMRES it is the number
 98:    of GMRES iterations SINCE the last restart. Any extra space in r[] and c[]
 99:    will be ignored.

101:    KSPComputeEigenvalues() does not usually provide accurate estimates; it is
102:    intended only for assistance in understanding the convergence of iterative
103:    methods, not for eigenanalysis. For accurate computation of eigenvalues we recommend using
104:    the excellent package SLEPc.

106:    One must call KSPSetComputeEigenvalues() before calling KSPSetUp()
107:    in order for this routine to work correctly.

109:    Many users may just want to use the monitoring routine
110:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
111:    to print the singular values at each iteration of the linear solve.

113:    Level: advanced

115: .keywords: compute, extreme, singular, values

117: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeExtremeSingularValues(), KSP
118: @*/
119: PetscErrorCode  KSPComputeEigenvalues(KSP ksp,PetscInt n,PetscReal r[],PetscReal c[],PetscInt *neig)
120: {

127:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Requested < 0 Eigenvalues");
129:   if (!ksp->calc_sings) SETERRQ(PetscObjectComm((PetscObject)ksp),4,"Eigenvalues not requested before KSPSetUp()");

131:   if (n && ksp->ops->computeeigenvalues) {
132:     (*ksp->ops->computeeigenvalues)(ksp,n,r,c,neig);
133:   } else {
134:     *neig = 0;
135:   }
136:   return(0);
137: }

139: /*@
140:    KSPComputeRitz - Computes the Ritz or harmonic Ritz pairs associated to the
141:    smallest or largest in modulus, for the preconditioned operator.
142:    Called after KSPSolve().

144:    Not Collective

146:    Input Parameter:
147: +  ksp   - iterative context obtained from KSPCreate()
148: .  ritz  - PETSC_TRUE or PETSC_FALSE for ritz pairs or harmonic Ritz pairs, respectively
149: .  small - PETSC_TRUE or PETSC_FALSE for smallest or largest (harmonic) Ritz values, respectively
150: .  nrit  - number of (harmonic) Ritz pairs to compute

152:    Output Parameters:
153: +  nrit  - actual number of computed (harmonic) Ritz pairs 
154: .  S     - multidimensional vector with Ritz vectors
155: .  tetar - real part of the Ritz values        
156: .  tetai - imaginary part of the Ritz values

158:    Notes:
159:    -For GMRES, the (harmonic) Ritz pairs are computed from the Hessenberg matrix obtained during 
160:    the last complete cycle, or obtained at the end of the solution if the method is stopped before 
161:    a restart. Then, the number of actual (harmonic) Ritz pairs computed is less or equal to the restart
162:    parameter for GMRES if a complete cycle has been performed or less or equal to the number of GMRES 
163:    iterations.
164:    -Moreover, for real matrices, the (harmonic) Ritz pairs are possibly complex-valued. In such a case,
165:    the routine selects the complex (harmonic) Ritz value and its conjugate, and two successive columns of S 
166:    are equal to the real and the imaginary parts of the associated vectors. 
167:    -the (harmonic) Ritz pairs are given in order of increasing (harmonic) Ritz values in modulus
168:    -this is currently not implemented when PETSc is built with complex numbers

170:    One must call KSPSetComputeRitz() before calling KSPSetUp()
171:    in order for this routine to work correctly.

173:    Level: advanced

175: .keywords: compute, ritz, values

177: .seealso: KSPSetComputeRitz(), KSP
178: @*/
179: PetscErrorCode  KSPComputeRitz(KSP ksp,PetscBool ritz,PetscBool small,PetscInt *nrit,Vec S[],PetscReal tetar[],PetscReal tetai[])
180: {

185:   if (!ksp->calc_ritz) SETERRQ(PetscObjectComm((PetscObject)ksp),4,"Ritz pairs not requested before KSPSetUp()");
186:   if (ksp->ops->computeritz) {(*ksp->ops->computeritz)(ksp,ritz,small,nrit,S,tetar,tetai);}
187:   return(0);
188: }
189: /*@
190:    KSPSetUpOnBlocks - Sets up the preconditioner for each block in
191:    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
192:    methods.

194:    Collective on KSP

196:    Input Parameter:
197: .  ksp - the KSP context

199:    Notes:
200:    KSPSetUpOnBlocks() is a routine that the user can optinally call for
201:    more precise profiling (via -log_view) of the setup phase for these
202:    block preconditioners.  If the user does not call KSPSetUpOnBlocks(),
203:    it will automatically be called from within KSPSolve().

205:    Calling KSPSetUpOnBlocks() is the same as calling PCSetUpOnBlocks()
206:    on the PC context within the KSP context.

208:    Level: advanced

210: .keywords: setup, blocks

212: .seealso: PCSetUpOnBlocks(), KSPSetUp(), PCSetUp(), KSP
213: @*/
214: PetscErrorCode  KSPSetUpOnBlocks(KSP ksp)
215: {
216:   PC             pc;
218:   PCFailedReason pcreason;

222:   KSPGetPC(ksp,&pc);
223:   PCSetUpOnBlocks(pc);
224:   PCGetFailedReason(pc,&pcreason);
225:   if (pcreason) {
226:     ksp->reason = KSP_DIVERGED_PC_FAILED;
227:   }
228:   return(0);
229: }

231: /*@
232:    KSPSetReusePreconditioner - reuse the current preconditioner, do not construct a new one even if the operator changes

234:    Collective on KSP

236:    Input Parameters:
237: +  ksp   - iterative context obtained from KSPCreate()
238: -  flag - PETSC_TRUE to reuse the current preconditioner

240:    Level: intermediate

242: .keywords: setup

244: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner(), KSP
245: @*/
246: PetscErrorCode  KSPSetReusePreconditioner(KSP ksp,PetscBool flag)
247: {
248:   PC             pc;

253:   KSPGetPC(ksp,&pc);
254:   PCSetReusePreconditioner(pc,flag);
255:   return(0);
256: }

258: /*@
259:    KSPSetSkipPCSetFromOptions - prevents KSPSetFromOptions() from call PCSetFromOptions(). This is used if the same PC is shared by more than one KSP so its options are not resetable for each KSP

261:    Collective on KSP

263:    Input Parameters:
264: +  ksp   - iterative context obtained from KSPCreate()
265: -  flag - PETSC_TRUE to skip calling the PCSetFromOptions()

267:    Level: intermediate

269: .keywords: setup

271: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner(), KSP
272: @*/
273: PetscErrorCode  KSPSetSkipPCSetFromOptions(KSP ksp,PetscBool flag)
274: {
277:   ksp->skippcsetfromoptions = flag;
278:   return(0);
279: }

281: /*@
282:    KSPSetUp - Sets up the internal data structures for the
283:    later use of an iterative solver.

285:    Collective on KSP

287:    Input Parameter:
288: .  ksp   - iterative context obtained from KSPCreate()

290:    Level: developer

292: .keywords: setup

294: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), KSP
295: @*/
296: PetscErrorCode KSPSetUp(KSP ksp)
297: {
299:   Mat            A,B;
300:   Mat            mat,pmat;
301:   MatNullSpace   nullsp;
302:   PCFailedReason pcreason;
303: 

307:   /* reset the convergence flag from the previous solves */
308:   ksp->reason = KSP_CONVERGED_ITERATING;

310:   if (!((PetscObject)ksp)->type_name) {
311:     KSPSetType(ksp,KSPGMRES);
312:   }
313:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);

315:   if (ksp->dmActive && !ksp->setupstage) {
316:     /* first time in so build matrix and vector data structures using DM */
317:     if (!ksp->vec_rhs) {DMCreateGlobalVector(ksp->dm,&ksp->vec_rhs);}
318:     if (!ksp->vec_sol) {DMCreateGlobalVector(ksp->dm,&ksp->vec_sol);}
319:     DMCreateMatrix(ksp->dm,&A);
320:     KSPSetOperators(ksp,A,A);
321:     PetscObjectDereference((PetscObject)A);
322:   }

324:   if (ksp->dmActive) {
325:     DMKSP kdm;
326:     DMGetDMKSP(ksp->dm,&kdm);

328:     if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
329:       /* only computes initial guess the first time through */
330:       (*kdm->ops->computeinitialguess)(ksp,ksp->vec_sol,kdm->initialguessctx);
331:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
332:     }
333:     if (kdm->ops->computerhs) {
334:       (*kdm->ops->computerhs)(ksp,ksp->vec_rhs,kdm->rhsctx);
335:     }

337:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
338:       if (kdm->ops->computeoperators) {
339:         KSPGetOperators(ksp,&A,&B);
340:         (*kdm->ops->computeoperators)(ksp,A,B,kdm->operatorsctx);
341:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(ksp,PETSC_FALSE);");
342:     }
343:   }

345:   if (ksp->setupstage == KSP_SETUP_NEWRHS) return(0);
346:   PetscLogEventBegin(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);

348:   switch (ksp->setupstage) {
349:   case KSP_SETUP_NEW:
350:     (*ksp->ops->setup)(ksp);
351:     break;
352:   case KSP_SETUP_NEWMATRIX: {   /* This should be replaced with a more general mechanism */
353:     if (ksp->setupnewmatrix) {
354:       (*ksp->ops->setup)(ksp);
355:     }
356:   } break;
357:   default: break;
358:   }

360:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
361:   PCGetOperators(ksp->pc,&mat,&pmat);
362:   /* scale the matrix if requested */
363:   if (ksp->dscale) {
364:     PetscScalar *xx;
365:     PetscInt    i,n;
366:     PetscBool   zeroflag = PETSC_FALSE;
367:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
368:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
369:       MatCreateVecs(pmat,&ksp->diagonal,0);
370:     }
371:     MatGetDiagonal(pmat,ksp->diagonal);
372:     VecGetLocalSize(ksp->diagonal,&n);
373:     VecGetArray(ksp->diagonal,&xx);
374:     for (i=0; i<n; i++) {
375:       if (xx[i] != 0.0) xx[i] = 1.0/PetscSqrtReal(PetscAbsScalar(xx[i]));
376:       else {
377:         xx[i]    = 1.0;
378:         zeroflag = PETSC_TRUE;
379:       }
380:     }
381:     VecRestoreArray(ksp->diagonal,&xx);
382:     if (zeroflag) {
383:       PetscInfo(ksp,"Zero detected in diagonal of matrix, using 1 at those locations\n");
384:     }
385:     MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
386:     if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
387:     ksp->dscalefix2 = PETSC_FALSE;
388:   }
389:   PetscLogEventEnd(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
390:   PCSetErrorIfFailure(ksp->pc,ksp->errorifnotconverged);
391:   PCSetUp(ksp->pc);
392:   PCGetFailedReason(ksp->pc,&pcreason);
393:   if (pcreason) {
394:     ksp->reason = KSP_DIVERGED_PC_FAILED;
395:   }

397:   MatGetNullSpace(mat,&nullsp);
398:   if (nullsp) {
399:     PetscBool test = PETSC_FALSE;
400:     PetscOptionsGetBool(((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_test_null_space",&test,NULL);
401:     if (test) {
402:       MatNullSpaceTest(nullsp,mat,NULL);
403:     }
404:   }
405:   ksp->setupstage = KSP_SETUP_NEWRHS;
406:   return(0);
407: }

409: static PetscErrorCode KSPReasonView_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
410: {
412:   PetscBool      isAscii;

415:   if (format != PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,format);}
416:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
417:   if (isAscii) {
418:     PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
419:     if (ksp->reason > 0) {
420:       if (((PetscObject) ksp)->prefix) {
421:         PetscViewerASCIIPrintf(viewer,"Linear %s solve converged due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
422:       } else {
423:         PetscViewerASCIIPrintf(viewer,"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
424:       }
425:     } else {
426:       if (((PetscObject) ksp)->prefix) {
427:         PetscViewerASCIIPrintf(viewer,"Linear %s solve did not converge due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
428:       } else {
429:         PetscViewerASCIIPrintf(viewer,"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
430:       }
431:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
432:         PCFailedReason reason;
433:         PCGetFailedReason(ksp->pc,&reason);
434:         PetscViewerASCIIPrintf(viewer,"               PC_FAILED due to %s \n",PCFailedReasons[reason]);
435:       }
436:     }
437:     PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
438:   }
439:   if (format != PETSC_VIEWER_DEFAULT) {PetscViewerPopFormat(viewer);}
440:   return(0);
441: }

443: /*@
444:    KSPReasonView - Displays the reason a KSP solve converged or diverged to a viewer

446:    Collective on KSP

448:    Parameter:
449: +  ksp - iterative context obtained from KSPCreate()
450: -  viewer - the viewer to display the reason


453:    Options Database Keys:
454: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations

456:    Level: beginner

458: .keywords: solve, linear system

460: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
461:           KSPSolveTranspose(), KSPGetIterationNumber(), KSP
462: @*/
463: PetscErrorCode KSPReasonView(KSP ksp,PetscViewer viewer)
464: {

468:   KSPReasonView_Internal(ksp, viewer, PETSC_VIEWER_DEFAULT);
469:   return(0);
470: }

472: #if defined(PETSC_HAVE_THREADSAFETY)
473: #define KSPReasonViewFromOptions KSPReasonViewFromOptionsUnsafe
474: #else
475: #endif
476: /*@C
477:   KSPReasonViewFromOptions - Processes command line options to determine if/how a KSPReason is to be viewed.

479:   Collective on KSP

481:   Input Parameters:
482: . ksp   - the KSP object

484:   Level: intermediate

486: @*/
487: PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
488: {
489:   PetscViewer       viewer;
490:   PetscBool         flg;
491:   PetscViewerFormat format;
492:   PetscErrorCode    ierr;

495:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_converged_reason",&viewer,&format,&flg);
496:   if (flg) {
497:     KSPReasonView_Internal(ksp, viewer, format);
498:     PetscViewerDestroy(&viewer);
499:   }
500:   return(0);
501: }

503:  #include <petscdraw.h>

505: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
506: {
507:   PetscReal     *r, *c;
508:   PetscInt       n, i, neig;
509:   PetscBool      isascii, isdraw;
510:   PetscMPIInt    rank;

514:   MPI_Comm_rank(PetscObjectComm((PetscObject) ksp), &rank);
515:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
516:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW,  &isdraw);
517:   if (isExplicit) {
518:     VecGetSize(ksp->vec_sol,&n);
519:     PetscMalloc2(n, &r, n, &c);
520:     KSPComputeEigenvaluesExplicitly(ksp, n, r, c);
521:     neig = n;
522:   } else {
523:     PetscInt nits;

525:     KSPGetIterationNumber(ksp, &nits);
526:     n    = nits+2;
527:     if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n");return(0);}
528:     PetscMalloc2(n, &r, n, &c);
529:     KSPComputeEigenvalues(ksp, n, r, c, &neig);
530:   }
531:   if (isascii) {
532:     PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively");
533:     for (i = 0; i < neig; ++i) {
534:       if (c[i] >= 0.0) {PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double) r[i],  (double) c[i]);}
535:       else             {PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double) r[i], -(double) c[i]);}
536:     }
537:   } else if (isdraw && !rank) {
538:     PetscDraw   draw;
539:     PetscDrawSP drawsp;

541:     if (format == PETSC_VIEWER_DRAW_CONTOUR) {
542:       KSPPlotEigenContours_Private(ksp,neig,r,c);
543:     } else {
544:       if (!ksp->eigviewer) {PetscViewerDrawOpen(PETSC_COMM_SELF,0,isExplicit ? "Explicitly Computed Eigenvalues" : "Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);}
545:       PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
546:       PetscDrawSPCreate(draw,1,&drawsp);
547:       PetscDrawSPReset(drawsp);
548:       for (i = 0; i < neig; ++i) {PetscDrawSPAddPoint(drawsp,r+i,c+i);}
549:       PetscDrawSPDraw(drawsp,PETSC_TRUE);
550:       PetscDrawSPSave(drawsp);
551:       PetscDrawSPDestroy(&drawsp);
552:     }
553:   }
554:   PetscFree2(r, c);
555:   return(0);
556: }

558: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
559: {
560:   PetscReal      smax, smin;
561:   PetscInt       nits;
562:   PetscBool      isascii;

566:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
567:   KSPGetIterationNumber(ksp, &nits);
568:   if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n");return(0);}
569:   KSPComputeExtremeSingularValues(ksp, &smax, &smin);
570:   if (isascii) {PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)smax,(double)smin,(double)(smax/smin));}
571:   return(0);
572: }

574: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
575: {
576:   PetscBool      isascii;

580:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
581:   if (ksp->dscale && !ksp->dscalefix) SETERRQ(PetscObjectComm((PetscObject) ksp), PETSC_ERR_ARG_WRONGSTATE, "Cannot compute final scale with -ksp_diagonal_scale except also with -ksp_diagonal_scale_fix");
582:   if (isascii) {
583:     Mat       A;
584:     Vec       t;
585:     PetscReal norm;

587:     PCGetOperators(ksp->pc, &A, NULL);
588:     VecDuplicate(ksp->vec_rhs, &t);
589:     KSP_MatMult(ksp, A, ksp->vec_sol, t);
590:     VecAYPX(t, -1.0, ksp->vec_rhs);
591:     VecNorm(t, NORM_2, &norm);
592:     VecDestroy(&t);
593:     PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double) norm);
594:   }
595:   return(0);
596: }

598: /*@
599:    KSPSolve - Solves linear system.

601:    Collective on KSP

603:    Parameter:
604: +  ksp - iterative context obtained from KSPCreate()
605: .  b - the right hand side vector
606: -  x - the solution  (this may be the same vector as b, then b will be overwritten with answer)

608:    Options Database Keys:
609: +  -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
610: .  -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
611: .  -ksp_plot_eigencontours - plot the computed eigenvalues in an X-window with contours
612: .  -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and using LAPACK
613: .  -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
614: .  -ksp_view_mat binary - save matrix to the default binary viewer
615: .  -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
616: .  -ksp_view_rhs binary - save right hand side vector to the default binary viewer
617: .  -ksp_view_solution binary - save computed solution vector to the default binary viewer
618:            (can be read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
619: .  -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
620: .  -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
621: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
622: .  -ksp_final_residual - print 2-norm of true linear system residual at the end of the solution process
623: -  -ksp_view - print the ksp data structure at the end of the system solution

625:    Notes:

627:    If one uses KSPSetDM() then x or b need not be passed. Use KSPGetSolution() to access the solution in this case.

629:    The operator is specified with KSPSetOperators().

631:    Call KSPGetConvergedReason() to determine if the solver converged or failed and
632:    why. The number of iterations can be obtained from KSPGetIterationNumber().

634:    If you provide a matrix that has a MatSetNullSpace() and MatSetTransposeNullSpace() this will use that information to solve singular systems
635:    in the least squares sense with a norm minimizing solution.
636: $
637: $                   A x = b   where b = b_p + b_t where b_t is not in the range of A (and hence by the fundamental theorem of linear algebra is in the nullspace(A') see MatSetNullSpace()
638: $
639: $    KSP first removes b_t producing the linear system  A x = b_p (which has multiple solutions) and solves this to find the ||x|| minimizing solution (and hence
640: $    it finds the solution x orthogonal to the nullspace(A). The algorithm is simply in each iteration of the Krylov method we remove the nullspace(A) from the search
641: $    direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
642: $
643: $    We recommend always using GMRES for such singular systems.
644: $    If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
645: $    If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).

647:    Developer Note: The reason we cannot always solve  nullspace(A) != nullspace(A') systems with right preconditioning is because we need to remove at each iteration
648:        the nullspace(AB) from the search direction. While we know the nullspace(A) the nullspace(AB) equals B^-1 times the nullspace(A) but except for trivial preconditioners
649:        such as diagonal scaling we cannot apply the inverse of the preconditioner to a vector and thus cannot compute the nullspace(AB).


652:    If using a direct method (e.g., via the KSP solver
653:    KSPPREONLY and a preconditioner such as PCLU/PCILU),
654:    then its=1.  See KSPSetTolerances() and KSPConvergedDefault()
655:    for more details.

657:    Understanding Convergence:
658:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
659:    KSPComputeEigenvaluesExplicitly() provide information on additional
660:    options to monitor convergence and print eigenvalue information.

662:    Level: beginner

664: .keywords: solve, linear system

666: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
667:           KSPSolveTranspose(), KSPGetIterationNumber(), MatNullSpaceCreate(), MatSetNullSpace(), MatSetTransposeNullSpace(), KSP
668: @*/
669: PetscErrorCode KSPSolve(KSP ksp,Vec b,Vec x)
670: {
671:   PetscErrorCode    ierr;
672:   PetscBool         flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
673:   Mat               mat,pmat;
674:   MPI_Comm          comm;
675:   MatNullSpace      nullsp;
676:   Vec               btmp,vec_rhs=0;

682:   comm = PetscObjectComm((PetscObject)ksp);
683:   if (x && x == b) {
684:     if (!ksp->guess_zero) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
685:     VecDuplicate(b,&x);
686:     inXisinB = PETSC_TRUE;
687:   }
688:   if (b) {
689:     PetscObjectReference((PetscObject)b);
690:     VecDestroy(&ksp->vec_rhs);
691:     ksp->vec_rhs = b;
692:   }
693:   if (x) {
694:     PetscObjectReference((PetscObject)x);
695:     VecDestroy(&ksp->vec_sol);
696:     ksp->vec_sol = x;
697:   }
698:   if (ksp->viewPre) {ObjectView((PetscObject) ksp, ksp->viewerPre, ksp->formatPre);}

700:   ksp->transpose_solve = PETSC_FALSE;

702:   if (ksp->presolve) {(*ksp->presolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->prectx);}

704:   /* reset the residual history list if requested */
705:   if (ksp->res_hist_reset) ksp->res_hist_len = 0;

707:   PetscLogEventBegin(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);

709:   if (ksp->guess) {
710:     PetscObjectState ostate,state;

712:     KSPGuessSetUp(ksp->guess);
713:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
714:     KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
715:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
716:     if (state != ostate) {
717:       ksp->guess_zero = PETSC_FALSE;
718:     } else {
719:       PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
720:       ksp->guess_zero = PETSC_TRUE;
721:     }
722:   }

724:   /* KSPSetUp() scales the matrix if needed */
725:   KSPSetUp(ksp);
726:   KSPSetUpOnBlocks(ksp);

728:   VecSetErrorIfLocked(ksp->vec_sol,3);

730:   PCGetOperators(ksp->pc,&mat,&pmat);
731:   /* diagonal scale RHS if called for */
732:   if (ksp->dscale) {
733:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
734:     /* second time in, but matrix was scaled back to original */
735:     if (ksp->dscalefix && ksp->dscalefix2) {
736:       Mat mat,pmat;

738:       PCGetOperators(ksp->pc,&mat,&pmat);
739:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
740:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
741:     }

743:     /* scale initial guess */
744:     if (!ksp->guess_zero) {
745:       if (!ksp->truediagonal) {
746:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
747:         VecCopy(ksp->diagonal,ksp->truediagonal);
748:         VecReciprocal(ksp->truediagonal);
749:       }
750:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
751:     }
752:   }
753:   PCPreSolve(ksp->pc,ksp);

755:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
756:   if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
757:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
758:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
759:     ksp->guess_zero = PETSC_FALSE;
760:   }

762:   /* can we mark the initial guess as zero for this solve? */
763:   guess_zero = ksp->guess_zero;
764:   if (!ksp->guess_zero) {
765:     PetscReal norm;

767:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
768:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
769:   }
770:   MatGetTransposeNullSpace(pmat,&nullsp);
771:   if (nullsp) {
772:     VecDuplicate(ksp->vec_rhs,&btmp);
773:     VecCopy(ksp->vec_rhs,btmp);
774:     MatNullSpaceRemove(nullsp,btmp);
775:     vec_rhs      = ksp->vec_rhs;
776:     ksp->vec_rhs = btmp;
777:   }
778:   VecLockReadPush(ksp->vec_rhs);
779:   if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
780:     VecSetInf(ksp->vec_sol);
781:   }
782:   (*ksp->ops->solve)(ksp);

784:   VecLockReadPop(ksp->vec_rhs);
785:   if (nullsp) {
786:     ksp->vec_rhs = vec_rhs;
787:     VecDestroy(&btmp);
788:   }

790:   ksp->guess_zero = guess_zero;

792:   if (!ksp->reason) SETERRQ(comm,PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
793:   ksp->totalits += ksp->its;

795:   if (ksp->viewReason) {KSPReasonView_Internal(ksp, ksp->viewerReason, ksp->formatReason);}
796:   PCPostSolve(ksp->pc,ksp);

798:   /* diagonal scale solution if called for */
799:   if (ksp->dscale) {
800:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
801:     /* unscale right hand side and matrix */
802:     if (ksp->dscalefix) {
803:       Mat mat,pmat;

805:       VecReciprocal(ksp->diagonal);
806:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
807:       PCGetOperators(ksp->pc,&mat,&pmat);
808:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
809:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
810:       VecReciprocal(ksp->diagonal);
811:       ksp->dscalefix2 = PETSC_TRUE;
812:     }
813:   }
814:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
815:   if (ksp->guess) {
816:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
817:   }
818:   if (ksp->postsolve) {
819:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
820:   }

822:   PCGetOperators(ksp->pc,&mat,&pmat);
823:   if (ksp->viewEV)       {KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV,    ksp->formatEV);}
824:   if (ksp->viewEVExp)    {KSPViewEigenvalues_Internal(ksp, PETSC_TRUE,  ksp->viewerEVExp, ksp->formatEVExp);}
825:   if (ksp->viewSV)       {KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV);}
826:   if (ksp->viewFinalRes) {KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes);}
827:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,           ksp->viewerMat,    ksp->formatMat);}
828:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,          ksp->viewerPMat,   ksp->formatPMat);}
829:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs,  ksp->viewerRhs,    ksp->formatRhs);}
830:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol,  ksp->viewerSol,    ksp->formatSol);}
831:   if (ksp->view)         {ObjectView((PetscObject) ksp,           ksp->viewer,       ksp->format);}
832:   if (ksp->viewDScale)   {ObjectView((PetscObject) ksp->diagonal, ksp->viewerDScale, ksp->formatDScale);}
833:   if (ksp->viewMatExp)   {
834:     Mat A, B;

836:     PCGetOperators(ksp->pc, &A, NULL);
837:     MatComputeExplicitOperator(A, &B);
838:     ObjectView((PetscObject) B, ksp->viewerMatExp, ksp->formatMatExp);
839:     MatDestroy(&B);
840:   }
841:   if (ksp->viewPOpExp)   {
842:     Mat B;

844:     KSPComputeExplicitOperator(ksp, &B);
845:     ObjectView((PetscObject) B, ksp->viewerPOpExp, ksp->formatPOpExp);
846:     MatDestroy(&B);
847:   }

849:   if (inXisinB) {
850:     VecCopy(x,b);
851:     VecDestroy(&x);
852:   }
853:   PetscObjectSAWsBlock((PetscObject)ksp);
854:   if (ksp->errorifnotconverged && ksp->reason < 0 && ksp->reason != KSP_DIVERGED_ITS) SETERRQ1(comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged, reason %s",KSPConvergedReasons[ksp->reason]);
855:   return(0);
856: }

858: /*@
859:    KSPSolveTranspose - Solves the transpose of a linear system.

861:    Collective on KSP

863:    Input Parameter:
864: +  ksp - iterative context obtained from KSPCreate()
865: .  b - right hand side vector
866: -  x - solution vector

868:    Notes:
869:     For complex numbers this solve the non-Hermitian transpose system.

871:    This currently does NOT correctly use the null space of the operator and its transpose for solving singular systems.

873:    Developer Notes:
874:     We need to implement a KSPSolveHermitianTranspose()

876:    Level: developer

878: .keywords: solve, linear system

880: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
881:           KSPSolve(), KSP
882: @*/

884: PetscErrorCode  KSPSolveTranspose(KSP ksp,Vec b,Vec x)
885: {
887:   PetscBool      inXisinB=PETSC_FALSE;
888:   Vec            vec_rhs = 0,btmp;
889:   Mat            mat,pmat;
890:   MatNullSpace   nullsp;

896:   if (x == b) {
897:     VecDuplicate(b,&x);
898:     inXisinB = PETSC_TRUE;
899:   }
900:   PetscObjectReference((PetscObject)b);
901:   PetscObjectReference((PetscObject)x);
902:   VecDestroy(&ksp->vec_rhs);
903:   VecDestroy(&ksp->vec_sol);

905:   ksp->vec_rhs         = b;
906:   ksp->vec_sol         = x;
907:   ksp->transpose_solve = PETSC_TRUE;

909:   if (ksp->presolve) {(*ksp->presolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->prectx);}

911:   if (ksp->guess) {
912:     PetscObjectState ostate,state;

914:     KSPGuessSetUp(ksp->guess);
915:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
916:     KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
917:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
918:     if (state != ostate) {
919:       ksp->guess_zero = PETSC_FALSE;
920:     } else {
921:       PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
922:       ksp->guess_zero = PETSC_TRUE;
923:     }
924:   }

926:   KSPSetUp(ksp);
927:   KSPSetUpOnBlocks(ksp);
928:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}

930:   PCGetOperators(ksp->pc,&mat,&pmat);
931:   MatGetNullSpace(pmat,&nullsp);
932:   if (nullsp) {
933:     VecDuplicate(ksp->vec_rhs,&btmp);
934:     VecCopy(ksp->vec_rhs,btmp);
935:     MatNullSpaceRemove(nullsp,btmp);
936:     vec_rhs      = ksp->vec_rhs;
937:     ksp->vec_rhs = btmp;
938:   }

940:   (*ksp->ops->solve)(ksp);
941:   ksp->totalits += ksp->its;
942:   if (nullsp) {
943:     ksp->vec_rhs = vec_rhs;
944:     VecDestroy(&btmp);
945:   }
946:   if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
947:   if (ksp->viewReason) {KSPReasonView_Internal(ksp, ksp->viewerReason, ksp->formatReason);}
948:   if (ksp->guess) {
949:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
950:   }
951:   if (ksp->postsolve) {
952:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
953:   }

955:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,          ksp->viewerMat,  ksp->formatMat);}
956:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,         ksp->viewerPMat, ksp->formatPMat);}
957:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs, ksp->viewerRhs,  ksp->formatRhs);}
958:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol, ksp->viewerSol,  ksp->formatSol);}
959:   if (ksp->view)         {ObjectView((PetscObject) ksp,          ksp->viewer,     ksp->format);}

961:   if (inXisinB) {
962:     VecCopy(x,b);
963:     VecDestroy(&x);
964:   }
965:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
966:   return(0);
967: }

969: /*@
970:    KSPResetViewers - Resets all the viewers set from the options database during KSPSetFromOptions()

972:    Collective on KSP

974:    Input Parameter:
975: .  ksp - iterative context obtained from KSPCreate()

977:    Level: beginner

979: .keywords: destroy

981: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSPSetFromOptions(), KSP
982: @*/
983: PetscErrorCode  KSPResetViewers(KSP ksp)
984: {

989:   if (!ksp) return(0);
990:   PetscViewerDestroy(&ksp->viewer);
991:   PetscViewerDestroy(&ksp->viewerPre);
992:   PetscViewerDestroy(&ksp->viewerReason);
993:   PetscViewerDestroy(&ksp->viewerMat);
994:   PetscViewerDestroy(&ksp->viewerPMat);
995:   PetscViewerDestroy(&ksp->viewerRhs);
996:   PetscViewerDestroy(&ksp->viewerSol);
997:   PetscViewerDestroy(&ksp->viewerMatExp);
998:   PetscViewerDestroy(&ksp->viewerEV);
999:   PetscViewerDestroy(&ksp->viewerSV);
1000:   PetscViewerDestroy(&ksp->viewerEVExp);
1001:   PetscViewerDestroy(&ksp->viewerFinalRes);
1002:   PetscViewerDestroy(&ksp->viewerPOpExp);
1003:   PetscViewerDestroy(&ksp->viewerDScale);
1004:   ksp->view         = PETSC_FALSE;
1005:   ksp->viewPre      = PETSC_FALSE;
1006:   ksp->viewReason   = PETSC_FALSE;
1007:   ksp->viewMat      = PETSC_FALSE;
1008:   ksp->viewPMat     = PETSC_FALSE;
1009:   ksp->viewRhs      = PETSC_FALSE;
1010:   ksp->viewSol      = PETSC_FALSE;
1011:   ksp->viewMatExp   = PETSC_FALSE;
1012:   ksp->viewEV       = PETSC_FALSE;
1013:   ksp->viewSV       = PETSC_FALSE;
1014:   ksp->viewEVExp    = PETSC_FALSE;
1015:   ksp->viewFinalRes = PETSC_FALSE;
1016:   ksp->viewPOpExp   = PETSC_FALSE;
1017:   ksp->viewDScale   = PETSC_FALSE;
1018:   return(0);
1019: }

1021: /*@
1022:    KSPReset - Resets a KSP context to the kspsetupcalled = 0 state and removes any allocated Vecs and Mats

1024:    Collective on KSP

1026:    Input Parameter:
1027: .  ksp - iterative context obtained from KSPCreate()

1029:    Level: beginner

1031: .keywords: destroy

1033: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1034: @*/
1035: PetscErrorCode  KSPReset(KSP ksp)
1036: {

1041:   if (!ksp) return(0);
1042:   if (ksp->ops->reset) {
1043:     (*ksp->ops->reset)(ksp);
1044:   }
1045:   if (ksp->pc) {PCReset(ksp->pc);}
1046:   if (ksp->guess) {
1047:     KSPGuess guess = ksp->guess;
1048:     if (guess->ops->reset) { (*guess->ops->reset)(guess); }
1049:   }
1050:   VecDestroyVecs(ksp->nwork,&ksp->work);
1051:   VecDestroy(&ksp->vec_rhs);
1052:   VecDestroy(&ksp->vec_sol);
1053:   VecDestroy(&ksp->diagonal);
1054:   VecDestroy(&ksp->truediagonal);

1056:   KSPResetViewers(ksp);

1058:   ksp->setupstage = KSP_SETUP_NEW;
1059:   return(0);
1060: }

1062: /*@
1063:    KSPDestroy - Destroys KSP context.

1065:    Collective on KSP

1067:    Input Parameter:
1068: .  ksp - iterative context obtained from KSPCreate()

1070:    Level: beginner

1072: .keywords: destroy

1074: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1075: @*/
1076: PetscErrorCode  KSPDestroy(KSP *ksp)
1077: {
1079:   PC             pc;

1082:   if (!*ksp) return(0);
1084:   if (--((PetscObject)(*ksp))->refct > 0) {*ksp = 0; return(0);}

1086:   PetscObjectSAWsViewOff((PetscObject)*ksp);

1088:   /*
1089:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1090:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1091:    refcount (and may be shared, e.g., by other ksps).
1092:    */
1093:   pc         = (*ksp)->pc;
1094:   (*ksp)->pc = NULL;
1095:   KSPReset((*ksp));
1096:   (*ksp)->pc = pc;
1097:   if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

1099:   KSPGuessDestroy(&(*ksp)->guess);
1100:   DMDestroy(&(*ksp)->dm);
1101:   PCDestroy(&(*ksp)->pc);
1102:   PetscFree((*ksp)->res_hist_alloc);
1103:   if ((*ksp)->convergeddestroy) {
1104:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
1105:   }
1106:   KSPMonitorCancel((*ksp));
1107:   PetscViewerDestroy(&(*ksp)->eigviewer);
1108:   PetscHeaderDestroy(ksp);
1109:   return(0);
1110: }

1112: /*@
1113:     KSPSetPCSide - Sets the preconditioning side.

1115:     Logically Collective on KSP

1117:     Input Parameter:
1118: .   ksp - iterative context obtained from KSPCreate()

1120:     Output Parameter:
1121: .   side - the preconditioning side, where side is one of
1122: .vb
1123:       PC_LEFT - left preconditioning (default)
1124:       PC_RIGHT - right preconditioning
1125:       PC_SYMMETRIC - symmetric preconditioning
1126: .ve

1128:     Options Database Keys:
1129: .   -ksp_pc_side <right,left,symmetric>

1131:     Notes:
1132:     Left preconditioning is used by default for most Krylov methods except KSPFGMRES which only supports right preconditioning.

1134:     For methods changing the side of the preconditioner changes the norm type that is used, see KSPSetNormType().

1136:     Symmetric preconditioning is currently available only for the KSPQCG method. Note, however, that
1137:     symmetric preconditioning can be emulated by using either right or left
1138:     preconditioning and a pre or post processing step.

1140:     Setting the PC side often affects the default norm type.  See KSPSetNormType() for details.

1142:     Level: intermediate

1144: .keywords: set, right, left, symmetric, side, preconditioner, flag

1146: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType(), KSP
1147: @*/
1148: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
1149: {
1153:   ksp->pc_side = ksp->pc_side_set = side;
1154:   return(0);
1155: }

1157: /*@
1158:     KSPGetPCSide - Gets the preconditioning side.

1160:     Not Collective

1162:     Input Parameter:
1163: .   ksp - iterative context obtained from KSPCreate()

1165:     Output Parameter:
1166: .   side - the preconditioning side, where side is one of
1167: .vb
1168:       PC_LEFT - left preconditioning (default)
1169:       PC_RIGHT - right preconditioning
1170:       PC_SYMMETRIC - symmetric preconditioning
1171: .ve

1173:     Level: intermediate

1175: .keywords: get, right, left, symmetric, side, preconditioner, flag

1177: .seealso: KSPSetPCSide(), KSP
1178: @*/
1179: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
1180: {

1186:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
1187:   *side = ksp->pc_side;
1188:   return(0);
1189: }

1191: /*@
1192:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1193:    iteration tolerances used by the default KSP convergence tests.

1195:    Not Collective

1197:    Input Parameter:
1198: .  ksp - the Krylov subspace context

1200:    Output Parameters:
1201: +  rtol - the relative convergence tolerance
1202: .  abstol - the absolute convergence tolerance
1203: .  dtol - the divergence tolerance
1204: -  maxits - maximum number of iterations

1206:    Notes:
1207:    The user can specify NULL for any parameter that is not needed.

1209:    Level: intermediate

1211: .keywords: get, tolerance, absolute, relative, divergence, convergence,
1212:            maximum, iterations

1214: .seealso: KSPSetTolerances(), KSP
1215: @*/
1216: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1217: {
1220:   if (abstol) *abstol = ksp->abstol;
1221:   if (rtol) *rtol = ksp->rtol;
1222:   if (dtol) *dtol = ksp->divtol;
1223:   if (maxits) *maxits = ksp->max_it;
1224:   return(0);
1225: }

1227: /*@
1228:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1229:    iteration tolerances used by the default KSP convergence testers.

1231:    Logically Collective on KSP

1233:    Input Parameters:
1234: +  ksp - the Krylov subspace context
1235: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1236: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1237: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1238: -  maxits - maximum number of iterations to use

1240:    Options Database Keys:
1241: +  -ksp_atol <abstol> - Sets abstol
1242: .  -ksp_rtol <rtol> - Sets rtol
1243: .  -ksp_divtol <dtol> - Sets dtol
1244: -  -ksp_max_it <maxits> - Sets maxits

1246:    Notes:
1247:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

1249:    See KSPConvergedDefault() for details how these parameters are used in the default convergence test.  See also KSPSetConvergenceTest()
1250:    for setting user-defined stopping criteria.

1252:    Level: intermediate

1254: .keywords: set, tolerance, absolute, relative, divergence,
1255:            convergence, maximum, iterations

1257: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
1258: @*/
1259: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1260: {

1268:   if (rtol != PETSC_DEFAULT) {
1269:     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol);
1270:     ksp->rtol = rtol;
1271:   }
1272:   if (abstol != PETSC_DEFAULT) {
1273:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1274:     ksp->abstol = abstol;
1275:   }
1276:   if (dtol != PETSC_DEFAULT) {
1277:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1278:     ksp->divtol = dtol;
1279:   }
1280:   if (maxits != PETSC_DEFAULT) {
1281:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1282:     ksp->max_it = maxits;
1283:   }
1284:   return(0);
1285: }

1287: /*@
1288:    KSPSetInitialGuessNonzero - Tells the iterative solver that the
1289:    initial guess is nonzero; otherwise KSP assumes the initial guess
1290:    is to be zero (and thus zeros it out before solving).

1292:    Logically Collective on KSP

1294:    Input Parameters:
1295: +  ksp - iterative context obtained from KSPCreate()
1296: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

1298:    Options database keys:
1299: .  -ksp_initial_guess_nonzero : use nonzero initial guess; this takes an optional truth value (0/1/no/yes/true/false)

1301:    Level: beginner

1303:    Notes:
1304:     If this is not called the X vector is zeroed in the call to KSPSolve().

1306: .keywords: set, initial guess, nonzero

1308: .seealso: KSPGetInitialGuessNonzero(), KSPSetGuessType(), KSPGuessType, KSP
1309: @*/
1310: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1311: {
1315:   ksp->guess_zero = (PetscBool) !(int)flg;
1316:   return(0);
1317: }

1319: /*@
1320:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1321:    a zero initial guess.

1323:    Not Collective

1325:    Input Parameter:
1326: .  ksp - iterative context obtained from KSPCreate()

1328:    Output Parameter:
1329: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1331:    Level: intermediate

1333: .keywords: set, initial guess, nonzero

1335: .seealso: KSPSetInitialGuessNonzero(), KSP
1336: @*/
1337: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1338: {
1342:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1343:   else *flag = PETSC_TRUE;
1344:   return(0);
1345: }

1347: /*@
1348:    KSPSetErrorIfNotConverged - Causes KSPSolve() to generate an error if the solver has not converged.

1350:    Logically Collective on KSP

1352:    Input Parameters:
1353: +  ksp - iterative context obtained from KSPCreate()
1354: -  flg - PETSC_TRUE indicates you want the error generated

1356:    Options database keys:
1357: .  -ksp_error_if_not_converged : this takes an optional truth value (0/1/no/yes/true/false)

1359:    Level: intermediate

1361:    Notes:
1362:     Normally PETSc continues if a linear solver fails to converge, you can call KSPGetConvergedReason() after a KSPSolve()
1363:     to determine if it has converged.


1366: .seealso: KSPGetErrorIfNotConverged(), KSP
1367: @*/
1368: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1369: {
1373:   ksp->errorifnotconverged = flg;
1374:   return(0);
1375: }

1377: /*@
1378:    KSPGetErrorIfNotConverged - Will KSPSolve() generate an error if the solver does not converge?

1380:    Not Collective

1382:    Input Parameter:
1383: .  ksp - iterative context obtained from KSPCreate()

1385:    Output Parameter:
1386: .  flag - PETSC_TRUE if it will generate an error, else PETSC_FALSE

1388:    Level: intermediate

1390: .seealso: KSPSetErrorIfNotConverged(), KSP
1391: @*/
1392: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1393: {
1397:   *flag = ksp->errorifnotconverged;
1398:   return(0);
1399: }

1401: /*@
1402:    KSPSetInitialGuessKnoll - Tells the iterative solver to use PCApply(pc,b,..) to compute the initial guess (The Knoll trick)

1404:    Logically Collective on KSP

1406:    Input Parameters:
1407: +  ksp - iterative context obtained from KSPCreate()
1408: -  flg - PETSC_TRUE or PETSC_FALSE

1410:    Level: advanced

1412:    Developer Note: the Knoll trick is not currently implemented using the KSPGuess class

1414: .keywords: set, initial guess, nonzero

1416: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1417: @*/
1418: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1419: {
1423:   ksp->guess_knoll = flg;
1424:   return(0);
1425: }

1427: /*@
1428:    KSPGetInitialGuessKnoll - Determines whether the KSP solver is using the Knoll trick (using PCApply(pc,b,...) to compute
1429:      the initial guess

1431:    Not Collective

1433:    Input Parameter:
1434: .  ksp - iterative context obtained from KSPCreate()

1436:    Output Parameter:
1437: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1439:    Level: advanced

1441: .keywords: set, initial guess, nonzero

1443: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1444: @*/
1445: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1446: {
1450:   *flag = ksp->guess_knoll;
1451:   return(0);
1452: }

1454: /*@
1455:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1456:    values will be calculated via a Lanczos or Arnoldi process as the linear
1457:    system is solved.

1459:    Not Collective

1461:    Input Parameter:
1462: .  ksp - iterative context obtained from KSPCreate()

1464:    Output Parameter:
1465: .  flg - PETSC_TRUE or PETSC_FALSE

1467:    Options Database Key:
1468: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1470:    Notes:
1471:    Currently this option is not valid for all iterative methods.

1473:    Many users may just want to use the monitoring routine
1474:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1475:    to print the singular values at each iteration of the linear solve.

1477:    Level: advanced

1479: .keywords: set, compute, singular values

1481: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1482: @*/
1483: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1484: {
1488:   *flg = ksp->calc_sings;
1489:   return(0);
1490: }

1492: /*@
1493:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1494:    values will be calculated via a Lanczos or Arnoldi process as the linear
1495:    system is solved.

1497:    Logically Collective on KSP

1499:    Input Parameters:
1500: +  ksp - iterative context obtained from KSPCreate()
1501: -  flg - PETSC_TRUE or PETSC_FALSE

1503:    Options Database Key:
1504: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1506:    Notes:
1507:    Currently this option is not valid for all iterative methods.

1509:    Many users may just want to use the monitoring routine
1510:    KSPMonitorSingularValue() (which can be set with option -ksp_monitor_singular_value)
1511:    to print the singular values at each iteration of the linear solve.

1513:    Level: advanced

1515: .keywords: set, compute, singular values

1517: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1518: @*/
1519: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1520: {
1524:   ksp->calc_sings = flg;
1525:   return(0);
1526: }

1528: /*@
1529:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1530:    values will be calculated via a Lanczos or Arnoldi process as the linear
1531:    system is solved.

1533:    Not Collective

1535:    Input Parameter:
1536: .  ksp - iterative context obtained from KSPCreate()

1538:    Output Parameter:
1539: .  flg - PETSC_TRUE or PETSC_FALSE

1541:    Notes:
1542:    Currently this option is not valid for all iterative methods.

1544:    Level: advanced

1546: .keywords: set, compute, eigenvalues

1548: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1549: @*/
1550: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1551: {
1555:   *flg = ksp->calc_sings;
1556:   return(0);
1557: }

1559: /*@
1560:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1561:    values will be calculated via a Lanczos or Arnoldi process as the linear
1562:    system is solved.

1564:    Logically Collective on KSP

1566:    Input Parameters:
1567: +  ksp - iterative context obtained from KSPCreate()
1568: -  flg - PETSC_TRUE or PETSC_FALSE

1570:    Notes:
1571:    Currently this option is not valid for all iterative methods.

1573:    Level: advanced

1575: .keywords: set, compute, eigenvalues

1577: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1578: @*/
1579: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1580: {
1584:   ksp->calc_sings = flg;
1585:   return(0);
1586: }

1588: /*@
1589:    KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1590:    will be calculated via a Lanczos or Arnoldi process as the linear
1591:    system is solved.

1593:    Logically Collective on KSP

1595:    Input Parameters:
1596: +  ksp - iterative context obtained from KSPCreate()
1597: -  flg - PETSC_TRUE or PETSC_FALSE

1599:    Notes:
1600:    Currently this option is only valid for the GMRES method.

1602:    Level: advanced

1604: .keywords: set, compute, ritz

1606: .seealso: KSPComputeRitz(), KSP
1607: @*/
1608: PetscErrorCode  KSPSetComputeRitz(KSP ksp, PetscBool flg)
1609: {
1613:   ksp->calc_ritz = flg;
1614:   return(0);
1615: }

1617: /*@
1618:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1619:    be solved.

1621:    Not Collective

1623:    Input Parameter:
1624: .  ksp - iterative context obtained from KSPCreate()

1626:    Output Parameter:
1627: .  r - right-hand-side vector

1629:    Level: developer

1631: .keywords: get, right-hand-side, rhs

1633: .seealso: KSPGetSolution(), KSPSolve(), KSP
1634: @*/
1635: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1636: {
1640:   *r = ksp->vec_rhs;
1641:   return(0);
1642: }

1644: /*@
1645:    KSPGetSolution - Gets the location of the solution for the
1646:    linear system to be solved.  Note that this may not be where the solution
1647:    is stored during the iterative process; see KSPBuildSolution().

1649:    Not Collective

1651:    Input Parameters:
1652: .  ksp - iterative context obtained from KSPCreate()

1654:    Output Parameters:
1655: .  v - solution vector

1657:    Level: developer

1659: .keywords: get, solution

1661: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve(), KSP
1662: @*/
1663: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1664: {
1668:   *v = ksp->vec_sol;
1669:   return(0);
1670: }

1672: /*@
1673:    KSPSetPC - Sets the preconditioner to be used to calculate the
1674:    application of the preconditioner on a vector.

1676:    Collective on KSP

1678:    Input Parameters:
1679: +  ksp - iterative context obtained from KSPCreate()
1680: -  pc   - the preconditioner object

1682:    Notes:
1683:    Use KSPGetPC() to retrieve the preconditioner context (for example,
1684:    to free it at the end of the computations).

1686:    Level: developer

1688: .keywords: set, precondition, Binv

1690: .seealso: KSPGetPC(), KSP
1691: @*/
1692: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1693: {

1700:   PetscObjectReference((PetscObject)pc);
1701:   PCDestroy(&ksp->pc);
1702:   ksp->pc = pc;
1703:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1704:   return(0);
1705: }

1707: /*@
1708:    KSPGetPC - Returns a pointer to the preconditioner context
1709:    set with KSPSetPC().

1711:    Not Collective

1713:    Input Parameters:
1714: .  ksp - iterative context obtained from KSPCreate()

1716:    Output Parameter:
1717: .  pc - preconditioner context

1719:    Level: developer

1721: .keywords: get, preconditioner, Binv

1723: .seealso: KSPSetPC(), KSP
1724: @*/
1725: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1726: {

1732:   if (!ksp->pc) {
1733:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1734:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1735:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1736:     PetscObjectSetOptions((PetscObject)ksp->pc,((PetscObject)ksp)->options);
1737:   }
1738:   *pc = ksp->pc;
1739:   return(0);
1740: }

1742: /*@
1743:    KSPMonitor - runs the user provided monitor routines, if they exist

1745:    Collective on KSP

1747:    Input Parameters:
1748: +  ksp - iterative context obtained from KSPCreate()
1749: .  it - iteration number
1750: -  rnorm - relative norm of the residual

1752:    Notes:
1753:    This routine is called by the KSP implementations.
1754:    It does not typically need to be called by the user.

1756:    Level: developer

1758: .seealso: KSPMonitorSet()
1759: @*/
1760: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1761: {
1762:   PetscInt       i, n = ksp->numbermonitors;

1766:   for (i=0; i<n; i++) {
1767:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1768:   }
1769:   return(0);
1770: }

1772: /*

1774:     Checks if two monitors are identical; if they are then it destroys the new one
1775: */
1776: PetscErrorCode PetscMonitorCompare(PetscErrorCode (*nmon)(void),void *nmctx,PetscErrorCode (*nmdestroy)(void**),PetscErrorCode (*mon)(void),void *mctx,PetscErrorCode (*mdestroy)(void**),PetscBool *identical)
1777: {
1778:   *identical = PETSC_FALSE;
1779:   if (nmon == mon && nmdestroy == mdestroy) {
1780:     if (nmctx == mctx) *identical = PETSC_TRUE;
1781:     else if (nmdestroy == (PetscErrorCode (*)(void**)) PetscViewerAndFormatDestroy) {
1782:       PetscViewerAndFormat *old = (PetscViewerAndFormat*)mctx, *newo = (PetscViewerAndFormat*)nmctx;
1783:       if (old->viewer == newo->viewer && old->format == newo->format) *identical = PETSC_TRUE;
1784:     }
1785:     if (*identical) {
1786:       if (mdestroy) {
1788:         (*mdestroy)(&nmctx);
1789:       }
1790:     }
1791:   }
1792:   return(0);
1793: }

1795: /*@C
1796:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1797:    the residual/error etc.

1799:    Logically Collective on KSP

1801:    Input Parameters:
1802: +  ksp - iterative context obtained from KSPCreate()
1803: .  monitor - pointer to function (if this is NULL, it turns off monitoring
1804: .  mctx    - [optional] context for private data for the
1805:              monitor routine (use NULL if no context is desired)
1806: -  monitordestroy - [optional] routine that frees monitor context
1807:           (may be NULL)

1809:    Calling Sequence of monitor:
1810: $     monitor (KSP ksp, int it, PetscReal rnorm, void *mctx)

1812: +  ksp - iterative context obtained from KSPCreate()
1813: .  it - iteration number
1814: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1815: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

1817:    Options Database Keys:
1818: +    -ksp_monitor        - sets KSPMonitorDefault()
1819: .    -ksp_monitor_true_residual    - sets KSPMonitorTrueResidualNorm()
1820: .    -ksp_monitor_max    - sets KSPMonitorTrueResidualMaxNorm()
1821: .    -ksp_monitor_lg_residualnorm    - sets line graph monitor,
1822:                            uses KSPMonitorLGResidualNormCreate()
1823: .    -ksp_monitor_lg_true_residualnorm   - sets line graph monitor,
1824:                            uses KSPMonitorLGResidualNormCreate()
1825: .    -ksp_monitor_singular_value    - sets KSPMonitorSingularValue()
1826: -    -ksp_monitor_cancel - cancels all monitors that have
1827:                           been hardwired into a code by
1828:                           calls to KSPMonitorSet(), but
1829:                           does not cancel those set via
1830:                           the options database.

1832:    Notes:
1833:    The default is to do nothing.  To print the residual, or preconditioned
1834:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1835:    KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1836:    context.

1838:    Several different monitoring routines may be set by calling
1839:    KSPMonitorSet() multiple times; all will be called in the
1840:    order in which they were set.

1842:    Fortran Notes:
1843:     Only a single monitor function can be set for each KSP object

1845:    Level: beginner

1847: .keywords: set, monitor

1849: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel(), KSP
1850: @*/
1851: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1852: {
1853:   PetscInt       i;
1855:   PetscBool      identical;

1859:   for (i=0; i<ksp->numbermonitors;i++) {
1860:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
1861:     if (identical) return(0);
1862:   }
1863:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1864:   ksp->monitor[ksp->numbermonitors]          = monitor;
1865:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
1866:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1867:   return(0);
1868: }

1870: /*@
1871:    KSPMonitorCancel - Clears all monitors for a KSP object.

1873:    Logically Collective on KSP

1875:    Input Parameters:
1876: .  ksp - iterative context obtained from KSPCreate()

1878:    Options Database Key:
1879: .  -ksp_monitor_cancel - Cancels all monitors that have
1880:     been hardwired into a code by calls to KSPMonitorSet(),
1881:     but does not cancel those set via the options database.

1883:    Level: intermediate

1885: .keywords: set, monitor

1887: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet(), KSP
1888: @*/
1889: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1890: {
1892:   PetscInt       i;

1896:   for (i=0; i<ksp->numbermonitors; i++) {
1897:     if (ksp->monitordestroy[i]) {
1898:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1899:     }
1900:   }
1901:   ksp->numbermonitors = 0;
1902:   return(0);
1903: }

1905: /*@C
1906:    KSPGetMonitorContext - Gets the monitoring context, as set by
1907:    KSPMonitorSet() for the FIRST monitor only.

1909:    Not Collective

1911:    Input Parameter:
1912: .  ksp - iterative context obtained from KSPCreate()

1914:    Output Parameter:
1915: .  ctx - monitoring context

1917:    Level: intermediate

1919: .keywords: get, monitor, context

1921: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSP
1922: @*/
1923: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1924: {
1927:   *ctx =      (ksp->monitorcontext[0]);
1928:   return(0);
1929: }

1931: /*@
1932:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1933:    If set, this array will contain the residual norms computed at each
1934:    iteration of the solver.

1936:    Not Collective

1938:    Input Parameters:
1939: +  ksp - iterative context obtained from KSPCreate()
1940: .  a   - array to hold history
1941: .  na  - size of a
1942: -  reset - PETSC_TRUE indicates the history counter is reset to zero
1943:            for each new linear solve

1945:    Level: advanced

1947:    Notes:
1948:     The array is NOT freed by PETSc so the user needs to keep track of
1949:            it and destroy once the KSP object is destroyed.

1951:    If 'a' is NULL then space is allocated for the history. If 'na' PETSC_DECIDE or PETSC_DEFAULT then a
1952:    default array of length 10000 is allocated.

1954: .keywords: set, residual, history, norm

1956: .seealso: KSPGetResidualHistory(), KSP

1958: @*/
1959: PetscErrorCode  KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
1960: {


1966:   PetscFree(ksp->res_hist_alloc);
1967:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
1968:     ksp->res_hist     = a;
1969:     ksp->res_hist_max = na;
1970:   } else {
1971:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
1972:     else                                           ksp->res_hist_max = 10000; /* like default ksp->max_it */
1973:     PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);

1975:     ksp->res_hist = ksp->res_hist_alloc;
1976:   }
1977:   ksp->res_hist_len   = 0;
1978:   ksp->res_hist_reset = reset;
1979:   return(0);
1980: }

1982: /*@C
1983:    KSPGetResidualHistory - Gets the array used to hold the residual history
1984:    and the number of residuals it contains.

1986:    Not Collective

1988:    Input Parameter:
1989: .  ksp - iterative context obtained from KSPCreate()

1991:    Output Parameters:
1992: +  a   - pointer to array to hold history (or NULL)
1993: -  na  - number of used entries in a (or NULL)

1995:    Level: advanced

1997:    Notes:
1998:      Can only be called after a KSPSetResidualHistory() otherwise a and na are set to zero

2000:      The Fortran version of this routine has a calling sequence
2001: $   call KSPGetResidualHistory(KSP ksp, integer na, integer ierr)
2002:     note that you have passed a Fortran array into KSPSetResidualHistory() and you need
2003:     to access the residual values from this Fortran array you provided. Only the na (number of
2004:     residual norms currently held) is set.

2006: .keywords: get, residual, history, norm

2008: .seealso: KSPGetResidualHistory(), KSP

2010: @*/
2011: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
2012: {
2015:   if (a) *a = ksp->res_hist;
2016:   if (na) *na = ksp->res_hist_len;
2017:   return(0);
2018: }

2020: /*@C
2021:    KSPSetConvergenceTest - Sets the function to be used to determine
2022:    convergence.

2024:    Logically Collective on KSP

2026:    Input Parameters:
2027: +  ksp - iterative context obtained from KSPCreate()
2028: .  converge - pointer to int function
2029: .  cctx    - context for private data for the convergence routine (may be null)
2030: -  destroy - a routine for destroying the context (may be null)

2032:    Calling sequence of converge:
2033: $     converge (KSP ksp, int it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)

2035: +  ksp - iterative context obtained from KSPCreate()
2036: .  it - iteration number
2037: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2038: .  reason - the reason why it has converged or diverged
2039: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


2042:    Notes:
2043:    Must be called after the KSP type has been set so put this after
2044:    a call to KSPSetType(), or KSPSetFromOptions().

2046:    The default convergence test, KSPConvergedDefault(), aborts if the
2047:    residual grows to more than 10000 times the initial residual.

2049:    The default is a combination of relative and absolute tolerances.
2050:    The residual value that is tested may be an approximation; routines
2051:    that need exact values should compute them.

2053:    In the default PETSc convergence test, the precise values of reason
2054:    are macros such as KSP_CONVERGED_RTOL, which are defined in petscksp.h.

2056:    Level: advanced

2058: .keywords: set, convergence, test, context

2060: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPGetConvergenceTest(), KSPGetAndClearConvergenceTest()
2061: @*/
2062: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2063: {

2068:   if (ksp->convergeddestroy) {
2069:     (*ksp->convergeddestroy)(ksp->cnvP);
2070:   }
2071:   ksp->converged        = converge;
2072:   ksp->convergeddestroy = destroy;
2073:   ksp->cnvP             = (void*)cctx;
2074:   return(0);
2075: }

2077: /*@C
2078:    KSPGetConvergenceTest - Gets the function to be used to determine
2079:    convergence.

2081:    Logically Collective on KSP

2083:    Input Parameter:
2084: .   ksp - iterative context obtained from KSPCreate()

2086:    Output Parameter:
2087: +  converge - pointer to convergence test function
2088: .  cctx    - context for private data for the convergence routine (may be null)
2089: -  destroy - a routine for destroying the context (may be null)

2091:    Calling sequence of converge:
2092: $     converge (KSP ksp, int it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)

2094: +  ksp - iterative context obtained from KSPCreate()
2095: .  it - iteration number
2096: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2097: .  reason - the reason why it has converged or diverged
2098: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2100:    Level: advanced

2102: .keywords: set, convergence, test, context

2104: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetAndClearConvergenceTest()
2105: @*/
2106: PetscErrorCode  KSPGetConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2107: {
2110:   if (converge) *converge = ksp->converged;
2111:   if (destroy)  *destroy  = ksp->convergeddestroy;
2112:   if (cctx)     *cctx     = ksp->cnvP;
2113:   return(0);
2114: }

2116: /*@C
2117:    KSPGetAndClearConvergenceTest - Gets the function to be used to determine convergence. Removes the current test without calling destroy on the test context

2119:    Logically Collective on KSP

2121:    Input Parameter:
2122: .   ksp - iterative context obtained from KSPCreate()

2124:    Output Parameter:
2125: +  converge - pointer to convergence test function
2126: .  cctx    - context for private data for the convergence routine
2127: -  destroy - a routine for destroying the context

2129:    Calling sequence of converge:
2130: $     converge (KSP ksp, int it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)

2132: +  ksp - iterative context obtained from KSPCreate()
2133: .  it - iteration number
2134: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2135: .  reason - the reason why it has converged or diverged
2136: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2138:    Level: advanced

2140:    Notes: This is intended to be used to allow transfering the convergence test (and its context) to another testing object (for example another KSP) and then calling
2141:           KSPSetConvergenceTest() on this original KSP. If you just called KSPGetConvergenceTest() followed by KSPSetConvergenceTest() the original context information
2142:           would be destroyed and hence the transfered context would be invalid and trigger a crash on use

2144: .keywords: set, convergence, test, context

2146: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetConvergenceTest()
2147: @*/
2148: PetscErrorCode  KSPGetAndClearConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2149: {
2152:   *converge             = ksp->converged;
2153:   *destroy              = ksp->convergeddestroy;
2154:   *cctx                 = ksp->cnvP;
2155:   ksp->converged        = NULL;
2156:   ksp->cnvP             = NULL;
2157:   ksp->convergeddestroy = NULL;
2158:   return(0);
2159: }

2161: /*@C
2162:    KSPGetConvergenceContext - Gets the convergence context set with
2163:    KSPSetConvergenceTest().

2165:    Not Collective

2167:    Input Parameter:
2168: .  ksp - iterative context obtained from KSPCreate()

2170:    Output Parameter:
2171: .  ctx - monitoring context

2173:    Level: advanced

2175: .keywords: get, convergence, test, context

2177: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
2178: @*/
2179: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
2180: {
2183:   *ctx = ksp->cnvP;
2184:   return(0);
2185: }

2187: /*@C
2188:    KSPBuildSolution - Builds the approximate solution in a vector provided.
2189:    This routine is NOT commonly needed (see KSPSolve()).

2191:    Collective on KSP

2193:    Input Parameter:
2194: .  ctx - iterative context obtained from KSPCreate()

2196:    Output Parameter:
2197:    Provide exactly one of
2198: +  v - location to stash solution.
2199: -  V - the solution is returned in this location. This vector is created
2200:        internally. This vector should NOT be destroyed by the user with
2201:        VecDestroy().

2203:    Notes:
2204:    This routine can be used in one of two ways
2205: .vb
2206:       KSPBuildSolution(ksp,NULL,&V);
2207:    or
2208:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2209: .ve
2210:    In the first case an internal vector is allocated to store the solution
2211:    (the user cannot destroy this vector). In the second case the solution
2212:    is generated in the vector that the user provides. Note that for certain
2213:    methods, such as KSPCG, the second case requires a copy of the solution,
2214:    while in the first case the call is essentially free since it simply
2215:    returns the vector where the solution already is stored. For some methods
2216:    like GMRES this is a reasonably expensive operation and should only be
2217:    used in truly needed.

2219:    Level: advanced

2221: .keywords: build, solution

2223: .seealso: KSPGetSolution(), KSPBuildResidual(), KSP
2224: @*/
2225: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2226: {

2231:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2232:   if (!V) V = &v;
2233:   (*ksp->ops->buildsolution)(ksp,v,V);
2234:   return(0);
2235: }

2237: /*@C
2238:    KSPBuildResidual - Builds the residual in a vector provided.

2240:    Collective on KSP

2242:    Input Parameter:
2243: .  ksp - iterative context obtained from KSPCreate()

2245:    Output Parameters:
2246: +  v - optional location to stash residual.  If v is not provided,
2247:        then a location is generated.
2248: .  t - work vector.  If not provided then one is generated.
2249: -  V - the residual

2251:    Notes:
2252:    Regardless of whether or not v is provided, the residual is
2253:    returned in V.

2255:    Level: advanced

2257: .keywords: KSP, build, residual

2259: .seealso: KSPBuildSolution()
2260: @*/
2261: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2262: {
2264:   PetscBool      flag = PETSC_FALSE;
2265:   Vec            w    = v,tt = t;

2269:   if (!w) {
2270:     VecDuplicate(ksp->vec_rhs,&w);
2271:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2272:   }
2273:   if (!tt) {
2274:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2275:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2276:   }
2277:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
2278:   if (flag) {VecDestroy(&tt);}
2279:   return(0);
2280: }

2282: /*@
2283:    KSPSetDiagonalScale - Tells KSP to symmetrically diagonally scale the system
2284:      before solving. This actually CHANGES the matrix (and right hand side).

2286:    Logically Collective on KSP

2288:    Input Parameter:
2289: +  ksp - the KSP context
2290: -  scale - PETSC_TRUE or PETSC_FALSE

2292:    Options Database Key:
2293: +   -ksp_diagonal_scale -
2294: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


2297:     Notes:
2298:     Scales the matrix by  D^(-1/2)  A  D^(-1/2)  [D^(1/2) x ] = D^(-1/2) b
2299:        where D_{ii} is 1/abs(A_{ii}) unless A_{ii} is zero and then it is 1.

2301:     BE CAREFUL with this routine: it actually scales the matrix and right
2302:     hand side that define the system. After the system is solved the matrix
2303:     and right hand side remain scaled unless you use KSPSetDiagonalScaleFix()

2305:     This should NOT be used within the SNES solves if you are using a line
2306:     search.

2308:     If you use this with the PCType Eisenstat preconditioner than you can
2309:     use the PCEisenstatSetNoDiagonalScaling() option, or -pc_eisenstat_no_diagonal_scaling
2310:     to save some unneeded, redundant flops.

2312:    Level: intermediate

2314: .keywords: set, options, prefix, database

2316: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2317: @*/
2318: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2319: {
2323:   ksp->dscale = scale;
2324:   return(0);
2325: }

2327: /*@
2328:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2329:                           right hand side

2331:    Not Collective

2333:    Input Parameter:
2334: .  ksp - the KSP context

2336:    Output Parameter:
2337: .  scale - PETSC_TRUE or PETSC_FALSE

2339:    Notes:
2340:     BE CAREFUL with this routine: it actually scales the matrix and right
2341:     hand side that define the system. After the system is solved the matrix
2342:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

2344:    Level: intermediate

2346: .keywords: set, options, prefix, database

2348: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2349: @*/
2350: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
2351: {
2355:   *scale = ksp->dscale;
2356:   return(0);
2357: }

2359: /*@
2360:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2361:      back after solving.

2363:    Logically Collective on KSP

2365:    Input Parameter:
2366: +  ksp - the KSP context
2367: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2368:          rescale (default)

2370:    Notes:
2371:      Must be called after KSPSetDiagonalScale()

2373:      Using this will slow things down, because it rescales the matrix before and
2374:      after each linear solve. This is intended mainly for testing to allow one
2375:      to easily get back the original system to make sure the solution computed is
2376:      accurate enough.

2378:    Level: intermediate

2380: .keywords: set, options, prefix, database

2382: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix(), KSP
2383: @*/
2384: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2385: {
2389:   ksp->dscalefix = fix;
2390:   return(0);
2391: }

2393: /*@
2394:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2395:      back after solving.

2397:    Not Collective

2399:    Input Parameter:
2400: .  ksp - the KSP context

2402:    Output Parameter:
2403: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2404:          rescale (default)

2406:    Notes:
2407:      Must be called after KSPSetDiagonalScale()

2409:      If PETSC_TRUE will slow things down, because it rescales the matrix before and
2410:      after each linear solve. This is intended mainly for testing to allow one
2411:      to easily get back the original system to make sure the solution computed is
2412:      accurate enough.

2414:    Level: intermediate

2416: .keywords: set, options, prefix, database

2418: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2419: @*/
2420: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2421: {
2425:   *fix = ksp->dscalefix;
2426:   return(0);
2427: }

2429: /*@C
2430:    KSPSetComputeOperators - set routine to compute the linear operators

2432:    Logically Collective

2434:    Input Arguments:
2435: +  ksp - the KSP context
2436: .  func - function to compute the operators
2437: -  ctx - optional context

2439:    Calling sequence of func:
2440: $  func(KSP ksp,Mat A,Mat B,void *ctx)

2442: +  ksp - the KSP context
2443: .  A - the linear operator
2444: .  B - preconditioning matrix
2445: -  ctx - optional user-provided context

2447:    Notes:
2448:     The user provided func() will be called automatically at the very next call to KSPSolve(). It will not be called at future KSPSolve() calls
2449:           unless either KSPSetComputeOperators() or KSPSetOperators() is called before that KSPSolve() is called.

2451:           To reuse the same preconditioner for the next KSPSolve() and not compute a new one based on the most recently computed matrix call KSPSetReusePreconditioner()

2453:    Level: beginner

2455: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2456: @*/
2457: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2458: {
2460:   DM             dm;

2464:   KSPGetDM(ksp,&dm);
2465:   DMKSPSetComputeOperators(dm,func,ctx);
2466:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2467:   return(0);
2468: }

2470: /*@C
2471:    KSPSetComputeRHS - set routine to compute the right hand side of the linear system

2473:    Logically Collective

2475:    Input Arguments:
2476: +  ksp - the KSP context
2477: .  func - function to compute the right hand side
2478: -  ctx - optional context

2480:    Calling sequence of func:
2481: $  func(KSP ksp,Vec b,void *ctx)

2483: +  ksp - the KSP context
2484: .  b - right hand side of linear system
2485: -  ctx - optional user-provided context

2487:    Notes:
2488:     The routine you provide will be called EACH you call KSPSolve() to prepare the new right hand side for that solve

2490:    Level: beginner

2492: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2493: @*/
2494: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2495: {
2497:   DM             dm;

2501:   KSPGetDM(ksp,&dm);
2502:   DMKSPSetComputeRHS(dm,func,ctx);
2503:   return(0);
2504: }

2506: /*@C
2507:    KSPSetComputeInitialGuess - set routine to compute the initial guess of the linear system

2509:    Logically Collective

2511:    Input Arguments:
2512: +  ksp - the KSP context
2513: .  func - function to compute the initial guess
2514: -  ctx - optional context

2516:    Calling sequence of func:
2517: $  func(KSP ksp,Vec x,void *ctx)

2519: +  ksp - the KSP context
2520: .  x - solution vector
2521: -  ctx - optional user-provided context

2523:    Notes: This should only be used in conjunction with KSPSetComputeRHS(), KSPSetComputeOperators(), otherwise
2524:    call KSPSetInitialGuessNonzero() and set the initial guess values in the solution vector passed to KSPSolve().

2526:    Level: beginner

2528: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2529: @*/
2530: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2531: {
2533:   DM             dm;

2537:   KSPGetDM(ksp,&dm);
2538:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2539:   return(0);
2540: }