Actual source code: itfunc.c

petsc-master 2020-05-26
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_view_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_view_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: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeEigenvalues(), KSP
 52: @*/
 53: PetscErrorCode  KSPComputeExtremeSingularValues(KSP ksp,PetscReal *emax,PetscReal *emin)
 54: {

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

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

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

 76:    Not Collective

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

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

 88:    Options Database Keys:
 89: .  -ksp_view_eigenvalues - Prints eigenvalues to stdout

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

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

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

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

110:    Level: advanced

112: .seealso: KSPSetComputeSingularValues(), KSPMonitorSingularValue(), KSPComputeExtremeSingularValues(), KSP
113: @*/
114: PetscErrorCode  KSPComputeEigenvalues(KSP ksp,PetscInt n,PetscReal r[],PetscReal c[],PetscInt *neig)
115: {

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

126:   if (n && ksp->ops->computeeigenvalues) {
127:     (*ksp->ops->computeeigenvalues)(ksp,n,r,c,neig);
128:   } else {
129:     *neig = 0;
130:   }
131:   return(0);
132: }

134: /*@
135:    KSPComputeRitz - Computes the Ritz or harmonic Ritz pairs associated to the
136:    smallest or largest in modulus, for the preconditioned operator.
137:    Called after KSPSolve().

139:    Not Collective

141:    Input Parameter:
142: +  ksp   - iterative context obtained from KSPCreate()
143: .  ritz  - PETSC_TRUE or PETSC_FALSE for ritz pairs or harmonic Ritz pairs, respectively
144: .  small - PETSC_TRUE or PETSC_FALSE for smallest or largest (harmonic) Ritz values, respectively
145: -  nrit  - number of (harmonic) Ritz pairs to compute

147:    Output Parameters:
148: +  nrit  - actual number of computed (harmonic) Ritz pairs
149: .  S     - multidimensional vector with Ritz vectors
150: .  tetar - real part of the Ritz values
151: -  tetai - imaginary part of the Ritz values

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

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

168:    Level: advanced

170: .seealso: KSPSetComputeRitz(), KSP
171: @*/
172: PetscErrorCode  KSPComputeRitz(KSP ksp,PetscBool ritz,PetscBool small,PetscInt *nrit,Vec S[],PetscReal tetar[],PetscReal tetai[])
173: {

178:   if (!ksp->calc_ritz) SETERRQ(PetscObjectComm((PetscObject)ksp),4,"Ritz pairs not requested before KSPSetUp()");
179:   if (ksp->ops->computeritz) {(*ksp->ops->computeritz)(ksp,ritz,small,nrit,S,tetar,tetai);}
180:   return(0);
181: }
182: /*@
183:    KSPSetUpOnBlocks - Sets up the preconditioner for each block in
184:    the block Jacobi, block Gauss-Seidel, and overlapping Schwarz
185:    methods.

187:    Collective on ksp

189:    Input Parameter:
190: .  ksp - the KSP context

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

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

201:    Level: advanced

203: .seealso: PCSetUpOnBlocks(), KSPSetUp(), PCSetUp(), KSP
204: @*/
205: PetscErrorCode  KSPSetUpOnBlocks(KSP ksp)
206: {
207:   PC             pc;
209:   PCFailedReason pcreason;

213:   KSPGetPC(ksp,&pc);
214:   PCSetUpOnBlocks(pc);
215:   PCGetFailedReason(pc,&pcreason);
216:   if (pcreason) {
217:     ksp->reason = KSP_DIVERGED_PC_FAILED;
218:   }
219:   return(0);
220: }

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

225:    Collective on ksp

227:    Input Parameters:
228: +  ksp   - iterative context obtained from KSPCreate()
229: -  flag - PETSC_TRUE to reuse the current preconditioner

231:    Level: intermediate

233: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner(), KSP
234: @*/
235: PetscErrorCode  KSPSetReusePreconditioner(KSP ksp,PetscBool flag)
236: {
237:   PC             pc;

242:   KSPGetPC(ksp,&pc);
243:   PCSetReusePreconditioner(pc,flag);
244:   return(0);
245: }

247: /*@
248:    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

250:    Collective on ksp

252:    Input Parameters:
253: +  ksp   - iterative context obtained from KSPCreate()
254: -  flag - PETSC_TRUE to skip calling the PCSetFromOptions()

256:    Level: intermediate

258: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner(), KSP
259: @*/
260: PetscErrorCode  KSPSetSkipPCSetFromOptions(KSP ksp,PetscBool flag)
261: {
264:   ksp->skippcsetfromoptions = flag;
265:   return(0);
266: }

268: /*@
269:    KSPSetUp - Sets up the internal data structures for the
270:    later use of an iterative solver.

272:    Collective on ksp

274:    Input Parameter:
275: .  ksp   - iterative context obtained from KSPCreate()

277:    Level: developer

279: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), KSP
280: @*/
281: PetscErrorCode KSPSetUp(KSP ksp)
282: {
284:   Mat            A,B;
285:   Mat            mat,pmat;
286:   MatNullSpace   nullsp;
287:   PCFailedReason pcreason;


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

295:   if (!((PetscObject)ksp)->type_name) {
296:     KSPSetType(ksp,KSPGMRES);
297:   }
298:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);

300:   if (ksp->dmActive && !ksp->setupstage) {
301:     /* first time in so build matrix and vector data structures using DM */
302:     if (!ksp->vec_rhs) {DMCreateGlobalVector(ksp->dm,&ksp->vec_rhs);}
303:     if (!ksp->vec_sol) {DMCreateGlobalVector(ksp->dm,&ksp->vec_sol);}
304:     DMCreateMatrix(ksp->dm,&A);
305:     KSPSetOperators(ksp,A,A);
306:     PetscObjectDereference((PetscObject)A);
307:   }

309:   if (ksp->dmActive) {
310:     DMKSP kdm;
311:     DMGetDMKSP(ksp->dm,&kdm);

313:     if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
314:       /* only computes initial guess the first time through */
315:       (*kdm->ops->computeinitialguess)(ksp,ksp->vec_sol,kdm->initialguessctx);
316:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
317:     }
318:     if (kdm->ops->computerhs) {
319:       (*kdm->ops->computerhs)(ksp,ksp->vec_rhs,kdm->rhsctx);
320:     }

322:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
323:       if (kdm->ops->computeoperators) {
324:         KSPGetOperators(ksp,&A,&B);
325:         (*kdm->ops->computeoperators)(ksp,A,B,kdm->operatorsctx);
326:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(ksp,PETSC_FALSE);");
327:     }
328:   }

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

333:   switch (ksp->setupstage) {
334:   case KSP_SETUP_NEW:
335:     (*ksp->ops->setup)(ksp);
336:     break;
337:   case KSP_SETUP_NEWMATRIX: {   /* This should be replaced with a more general mechanism */
338:     if (ksp->setupnewmatrix) {
339:       (*ksp->ops->setup)(ksp);
340:     }
341:   } break;
342:   default: break;
343:   }

345:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
346:   PCGetOperators(ksp->pc,&mat,&pmat);
347:   /* scale the matrix if requested */
348:   if (ksp->dscale) {
349:     PetscScalar *xx;
350:     PetscInt    i,n;
351:     PetscBool   zeroflag = PETSC_FALSE;
352:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
353:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
354:       MatCreateVecs(pmat,&ksp->diagonal,NULL);
355:     }
356:     MatGetDiagonal(pmat,ksp->diagonal);
357:     VecGetLocalSize(ksp->diagonal,&n);
358:     VecGetArray(ksp->diagonal,&xx);
359:     for (i=0; i<n; i++) {
360:       if (xx[i] != 0.0) xx[i] = 1.0/PetscSqrtReal(PetscAbsScalar(xx[i]));
361:       else {
362:         xx[i]    = 1.0;
363:         zeroflag = PETSC_TRUE;
364:       }
365:     }
366:     VecRestoreArray(ksp->diagonal,&xx);
367:     if (zeroflag) {
368:       PetscInfo(ksp,"Zero detected in diagonal of matrix, using 1 at those locations\n");
369:     }
370:     MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
371:     if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
372:     ksp->dscalefix2 = PETSC_FALSE;
373:   }
374:   PetscLogEventEnd(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
375:   PCSetErrorIfFailure(ksp->pc,ksp->errorifnotconverged);
376:   PCSetUp(ksp->pc);
377:   PCGetFailedReason(ksp->pc,&pcreason);
378:   if (pcreason) {
379:     ksp->reason = KSP_DIVERGED_PC_FAILED;
380:   }

382:   MatGetNullSpace(mat,&nullsp);
383:   if (nullsp) {
384:     PetscBool test = PETSC_FALSE;
385:     PetscOptionsGetBool(((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_test_null_space",&test,NULL);
386:     if (test) {
387:       MatNullSpaceTest(nullsp,mat,NULL);
388:     }
389:   }
390:   ksp->setupstage = KSP_SETUP_NEWRHS;
391:   return(0);
392: }

394: static PetscErrorCode KSPReasonView_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
395: {
397:   PetscBool      isAscii;

400:   if (format != PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,format);}
401:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
402:   if (isAscii) {
403:     PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
404:     if (ksp->reason > 0) {
405:       if (((PetscObject) ksp)->prefix) {
406:         PetscViewerASCIIPrintf(viewer,"Linear %s solve converged due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
407:       } else {
408:         PetscViewerASCIIPrintf(viewer,"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
409:       }
410:     } else {
411:       if (((PetscObject) ksp)->prefix) {
412:         PetscViewerASCIIPrintf(viewer,"Linear %s solve did not converge due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
413:       } else {
414:         PetscViewerASCIIPrintf(viewer,"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
415:       }
416:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
417:         PCFailedReason reason;
418:         PCGetFailedReason(ksp->pc,&reason);
419:         PetscViewerASCIIPrintf(viewer,"               PC_FAILED due to %s \n",PCFailedReasons[reason]);
420:       }
421:     }
422:     PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
423:   }
424:   if (format != PETSC_VIEWER_DEFAULT) {PetscViewerPopFormat(viewer);}
425:   return(0);
426: }

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

431:    Collective on ksp

433:    Parameter:
434: +  ksp - iterative context obtained from KSPCreate()
435: -  viewer - the viewer to display the reason


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

441:    Level: beginner

443: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
444:           KSPSolveTranspose(), KSPGetIterationNumber(), KSP
445: @*/
446: PetscErrorCode KSPReasonView(KSP ksp,PetscViewer viewer)
447: {

451:   KSPReasonView_Internal(ksp, viewer, PETSC_VIEWER_DEFAULT);
452:   return(0);
453: }

455: #if defined(PETSC_HAVE_THREADSAFETY)
456: #define KSPReasonViewFromOptions KSPReasonViewFromOptionsUnsafe
457: #else
458: #endif
459: /*@C
460:   KSPReasonViewFromOptions - Processes command line options to determine if/how a KSPReason is to be viewed.

462:   Collective on ksp

464:   Input Parameters:
465: . ksp   - the KSP object

467:   Level: intermediate

469: @*/
470: PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
471: {
472:   PetscViewer       viewer;
473:   PetscBool         flg;
474:   PetscViewerFormat format;
475:   PetscErrorCode    ierr;

478:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_converged_reason",&viewer,&format,&flg);
479:   if (flg) {
480:     KSPReasonView_Internal(ksp, viewer, format);
481:     PetscViewerDestroy(&viewer);
482:   }
483:   return(0);
484: }

486:  #include <petscdraw.h>

488: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
489: {
490:   PetscReal     *r, *c;
491:   PetscInt       n, i, neig;
492:   PetscBool      isascii, isdraw;
493:   PetscMPIInt    rank;

497:   MPI_Comm_rank(PetscObjectComm((PetscObject) ksp), &rank);
498:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
499:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW,  &isdraw);
500:   if (isExplicit) {
501:     VecGetSize(ksp->vec_sol,&n);
502:     PetscMalloc2(n, &r, n, &c);
503:     KSPComputeEigenvaluesExplicitly(ksp, n, r, c);
504:     neig = n;
505:   } else {
506:     PetscInt nits;

508:     KSPGetIterationNumber(ksp, &nits);
509:     n    = nits+2;
510:     if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n");return(0);}
511:     PetscMalloc2(n, &r, n, &c);
512:     KSPComputeEigenvalues(ksp, n, r, c, &neig);
513:   }
514:   if (isascii) {
515:     PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively");
516:     for (i = 0; i < neig; ++i) {
517:       if (c[i] >= 0.0) {PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double) r[i],  (double) c[i]);}
518:       else             {PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double) r[i], -(double) c[i]);}
519:     }
520:   } else if (isdraw && !rank) {
521:     PetscDraw   draw;
522:     PetscDrawSP drawsp;

524:     if (format == PETSC_VIEWER_DRAW_CONTOUR) {
525:       KSPPlotEigenContours_Private(ksp,neig,r,c);
526:     } else {
527:       if (!ksp->eigviewer) {PetscViewerDrawOpen(PETSC_COMM_SELF,NULL,isExplicit ? "Explicitly Computed Eigenvalues" : "Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);}
528:       PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
529:       PetscDrawSPCreate(draw,1,&drawsp);
530:       PetscDrawSPReset(drawsp);
531:       for (i = 0; i < neig; ++i) {PetscDrawSPAddPoint(drawsp,r+i,c+i);}
532:       PetscDrawSPDraw(drawsp,PETSC_TRUE);
533:       PetscDrawSPSave(drawsp);
534:       PetscDrawSPDestroy(&drawsp);
535:     }
536:   }
537:   PetscFree2(r, c);
538:   return(0);
539: }

541: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
542: {
543:   PetscReal      smax, smin;
544:   PetscInt       nits;
545:   PetscBool      isascii;

549:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
550:   KSPGetIterationNumber(ksp, &nits);
551:   if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n");return(0);}
552:   KSPComputeExtremeSingularValues(ksp, &smax, &smin);
553:   if (isascii) {PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)smax,(double)smin,(double)(smax/smin));}
554:   return(0);
555: }

557: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
558: {
559:   PetscBool      isascii;

563:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
564:   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");
565:   if (isascii) {
566:     Mat       A;
567:     Vec       t;
568:     PetscReal norm;

570:     PCGetOperators(ksp->pc, &A, NULL);
571:     VecDuplicate(ksp->vec_rhs, &t);
572:     KSP_MatMult(ksp, A, ksp->vec_sol, t);
573:     VecAYPX(t, -1.0, ksp->vec_rhs);
574:     VecNorm(t, NORM_2, &norm);
575:     VecDestroy(&t);
576:     PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double) norm);
577:   }
578:   return(0);
579: }

581: static PetscErrorCode KSPSolve_Private(KSP ksp,Vec b,Vec x)
582: {
584:   PetscBool      flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
585:   Mat            mat,pmat;
586:   MPI_Comm       comm;
587:   MatNullSpace   nullsp;
588:   Vec            btmp,vec_rhs=NULL;

591:   comm = PetscObjectComm((PetscObject)ksp);
592:   if (x && x == b) {
593:     if (!ksp->guess_zero) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
594:     VecDuplicate(b,&x);
595:     inXisinB = PETSC_TRUE;
596:   }
597:   if (b) {
598:     PetscObjectReference((PetscObject)b);
599:     VecDestroy(&ksp->vec_rhs);
600:     ksp->vec_rhs = b;
601:   }
602:   if (x) {
603:     PetscObjectReference((PetscObject)x);
604:     VecDestroy(&ksp->vec_sol);
605:     ksp->vec_sol = x;
606:   }

608:   if (ksp->viewPre) {ObjectView((PetscObject) ksp, ksp->viewerPre, ksp->formatPre);}

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

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

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

617:   if (ksp->guess) {
618:     PetscObjectState ostate,state;

620:     KSPGuessSetUp(ksp->guess);
621:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
622:     KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
623:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
624:     if (state != ostate) {
625:       ksp->guess_zero = PETSC_FALSE;
626:     } else {
627:       PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
628:       ksp->guess_zero = PETSC_TRUE;
629:     }
630:   }

632:   /* KSPSetUp() scales the matrix if needed */
633:   KSPSetUp(ksp);
634:   KSPSetUpOnBlocks(ksp);

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

638:   PCGetOperators(ksp->pc,&mat,&pmat);
639:   /* diagonal scale RHS if called for */
640:   if (ksp->dscale) {
641:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
642:     /* second time in, but matrix was scaled back to original */
643:     if (ksp->dscalefix && ksp->dscalefix2) {
644:       Mat mat,pmat;

646:       PCGetOperators(ksp->pc,&mat,&pmat);
647:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
648:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
649:     }

651:     /* scale initial guess */
652:     if (!ksp->guess_zero) {
653:       if (!ksp->truediagonal) {
654:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
655:         VecCopy(ksp->diagonal,ksp->truediagonal);
656:         VecReciprocal(ksp->truediagonal);
657:       }
658:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
659:     }
660:   }
661:   PCPreSolve(ksp->pc,ksp);

663:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
664:   if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
665:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
666:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
667:     ksp->guess_zero = PETSC_FALSE;
668:   }

670:   /* can we mark the initial guess as zero for this solve? */
671:   guess_zero = ksp->guess_zero;
672:   if (!ksp->guess_zero) {
673:     PetscReal norm;

675:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
676:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
677:   }
678:   if (ksp->transpose_solve) {
679:     MatGetNullSpace(pmat,&nullsp);
680:   } else {
681:     MatGetTransposeNullSpace(pmat,&nullsp);
682:   }
683:   if (nullsp) {
684:     VecDuplicate(ksp->vec_rhs,&btmp);
685:     VecCopy(ksp->vec_rhs,btmp);
686:     MatNullSpaceRemove(nullsp,btmp);
687:     vec_rhs      = ksp->vec_rhs;
688:     ksp->vec_rhs = btmp;
689:   }
690:   VecLockReadPush(ksp->vec_rhs);
691:   if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
692:     VecSetInf(ksp->vec_sol);
693:   }
694:   (*ksp->ops->solve)(ksp);

696:   VecLockReadPop(ksp->vec_rhs);
697:   if (nullsp) {
698:     ksp->vec_rhs = vec_rhs;
699:     VecDestroy(&btmp);
700:   }

702:   ksp->guess_zero = guess_zero;

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

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

710:   /* diagonal scale solution if called for */
711:   if (ksp->dscale) {
712:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
713:     /* unscale right hand side and matrix */
714:     if (ksp->dscalefix) {
715:       Mat mat,pmat;

717:       VecReciprocal(ksp->diagonal);
718:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
719:       PCGetOperators(ksp->pc,&mat,&pmat);
720:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
721:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
722:       VecReciprocal(ksp->diagonal);
723:       ksp->dscalefix2 = PETSC_TRUE;
724:     }
725:   }
726:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
727:   if (ksp->guess) {
728:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
729:   }
730:   if (ksp->postsolve) {
731:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
732:   }

734:   PCGetOperators(ksp->pc,&mat,&pmat);
735:   if (ksp->viewEV)       {KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV,    ksp->formatEV);}
736:   if (ksp->viewEVExp)    {KSPViewEigenvalues_Internal(ksp, PETSC_TRUE,  ksp->viewerEVExp, ksp->formatEVExp);}
737:   if (ksp->viewSV)       {KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV);}
738:   if (ksp->viewFinalRes) {KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes);}
739:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,           ksp->viewerMat,    ksp->formatMat);}
740:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,          ksp->viewerPMat,   ksp->formatPMat);}
741:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs,  ksp->viewerRhs,    ksp->formatRhs);}
742:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol,  ksp->viewerSol,    ksp->formatSol);}
743:   if (ksp->view)         {ObjectView((PetscObject) ksp,           ksp->viewer,       ksp->format);}
744:   if (ksp->viewDScale)   {ObjectView((PetscObject) ksp->diagonal, ksp->viewerDScale, ksp->formatDScale);}
745:   if (ksp->viewMatExp)   {
746:     Mat A, B;

748:     PCGetOperators(ksp->pc, &A, NULL);
749:     if (ksp->transpose_solve) {
750:       Mat AT;

752:       MatCreateTranspose(A, &AT);
753:       MatComputeOperator(AT, MATAIJ, &B);
754:       MatDestroy(&AT);
755:     } else {
756:       MatComputeOperator(A, MATAIJ, &B);
757:     }
758:     ObjectView((PetscObject) B, ksp->viewerMatExp, ksp->formatMatExp);
759:     MatDestroy(&B);
760:   }
761:   if (ksp->viewPOpExp)   {
762:     Mat B;

764:     KSPComputeOperator(ksp, MATAIJ, &B);
765:     ObjectView((PetscObject) B, ksp->viewerPOpExp, ksp->formatPOpExp);
766:     MatDestroy(&B);
767:   }

769:   if (inXisinB) {
770:     VecCopy(x,b);
771:     VecDestroy(&x);
772:   }
773:   PetscObjectSAWsBlock((PetscObject)ksp);
774:   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]);
775:   return(0);
776: }

778: /*@
779:    KSPSolve - Solves linear system.

781:    Collective on ksp

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

788:    Options Database Keys:
789: +  -ksp_view_eigenvalues - compute preconditioned operators eigenvalues
790: .  -ksp_view_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and using LAPACK
791: .  -ksp_view_mat binary - save matrix to the default binary viewer
792: .  -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
793: .  -ksp_view_rhs binary - save right hand side vector to the default binary viewer
794: .  -ksp_view_solution binary - save computed solution vector to the default binary viewer
795:            (can be read later with src/ksp/tutorials/ex10.c for testing solvers)
796: .  -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
797: .  -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
798: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
799: .  -ksp_view_final_residual - print 2-norm of true linear system residual at the end of the solution process
800: -  -ksp_view - print the ksp data structure at the end of the system solution

802:    Notes:

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

806:    The operator is specified with KSPSetOperators().

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

811:    If you provide a matrix that has a MatSetNullSpace() and MatSetTransposeNullSpace() this will use that information to solve singular systems
812:    in the least squares sense with a norm minimizing solution.
813: $
814: $                   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()
815: $
816: $    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
817: $    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
818: $    direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
819: $
820: $    We recommend always using GMRES for such singular systems.
821: $    If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
822: $    If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).

824:    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
825:        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
826:        such as diagonal scaling we cannot apply the inverse of the preconditioner to a vector and thus cannot compute the nullspace(AB).


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

834:    Understanding Convergence:
835:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
836:    KSPComputeEigenvaluesExplicitly() provide information on additional
837:    options to monitor convergence and print eigenvalue information.

839:    Level: beginner

841: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
842:           KSPSolveTranspose(), KSPGetIterationNumber(), MatNullSpaceCreate(), MatSetNullSpace(), MatSetTransposeNullSpace(), KSP
843: @*/
844: PetscErrorCode KSPSolve(KSP ksp,Vec b,Vec x)
845: {

852:   ksp->transpose_solve = PETSC_FALSE;
853:   KSPSolve_Private(ksp,b,x);
854:   return(0);
855: }

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

860:    Collective on ksp

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

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

870:    Developer Notes:
871:     We need to implement a KSPSolveHermitianTranspose()

873:    Level: developer

875: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
876:           KSPSolve(), KSP
877: @*/
878: PetscErrorCode KSPSolveTranspose(KSP ksp,Vec b,Vec x)
879: {

886:   ksp->transpose_solve = PETSC_TRUE;
887:   KSPSolve_Private(ksp,b,x);
888:   return(0);
889: }

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

894:    Collective on ksp

896:    Input Parameter:
897: .  ksp - iterative context obtained from KSPCreate()

899:    Level: beginner

901: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSPSetFromOptions(), KSP
902: @*/
903: PetscErrorCode  KSPResetViewers(KSP ksp)
904: {

909:   if (!ksp) return(0);
910:   PetscViewerDestroy(&ksp->viewer);
911:   PetscViewerDestroy(&ksp->viewerPre);
912:   PetscViewerDestroy(&ksp->viewerReason);
913:   PetscViewerDestroy(&ksp->viewerMat);
914:   PetscViewerDestroy(&ksp->viewerPMat);
915:   PetscViewerDestroy(&ksp->viewerRhs);
916:   PetscViewerDestroy(&ksp->viewerSol);
917:   PetscViewerDestroy(&ksp->viewerMatExp);
918:   PetscViewerDestroy(&ksp->viewerEV);
919:   PetscViewerDestroy(&ksp->viewerSV);
920:   PetscViewerDestroy(&ksp->viewerEVExp);
921:   PetscViewerDestroy(&ksp->viewerFinalRes);
922:   PetscViewerDestroy(&ksp->viewerPOpExp);
923:   PetscViewerDestroy(&ksp->viewerDScale);
924:   ksp->view         = PETSC_FALSE;
925:   ksp->viewPre      = PETSC_FALSE;
926:   ksp->viewReason   = PETSC_FALSE;
927:   ksp->viewMat      = PETSC_FALSE;
928:   ksp->viewPMat     = PETSC_FALSE;
929:   ksp->viewRhs      = PETSC_FALSE;
930:   ksp->viewSol      = PETSC_FALSE;
931:   ksp->viewMatExp   = PETSC_FALSE;
932:   ksp->viewEV       = PETSC_FALSE;
933:   ksp->viewSV       = PETSC_FALSE;
934:   ksp->viewEVExp    = PETSC_FALSE;
935:   ksp->viewFinalRes = PETSC_FALSE;
936:   ksp->viewPOpExp   = PETSC_FALSE;
937:   ksp->viewDScale   = PETSC_FALSE;
938:   return(0);
939: }

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

944:    Collective on ksp

946:    Input Parameter:
947: .  ksp - iterative context obtained from KSPCreate()

949:    Level: beginner

951: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
952: @*/
953: PetscErrorCode  KSPReset(KSP ksp)
954: {

959:   if (!ksp) return(0);
960:   if (ksp->ops->reset) {
961:     (*ksp->ops->reset)(ksp);
962:   }
963:   if (ksp->pc) {PCReset(ksp->pc);}
964:   if (ksp->guess) {
965:     KSPGuess guess = ksp->guess;
966:     if (guess->ops->reset) { (*guess->ops->reset)(guess); }
967:   }
968:   VecDestroyVecs(ksp->nwork,&ksp->work);
969:   VecDestroy(&ksp->vec_rhs);
970:   VecDestroy(&ksp->vec_sol);
971:   VecDestroy(&ksp->diagonal);
972:   VecDestroy(&ksp->truediagonal);

974:   KSPResetViewers(ksp);

976:   ksp->setupstage = KSP_SETUP_NEW;
977:   return(0);
978: }

980: /*@
981:    KSPDestroy - Destroys KSP context.

983:    Collective on ksp

985:    Input Parameter:
986: .  ksp - iterative context obtained from KSPCreate()

988:    Level: beginner

990: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
991: @*/
992: PetscErrorCode  KSPDestroy(KSP *ksp)
993: {
995:   PC             pc;

998:   if (!*ksp) return(0);
1000:   if (--((PetscObject)(*ksp))->refct > 0) {*ksp = NULL; return(0);}

1002:   PetscObjectSAWsViewOff((PetscObject)*ksp);

1004:   /*
1005:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1006:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1007:    refcount (and may be shared, e.g., by other ksps).
1008:    */
1009:   pc         = (*ksp)->pc;
1010:   (*ksp)->pc = NULL;
1011:   KSPReset((*ksp));
1012:   (*ksp)->pc = pc;
1013:   if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

1015:   KSPGuessDestroy(&(*ksp)->guess);
1016:   DMDestroy(&(*ksp)->dm);
1017:   PCDestroy(&(*ksp)->pc);
1018:   PetscFree((*ksp)->res_hist_alloc);
1019:   if ((*ksp)->convergeddestroy) {
1020:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
1021:   }
1022:   KSPMonitorCancel((*ksp));
1023:   PetscViewerDestroy(&(*ksp)->eigviewer);
1024:   PetscHeaderDestroy(ksp);
1025:   return(0);
1026: }

1028: /*@
1029:     KSPSetPCSide - Sets the preconditioning side.

1031:     Logically Collective on ksp

1033:     Input Parameter:
1034: .   ksp - iterative context obtained from KSPCreate()

1036:     Output Parameter:
1037: .   side - the preconditioning side, where side is one of
1038: .vb
1039:       PC_LEFT - left preconditioning (default)
1040:       PC_RIGHT - right preconditioning
1041:       PC_SYMMETRIC - symmetric preconditioning
1042: .ve

1044:     Options Database Keys:
1045: .   -ksp_pc_side <right,left,symmetric>

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

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

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

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

1058:     Level: intermediate

1060: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType(), KSP
1061: @*/
1062: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
1063: {
1067:   ksp->pc_side = ksp->pc_side_set = side;
1068:   return(0);
1069: }

1071: /*@
1072:     KSPGetPCSide - Gets the preconditioning side.

1074:     Not Collective

1076:     Input Parameter:
1077: .   ksp - iterative context obtained from KSPCreate()

1079:     Output Parameter:
1080: .   side - the preconditioning side, where side is one of
1081: .vb
1082:       PC_LEFT - left preconditioning (default)
1083:       PC_RIGHT - right preconditioning
1084:       PC_SYMMETRIC - symmetric preconditioning
1085: .ve

1087:     Level: intermediate

1089: .seealso: KSPSetPCSide(), KSP
1090: @*/
1091: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
1092: {

1098:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
1099:   *side = ksp->pc_side;
1100:   return(0);
1101: }

1103: /*@
1104:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1105:    iteration tolerances used by the default KSP convergence tests.

1107:    Not Collective

1109:    Input Parameter:
1110: .  ksp - the Krylov subspace context

1112:    Output Parameters:
1113: +  rtol - the relative convergence tolerance
1114: .  abstol - the absolute convergence tolerance
1115: .  dtol - the divergence tolerance
1116: -  maxits - maximum number of iterations

1118:    Notes:
1119:    The user can specify NULL for any parameter that is not needed.

1121:    Level: intermediate

1123:            maximum, iterations

1125: .seealso: KSPSetTolerances(), KSP
1126: @*/
1127: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1128: {
1131:   if (abstol) *abstol = ksp->abstol;
1132:   if (rtol) *rtol = ksp->rtol;
1133:   if (dtol) *dtol = ksp->divtol;
1134:   if (maxits) *maxits = ksp->max_it;
1135:   return(0);
1136: }

1138: /*@
1139:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1140:    iteration tolerances used by the default KSP convergence testers.

1142:    Logically Collective on ksp

1144:    Input Parameters:
1145: +  ksp - the Krylov subspace context
1146: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1147: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1148: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1149: -  maxits - maximum number of iterations to use

1151:    Options Database Keys:
1152: +  -ksp_atol <abstol> - Sets abstol
1153: .  -ksp_rtol <rtol> - Sets rtol
1154: .  -ksp_divtol <dtol> - Sets dtol
1155: -  -ksp_max_it <maxits> - Sets maxits

1157:    Notes:
1158:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

1163:    Level: intermediate

1165:            convergence, maximum, iterations

1167: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
1168: @*/
1169: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1170: {

1178:   if (rtol != PETSC_DEFAULT) {
1179:     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);
1180:     ksp->rtol = rtol;
1181:   }
1182:   if (abstol != PETSC_DEFAULT) {
1183:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1184:     ksp->abstol = abstol;
1185:   }
1186:   if (dtol != PETSC_DEFAULT) {
1187:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1188:     ksp->divtol = dtol;
1189:   }
1190:   if (maxits != PETSC_DEFAULT) {
1191:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1192:     ksp->max_it = maxits;
1193:   }
1194:   return(0);
1195: }

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

1202:    Logically Collective on ksp

1204:    Input Parameters:
1205: +  ksp - iterative context obtained from KSPCreate()
1206: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

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

1211:    Level: beginner

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

1216: .seealso: KSPGetInitialGuessNonzero(), KSPSetGuessType(), KSPGuessType, KSP
1217: @*/
1218: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1219: {
1223:   ksp->guess_zero = (PetscBool) !(int)flg;
1224:   return(0);
1225: }

1227: /*@
1228:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1229:    a zero initial guess.

1231:    Not Collective

1233:    Input Parameter:
1234: .  ksp - iterative context obtained from KSPCreate()

1236:    Output Parameter:
1237: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1239:    Level: intermediate

1241: .seealso: KSPSetInitialGuessNonzero(), KSP
1242: @*/
1243: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1244: {
1248:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1249:   else *flag = PETSC_TRUE;
1250:   return(0);
1251: }

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

1256:    Logically Collective on ksp

1258:    Input Parameters:
1259: +  ksp - iterative context obtained from KSPCreate()
1260: -  flg - PETSC_TRUE indicates you want the error generated

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

1265:    Level: intermediate

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


1272: .seealso: KSPGetErrorIfNotConverged(), KSP
1273: @*/
1274: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1275: {
1279:   ksp->errorifnotconverged = flg;
1280:   return(0);
1281: }

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

1286:    Not Collective

1288:    Input Parameter:
1289: .  ksp - iterative context obtained from KSPCreate()

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

1294:    Level: intermediate

1296: .seealso: KSPSetErrorIfNotConverged(), KSP
1297: @*/
1298: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1299: {
1303:   *flag = ksp->errorifnotconverged;
1304:   return(0);
1305: }

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

1310:    Logically Collective on ksp

1312:    Input Parameters:
1313: +  ksp - iterative context obtained from KSPCreate()
1314: -  flg - PETSC_TRUE or PETSC_FALSE

1316:    Level: advanced

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

1320: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1321: @*/
1322: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1323: {
1327:   ksp->guess_knoll = flg;
1328:   return(0);
1329: }

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

1335:    Not Collective

1337:    Input Parameter:
1338: .  ksp - iterative context obtained from KSPCreate()

1340:    Output Parameter:
1341: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1343:    Level: advanced

1345: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1346: @*/
1347: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1348: {
1352:   *flag = ksp->guess_knoll;
1353:   return(0);
1354: }

1356: /*@
1357:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1358:    values will be calculated via a Lanczos or Arnoldi process as the linear
1359:    system is solved.

1361:    Not Collective

1363:    Input Parameter:
1364: .  ksp - iterative context obtained from KSPCreate()

1366:    Output Parameter:
1367: .  flg - PETSC_TRUE or PETSC_FALSE

1369:    Options Database Key:
1370: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1372:    Notes:
1373:    Currently this option is not valid for all iterative methods.

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

1379:    Level: advanced

1381: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1382: @*/
1383: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1384: {
1388:   *flg = ksp->calc_sings;
1389:   return(0);
1390: }

1392: /*@
1393:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1394:    values will be calculated via a Lanczos or Arnoldi process as the linear
1395:    system is solved.

1397:    Logically Collective on ksp

1399:    Input Parameters:
1400: +  ksp - iterative context obtained from KSPCreate()
1401: -  flg - PETSC_TRUE or PETSC_FALSE

1403:    Options Database Key:
1404: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1406:    Notes:
1407:    Currently this option is not valid for all iterative methods.

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

1413:    Level: advanced

1415: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1416: @*/
1417: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1418: {
1422:   ksp->calc_sings = flg;
1423:   return(0);
1424: }

1426: /*@
1427:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1428:    values will be calculated via a Lanczos or Arnoldi process as the linear
1429:    system is solved.

1431:    Not Collective

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

1436:    Output Parameter:
1437: .  flg - PETSC_TRUE or PETSC_FALSE

1439:    Notes:
1440:    Currently this option is not valid for all iterative methods.

1442:    Level: advanced

1444: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1445: @*/
1446: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1447: {
1451:   *flg = ksp->calc_sings;
1452:   return(0);
1453: }

1455: /*@
1456:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1457:    values will be calculated via a Lanczos or Arnoldi process as the linear
1458:    system is solved.

1460:    Logically Collective on ksp

1462:    Input Parameters:
1463: +  ksp - iterative context obtained from KSPCreate()
1464: -  flg - PETSC_TRUE or PETSC_FALSE

1466:    Notes:
1467:    Currently this option is not valid for all iterative methods.

1469:    Level: advanced

1471: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1472: @*/
1473: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1474: {
1478:   ksp->calc_sings = flg;
1479:   return(0);
1480: }

1482: /*@
1483:    KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1484:    will be calculated via a Lanczos or Arnoldi process as the linear
1485:    system is solved.

1487:    Logically Collective on ksp

1489:    Input Parameters:
1490: +  ksp - iterative context obtained from KSPCreate()
1491: -  flg - PETSC_TRUE or PETSC_FALSE

1493:    Notes:
1494:    Currently this option is only valid for the GMRES method.

1496:    Level: advanced

1498: .seealso: KSPComputeRitz(), KSP
1499: @*/
1500: PetscErrorCode  KSPSetComputeRitz(KSP ksp, PetscBool flg)
1501: {
1505:   ksp->calc_ritz = flg;
1506:   return(0);
1507: }

1509: /*@
1510:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1511:    be solved.

1513:    Not Collective

1515:    Input Parameter:
1516: .  ksp - iterative context obtained from KSPCreate()

1518:    Output Parameter:
1519: .  r - right-hand-side vector

1521:    Level: developer

1523: .seealso: KSPGetSolution(), KSPSolve(), KSP
1524: @*/
1525: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1526: {
1530:   *r = ksp->vec_rhs;
1531:   return(0);
1532: }

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

1539:    Not Collective

1541:    Input Parameters:
1542: .  ksp - iterative context obtained from KSPCreate()

1544:    Output Parameters:
1545: .  v - solution vector

1547:    Level: developer

1549: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve(), KSP
1550: @*/
1551: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1552: {
1556:   *v = ksp->vec_sol;
1557:   return(0);
1558: }

1560: /*@
1561:    KSPSetPC - Sets the preconditioner to be used to calculate the
1562:    application of the preconditioner on a vector.

1564:    Collective on ksp

1566:    Input Parameters:
1567: +  ksp - iterative context obtained from KSPCreate()
1568: -  pc   - the preconditioner object

1570:    Notes:
1571:    Use KSPGetPC() to retrieve the preconditioner context (for example,
1572:    to free it at the end of the computations).

1574:    Level: developer

1576: .seealso: KSPGetPC(), KSP
1577: @*/
1578: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1579: {

1586:   PetscObjectReference((PetscObject)pc);
1587:   PCDestroy(&ksp->pc);
1588:   ksp->pc = pc;
1589:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1590:   return(0);
1591: }

1593: /*@
1594:    KSPGetPC - Returns a pointer to the preconditioner context
1595:    set with KSPSetPC().

1597:    Not Collective

1599:    Input Parameters:
1600: .  ksp - iterative context obtained from KSPCreate()

1602:    Output Parameter:
1603: .  pc - preconditioner context

1605:    Level: developer

1607: .seealso: KSPSetPC(), KSP
1608: @*/
1609: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1610: {

1616:   if (!ksp->pc) {
1617:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1618:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1619:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1620:     PetscObjectSetOptions((PetscObject)ksp->pc,((PetscObject)ksp)->options);
1621:   }
1622:   *pc = ksp->pc;
1623:   return(0);
1624: }

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

1629:    Collective on ksp

1631:    Input Parameters:
1632: +  ksp - iterative context obtained from KSPCreate()
1633: .  it - iteration number
1634: -  rnorm - relative norm of the residual

1636:    Notes:
1637:    This routine is called by the KSP implementations.
1638:    It does not typically need to be called by the user.

1640:    Level: developer

1642: .seealso: KSPMonitorSet()
1643: @*/
1644: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1645: {
1646:   PetscInt       i, n = ksp->numbermonitors;

1650:   for (i=0; i<n; i++) {
1651:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1652:   }
1653:   return(0);
1654: }

1656: /*@C
1657:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1658:    the residual/error etc.

1660:    Logically Collective on ksp

1662:    Input Parameters:
1663: +  ksp - iterative context obtained from KSPCreate()
1664: .  monitor - pointer to function (if this is NULL, it turns off monitoring
1665: .  mctx    - [optional] context for private data for the
1666:              monitor routine (use NULL if no context is desired)
1667: -  monitordestroy - [optional] routine that frees monitor context
1668:           (may be NULL)

1670:    Calling Sequence of monitor:
1671: $     monitor (KSP ksp, PetscInt it, PetscReal rnorm, void *mctx)

1673: +  ksp - iterative context obtained from KSPCreate()
1674: .  it - iteration number
1675: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1676: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

1678:    Options Database Keys:
1679: +    -ksp_monitor        - sets KSPMonitorDefault()
1680: .    -ksp_monitor_true_residual    - sets KSPMonitorTrueResidualNorm()
1681: .    -ksp_monitor_max    - sets KSPMonitorTrueResidualMaxNorm()
1682: .    -ksp_monitor_lg_residualnorm    - sets line graph monitor,
1683:                            uses KSPMonitorLGResidualNormCreate()
1684: .    -ksp_monitor_lg_true_residualnorm   - sets line graph monitor,
1685:                            uses KSPMonitorLGResidualNormCreate()
1686: .    -ksp_monitor_singular_value    - sets KSPMonitorSingularValue()
1687: -    -ksp_monitor_cancel - cancels all monitors that have
1688:                           been hardwired into a code by
1689:                           calls to KSPMonitorSet(), but
1690:                           does not cancel those set via
1691:                           the options database.

1693:    Notes:
1694:    The default is to do nothing.  To print the residual, or preconditioned
1695:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1696:    KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1697:    context.

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

1703:    Fortran Notes:
1704:     Only a single monitor function can be set for each KSP object

1706:    Level: beginner

1708: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel(), KSP
1709: @*/
1710: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1711: {
1712:   PetscInt       i;
1714:   PetscBool      identical;

1718:   for (i=0; i<ksp->numbermonitors;i++) {
1719:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
1720:     if (identical) return(0);
1721:   }
1722:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1723:   ksp->monitor[ksp->numbermonitors]          = monitor;
1724:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
1725:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1726:   return(0);
1727: }

1729: /*@
1730:    KSPMonitorCancel - Clears all monitors for a KSP object.

1732:    Logically Collective on ksp

1734:    Input Parameters:
1735: .  ksp - iterative context obtained from KSPCreate()

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

1742:    Level: intermediate

1744: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet(), KSP
1745: @*/
1746: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1747: {
1749:   PetscInt       i;

1753:   for (i=0; i<ksp->numbermonitors; i++) {
1754:     if (ksp->monitordestroy[i]) {
1755:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1756:     }
1757:   }
1758:   ksp->numbermonitors = 0;
1759:   return(0);
1760: }

1762: /*@C
1763:    KSPGetMonitorContext - Gets the monitoring context, as set by
1764:    KSPMonitorSet() for the FIRST monitor only.

1766:    Not Collective

1768:    Input Parameter:
1769: .  ksp - iterative context obtained from KSPCreate()

1771:    Output Parameter:
1772: .  ctx - monitoring context

1774:    Level: intermediate

1776: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSP
1777: @*/
1778: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1779: {
1782:   *ctx =      (ksp->monitorcontext[0]);
1783:   return(0);
1784: }

1786: /*@
1787:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1788:    If set, this array will contain the residual norms computed at each
1789:    iteration of the solver.

1791:    Not Collective

1793:    Input Parameters:
1794: +  ksp - iterative context obtained from KSPCreate()
1795: .  a   - array to hold history
1796: .  na  - size of a
1797: -  reset - PETSC_TRUE indicates the history counter is reset to zero
1798:            for each new linear solve

1800:    Level: advanced

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

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

1809: .seealso: KSPGetResidualHistory(), KSP

1811: @*/
1812: PetscErrorCode  KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
1813: {


1819:   PetscFree(ksp->res_hist_alloc);
1820:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
1821:     ksp->res_hist     = a;
1822:     ksp->res_hist_max = na;
1823:   } else {
1824:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
1825:     else                                           ksp->res_hist_max = 10000; /* like default ksp->max_it */
1826:     PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);

1828:     ksp->res_hist = ksp->res_hist_alloc;
1829:   }
1830:   ksp->res_hist_len   = 0;
1831:   ksp->res_hist_reset = reset;
1832:   return(0);
1833: }

1835: /*@C
1836:    KSPGetResidualHistory - Gets the array used to hold the residual history
1837:    and the number of residuals it contains.

1839:    Not Collective

1841:    Input Parameter:
1842: .  ksp - iterative context obtained from KSPCreate()

1844:    Output Parameters:
1845: +  a   - pointer to array to hold history (or NULL)
1846: -  na  - number of used entries in a (or NULL)

1848:    Level: advanced

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

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

1859: .seealso: KSPGetResidualHistory(), KSP

1861: @*/
1862: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1863: {
1866:   if (a) *a = ksp->res_hist;
1867:   if (na) *na = ksp->res_hist_len;
1868:   return(0);
1869: }

1871: /*@C
1872:    KSPSetConvergenceTest - Sets the function to be used to determine
1873:    convergence.

1875:    Logically Collective on ksp

1877:    Input Parameters:
1878: +  ksp - iterative context obtained from KSPCreate()
1879: .  converge - pointer to the function
1880: .  cctx    - context for private data for the convergence routine (may be null)
1881: -  destroy - a routine for destroying the context (may be null)

1883:    Calling sequence of converge:
1884: $     converge (KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)

1886: +  ksp - iterative context obtained from KSPCreate()
1887: .  it - iteration number
1888: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1889: .  reason - the reason why it has converged or diverged
1890: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


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

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

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

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

1907:    Level: advanced

1909: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPGetConvergenceTest(), KSPGetAndClearConvergenceTest()
1910: @*/
1911: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1912: {

1917:   if (ksp->convergeddestroy) {
1918:     (*ksp->convergeddestroy)(ksp->cnvP);
1919:   }
1920:   ksp->converged        = converge;
1921:   ksp->convergeddestroy = destroy;
1922:   ksp->cnvP             = (void*)cctx;
1923:   return(0);
1924: }

1926: /*@C
1927:    KSPGetConvergenceTest - Gets the function to be used to determine
1928:    convergence.

1930:    Logically Collective on ksp

1932:    Input Parameter:
1933: .   ksp - iterative context obtained from KSPCreate()

1935:    Output Parameter:
1936: +  converge - pointer to convergence test function
1937: .  cctx    - context for private data for the convergence routine (may be null)
1938: -  destroy - a routine for destroying the context (may be null)

1940:    Calling sequence of converge:
1941: $     converge (KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)

1943: +  ksp - iterative context obtained from KSPCreate()
1944: .  it - iteration number
1945: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1946: .  reason - the reason why it has converged or diverged
1947: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

1949:    Level: advanced

1951: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetAndClearConvergenceTest()
1952: @*/
1953: PetscErrorCode  KSPGetConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
1954: {
1957:   if (converge) *converge = ksp->converged;
1958:   if (destroy)  *destroy  = ksp->convergeddestroy;
1959:   if (cctx)     *cctx     = ksp->cnvP;
1960:   return(0);
1961: }

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

1966:    Logically Collective on ksp

1968:    Input Parameter:
1969: .   ksp - iterative context obtained from KSPCreate()

1971:    Output Parameter:
1972: +  converge - pointer to convergence test function
1973: .  cctx    - context for private data for the convergence routine
1974: -  destroy - a routine for destroying the context

1976:    Calling sequence of converge:
1977: $     converge (KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason,void *mctx)

1979: +  ksp - iterative context obtained from KSPCreate()
1980: .  it - iteration number
1981: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1982: .  reason - the reason why it has converged or diverged
1983: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

1985:    Level: advanced

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

1991: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetConvergenceTest()
1992: @*/
1993: PetscErrorCode  KSPGetAndClearConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
1994: {
1997:   *converge             = ksp->converged;
1998:   *destroy              = ksp->convergeddestroy;
1999:   *cctx                 = ksp->cnvP;
2000:   ksp->converged        = NULL;
2001:   ksp->cnvP             = NULL;
2002:   ksp->convergeddestroy = NULL;
2003:   return(0);
2004: }

2006: /*@C
2007:    KSPGetConvergenceContext - Gets the convergence context set with
2008:    KSPSetConvergenceTest().

2010:    Not Collective

2012:    Input Parameter:
2013: .  ksp - iterative context obtained from KSPCreate()

2015:    Output Parameter:
2016: .  ctx - monitoring context

2018:    Level: advanced

2020: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
2021: @*/
2022: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
2023: {
2026:   *ctx = ksp->cnvP;
2027:   return(0);
2028: }

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

2034:    Collective on ksp

2036:    Input Parameter:
2037: .  ctx - iterative context obtained from KSPCreate()

2039:    Output Parameter:
2040:    Provide exactly one of
2041: +  v - location to stash solution.
2042: -  V - the solution is returned in this location. This vector is created
2043:        internally. This vector should NOT be destroyed by the user with
2044:        VecDestroy().

2046:    Notes:
2047:    This routine can be used in one of two ways
2048: .vb
2049:       KSPBuildSolution(ksp,NULL,&V);
2050:    or
2051:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2052: .ve
2053:    In the first case an internal vector is allocated to store the solution
2054:    (the user cannot destroy this vector). In the second case the solution
2055:    is generated in the vector that the user provides. Note that for certain
2056:    methods, such as KSPCG, the second case requires a copy of the solution,
2057:    while in the first case the call is essentially free since it simply
2058:    returns the vector where the solution already is stored. For some methods
2059:    like GMRES this is a reasonably expensive operation and should only be
2060:    used in truly needed.

2062:    Level: advanced

2064: .seealso: KSPGetSolution(), KSPBuildResidual(), KSP
2065: @*/
2066: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2067: {

2072:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2073:   if (!V) V = &v;
2074:   (*ksp->ops->buildsolution)(ksp,v,V);
2075:   return(0);
2076: }

2078: /*@C
2079:    KSPBuildResidual - Builds the residual in a vector provided.

2081:    Collective on ksp

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

2086:    Output Parameters:
2087: +  v - optional location to stash residual.  If v is not provided,
2088:        then a location is generated.
2089: .  t - work vector.  If not provided then one is generated.
2090: -  V - the residual

2092:    Notes:
2093:    Regardless of whether or not v is provided, the residual is
2094:    returned in V.

2096:    Level: advanced

2098: .seealso: KSPBuildSolution()
2099: @*/
2100: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2101: {
2103:   PetscBool      flag = PETSC_FALSE;
2104:   Vec            w    = v,tt = t;

2108:   if (!w) {
2109:     VecDuplicate(ksp->vec_rhs,&w);
2110:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2111:   }
2112:   if (!tt) {
2113:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2114:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2115:   }
2116:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
2117:   if (flag) {VecDestroy(&tt);}
2118:   return(0);
2119: }

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

2125:    Logically Collective on ksp

2127:    Input Parameter:
2128: +  ksp - the KSP context
2129: -  scale - PETSC_TRUE or PETSC_FALSE

2131:    Options Database Key:
2132: +   -ksp_diagonal_scale -
2133: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


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

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

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

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

2151:    Level: intermediate

2153: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2154: @*/
2155: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2156: {
2160:   ksp->dscale = scale;
2161:   return(0);
2162: }

2164: /*@
2165:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2166:                           right hand side

2168:    Not Collective

2170:    Input Parameter:
2171: .  ksp - the KSP context

2173:    Output Parameter:
2174: .  scale - PETSC_TRUE or PETSC_FALSE

2176:    Notes:
2177:     BE CAREFUL with this routine: it actually scales the matrix and right
2178:     hand side that define the system. After the system is solved the matrix
2179:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

2181:    Level: intermediate

2183: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2184: @*/
2185: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
2186: {
2190:   *scale = ksp->dscale;
2191:   return(0);
2192: }

2194: /*@
2195:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2196:      back after solving.

2198:    Logically Collective on ksp

2200:    Input Parameter:
2201: +  ksp - the KSP context
2202: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2203:          rescale (default)

2205:    Notes:
2206:      Must be called after KSPSetDiagonalScale()

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

2213:    Level: intermediate

2215: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix(), KSP
2216: @*/
2217: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2218: {
2222:   ksp->dscalefix = fix;
2223:   return(0);
2224: }

2226: /*@
2227:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2228:      back after solving.

2230:    Not Collective

2232:    Input Parameter:
2233: .  ksp - the KSP context

2235:    Output Parameter:
2236: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2237:          rescale (default)

2239:    Notes:
2240:      Must be called after KSPSetDiagonalScale()

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

2247:    Level: intermediate

2249: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2250: @*/
2251: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2252: {
2256:   *fix = ksp->dscalefix;
2257:   return(0);
2258: }

2260: /*@C
2261:    KSPSetComputeOperators - set routine to compute the linear operators

2263:    Logically Collective

2265:    Input Arguments:
2266: +  ksp - the KSP context
2267: .  func - function to compute the operators
2268: -  ctx - optional context

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

2273: +  ksp - the KSP context
2274: .  A - the linear operator
2275: .  B - preconditioning matrix
2276: -  ctx - optional user-provided context

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

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

2284:    Level: beginner

2286: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2287: @*/
2288: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2289: {
2291:   DM             dm;

2295:   KSPGetDM(ksp,&dm);
2296:   DMKSPSetComputeOperators(dm,func,ctx);
2297:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2298:   return(0);
2299: }

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

2304:    Logically Collective

2306:    Input Arguments:
2307: +  ksp - the KSP context
2308: .  func - function to compute the right hand side
2309: -  ctx - optional context

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

2314: +  ksp - the KSP context
2315: .  b - right hand side of linear system
2316: -  ctx - optional user-provided context

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

2321:    Level: beginner

2323: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2324: @*/
2325: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2326: {
2328:   DM             dm;

2332:   KSPGetDM(ksp,&dm);
2333:   DMKSPSetComputeRHS(dm,func,ctx);
2334:   return(0);
2335: }

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

2340:    Logically Collective

2342:    Input Arguments:
2343: +  ksp - the KSP context
2344: .  func - function to compute the initial guess
2345: -  ctx - optional context

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

2350: +  ksp - the KSP context
2351: .  x - solution vector
2352: -  ctx - optional user-provided context

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

2357:    Level: beginner

2359: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2360: @*/
2361: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2362: {
2364:   DM             dm;

2368:   KSPGetDM(ksp,&dm);
2369:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2370:   return(0);
2371: }