Actual source code: itfunc.c

petsc-3.10.5 2019-03-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:   PCGetSetUpFailedReason(pc,&pcreason);
225:   if (pcreason) {
226:     ksp->reason = KSP_DIVERGED_PCSETUP_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:   PCGetSetUpFailedReason(ksp->pc,&pcreason);
393:   if (pcreason) {
394:     ksp->reason = KSP_DIVERGED_PCSETUP_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_PCSETUP_FAILED) {
432:         PCFailedReason reason;
433:         PCGetSetUpFailedReason(ksp->pc,&reason);
434:         PetscViewerASCIIPrintf(viewer,"               PCSETUP_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)->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:   if (ksp->presolve) {(*ksp->presolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->prectx);}
701:   PetscLogEventBegin(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);

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

707:   if (ksp->guess) {
708:     PetscObjectState ostate,state;

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

722:   /* KSPSetUp() scales the matrix if needed */
723:   KSPSetUp(ksp);
724:   KSPSetUpOnBlocks(ksp);

726:   VecLocked(ksp->vec_sol,3);

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

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

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

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

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

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

782:   VecLockPop(ksp->vec_rhs);
783:   if (nullsp) {
784:     ksp->vec_rhs = vec_rhs;
785:     VecDestroy(&btmp);
786:   }

788:   ksp->guess_zero = guess_zero;

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

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

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

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

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

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

841:     KSPComputeExplicitOperator(ksp, &B);
842:     ObjectView((PetscObject) B, ksp->viewerPOpExp, ksp->formatPOpExp);
843:     MatDestroy(&B);
844:   }

846:   if (inXisinB) {
847:     VecCopy(x,b);
848:     VecDestroy(&x);
849:   }
850:   PetscObjectSAWsBlock((PetscObject)ksp);
851:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(comm,PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
852:   return(0);
853: }

855: /*@
856:    KSPSolveTranspose - Solves the transpose of a linear system.

858:    Collective on KSP

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

865:    Notes:
866:     For complex numbers this solve the non-Hermitian transpose system.

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

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

873:    Level: developer

875: .keywords: solve, linear system

877: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
878:           KSPSolve(), KSP
879: @*/

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

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

902:   ksp->vec_rhs         = b;
903:   ksp->vec_sol         = x;
904:   ksp->transpose_solve = PETSC_TRUE;

906:   if (ksp->guess) {
907:     PetscObjectState ostate,state;

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

921:   KSPSetUp(ksp);
922:   KSPSetUpOnBlocks(ksp);
923:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}

925:   PCGetOperators(ksp->pc,&mat,&pmat);
926:   MatGetNullSpace(pmat,&nullsp);
927:   if (nullsp) {
928:     VecDuplicate(ksp->vec_rhs,&btmp);
929:     VecCopy(ksp->vec_rhs,btmp);
930:     MatNullSpaceRemove(nullsp,btmp);
931:     vec_rhs      = ksp->vec_rhs;
932:     ksp->vec_rhs = btmp;
933:   }

935:   (*ksp->ops->solve)(ksp);
936:   if (nullsp) {
937:     ksp->vec_rhs = vec_rhs;
938:     VecDestroy(&btmp);
939:   }
940:   if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
941:   if (ksp->viewReason) {KSPReasonView_Internal(ksp, ksp->viewerReason, ksp->formatReason);}
942:   if (ksp->guess) {
943:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
944:   }

946:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,          ksp->viewerMat,  ksp->formatMat);}
947:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,         ksp->viewerPMat, ksp->formatPMat);}
948:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs, ksp->viewerRhs,  ksp->formatRhs);}
949:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol, ksp->viewerSol,  ksp->formatSol);}

951:   if (inXisinB) {
952:     VecCopy(x,b);
953:     VecDestroy(&x);
954:   }
955:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
956:   return(0);
957: }

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

962:    Collective on KSP

964:    Input Parameter:
965: .  ksp - iterative context obtained from KSPCreate()

967:    Level: beginner

969: .keywords: destroy

971: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
972: @*/
973: PetscErrorCode  KSPReset(KSP ksp)
974: {

979:   if (!ksp) return(0);
980:   if (ksp->ops->reset) {
981:     (*ksp->ops->reset)(ksp);
982:   }
983:   if (ksp->pc) {PCReset(ksp->pc);}
984:   if (ksp->guess) {
985:     KSPGuess guess = ksp->guess;
986:     if (guess->ops->reset) { (*guess->ops->reset)(guess); }
987:   }
988:   VecDestroyVecs(ksp->nwork,&ksp->work);
989:   VecDestroy(&ksp->vec_rhs);
990:   VecDestroy(&ksp->vec_sol);
991:   VecDestroy(&ksp->diagonal);
992:   VecDestroy(&ksp->truediagonal);

994:   PetscViewerDestroy(&ksp->viewer);
995:   PetscViewerDestroy(&ksp->viewerPre);
996:   PetscViewerDestroy(&ksp->viewerReason);
997:   PetscViewerDestroy(&ksp->viewerMat);
998:   PetscViewerDestroy(&ksp->viewerPMat);
999:   PetscViewerDestroy(&ksp->viewerRhs);
1000:   PetscViewerDestroy(&ksp->viewerSol);
1001:   PetscViewerDestroy(&ksp->viewerMatExp);
1002:   PetscViewerDestroy(&ksp->viewerEV);
1003:   PetscViewerDestroy(&ksp->viewerSV);
1004:   PetscViewerDestroy(&ksp->viewerEVExp);
1005:   PetscViewerDestroy(&ksp->viewerFinalRes);
1006:   PetscViewerDestroy(&ksp->viewerPOpExp);

1008:   ksp->setupstage = KSP_SETUP_NEW;
1009:   return(0);
1010: }

1012: /*@
1013:    KSPDestroy - Destroys KSP context.

1015:    Collective on KSP

1017:    Input Parameter:
1018: .  ksp - iterative context obtained from KSPCreate()

1020:    Level: beginner

1022: .keywords: destroy

1024: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1025: @*/
1026: PetscErrorCode  KSPDestroy(KSP *ksp)
1027: {
1029:   PC             pc;

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

1036:   PetscObjectSAWsViewOff((PetscObject)*ksp);

1038:   /*
1039:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1040:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1041:    refcount (and may be shared, e.g., by other ksps).
1042:    */
1043:   pc         = (*ksp)->pc;
1044:   (*ksp)->pc = NULL;
1045:   KSPReset((*ksp));
1046:   (*ksp)->pc = pc;
1047:   if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

1049:   KSPGuessDestroy(&(*ksp)->guess);
1050:   DMDestroy(&(*ksp)->dm);
1051:   PCDestroy(&(*ksp)->pc);
1052:   PetscFree((*ksp)->res_hist_alloc);
1053:   if ((*ksp)->convergeddestroy) {
1054:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
1055:   }
1056:   KSPMonitorCancel((*ksp));
1057:   PetscViewerDestroy(&(*ksp)->eigviewer);
1058:   PetscHeaderDestroy(ksp);
1059:   return(0);
1060: }

1062: /*@
1063:     KSPSetPCSide - Sets the preconditioning side.

1065:     Logically Collective on KSP

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

1070:     Output Parameter:
1071: .   side - the preconditioning side, where side is one of
1072: .vb
1073:       PC_LEFT - left preconditioning (default)
1074:       PC_RIGHT - right preconditioning
1075:       PC_SYMMETRIC - symmetric preconditioning
1076: .ve

1078:     Options Database Keys:
1079: .   -ksp_pc_side <right,left,symmetric>

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

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

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

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

1092:     Level: intermediate

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

1096: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType(), KSP
1097: @*/
1098: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
1099: {
1103:   ksp->pc_side = ksp->pc_side_set = side;
1104:   return(0);
1105: }

1107: /*@
1108:     KSPGetPCSide - Gets the preconditioning side.

1110:     Not Collective

1112:     Input Parameter:
1113: .   ksp - iterative context obtained from KSPCreate()

1115:     Output Parameter:
1116: .   side - the preconditioning side, where side is one of
1117: .vb
1118:       PC_LEFT - left preconditioning (default)
1119:       PC_RIGHT - right preconditioning
1120:       PC_SYMMETRIC - symmetric preconditioning
1121: .ve

1123:     Level: intermediate

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

1127: .seealso: KSPSetPCSide(), KSP
1128: @*/
1129: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
1130: {

1136:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
1137:   *side = ksp->pc_side;
1138:   return(0);
1139: }

1141: /*@
1142:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1143:    iteration tolerances used by the default KSP convergence tests.

1145:    Not Collective

1147:    Input Parameter:
1148: .  ksp - the Krylov subspace context

1150:    Output Parameters:
1151: +  rtol - the relative convergence tolerance
1152: .  abstol - the absolute convergence tolerance
1153: .  dtol - the divergence tolerance
1154: -  maxits - maximum number of iterations

1156:    Notes:
1157:    The user can specify NULL for any parameter that is not needed.

1159:    Level: intermediate

1161: .keywords: get, tolerance, absolute, relative, divergence, convergence,
1162:            maximum, iterations

1164: .seealso: KSPSetTolerances(), KSP
1165: @*/
1166: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1167: {
1170:   if (abstol) *abstol = ksp->abstol;
1171:   if (rtol) *rtol = ksp->rtol;
1172:   if (dtol) *dtol = ksp->divtol;
1173:   if (maxits) *maxits = ksp->max_it;
1174:   return(0);
1175: }

1177: /*@
1178:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1179:    iteration tolerances used by the default KSP convergence testers.

1181:    Logically Collective on KSP

1183:    Input Parameters:
1184: +  ksp - the Krylov subspace context
1185: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1186: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1187: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1188: -  maxits - maximum number of iterations to use

1190:    Options Database Keys:
1191: +  -ksp_atol <abstol> - Sets abstol
1192: .  -ksp_rtol <rtol> - Sets rtol
1193: .  -ksp_divtol <dtol> - Sets dtol
1194: -  -ksp_max_it <maxits> - Sets maxits

1196:    Notes:
1197:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

1202:    Level: intermediate

1204: .keywords: set, tolerance, absolute, relative, divergence,
1205:            convergence, maximum, iterations

1207: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
1208: @*/
1209: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1210: {

1218:   if (rtol != PETSC_DEFAULT) {
1219:     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);
1220:     ksp->rtol = rtol;
1221:   }
1222:   if (abstol != PETSC_DEFAULT) {
1223:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1224:     ksp->abstol = abstol;
1225:   }
1226:   if (dtol != PETSC_DEFAULT) {
1227:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1228:     ksp->divtol = dtol;
1229:   }
1230:   if (maxits != PETSC_DEFAULT) {
1231:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1232:     ksp->max_it = maxits;
1233:   }
1234:   return(0);
1235: }

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

1242:    Logically Collective on KSP

1244:    Input Parameters:
1245: +  ksp - iterative context obtained from KSPCreate()
1246: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

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

1251:    Level: beginner

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

1256: .keywords: set, initial guess, nonzero

1258: .seealso: KSPGetInitialGuessNonzero(), KSPSetGuessType(), KSPGuessType, KSP
1259: @*/
1260: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1261: {
1265:   ksp->guess_zero = (PetscBool) !(int)flg;
1266:   return(0);
1267: }

1269: /*@
1270:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1271:    a zero initial guess.

1273:    Not Collective

1275:    Input Parameter:
1276: .  ksp - iterative context obtained from KSPCreate()

1278:    Output Parameter:
1279: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1281:    Level: intermediate

1283: .keywords: set, initial guess, nonzero

1285: .seealso: KSPSetInitialGuessNonzero(), KSP
1286: @*/
1287: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1288: {
1292:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1293:   else *flag = PETSC_TRUE;
1294:   return(0);
1295: }

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

1300:    Logically Collective on KSP

1302:    Input Parameters:
1303: +  ksp - iterative context obtained from KSPCreate()
1304: -  flg - PETSC_TRUE indicates you want the error generated

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

1309:    Level: intermediate

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


1316: .seealso: KSPGetErrorIfNotConverged(), KSP
1317: @*/
1318: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1319: {
1323:   ksp->errorifnotconverged = flg;
1324:   return(0);
1325: }

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

1330:    Not Collective

1332:    Input Parameter:
1333: .  ksp - iterative context obtained from KSPCreate()

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

1338:    Level: intermediate

1340: .seealso: KSPSetErrorIfNotConverged(), KSP
1341: @*/
1342: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1343: {
1347:   *flag = ksp->errorifnotconverged;
1348:   return(0);
1349: }

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

1354:    Logically Collective on KSP

1356:    Input Parameters:
1357: +  ksp - iterative context obtained from KSPCreate()
1358: -  flg - PETSC_TRUE or PETSC_FALSE

1360:    Level: advanced

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

1364: .keywords: set, initial guess, nonzero

1366: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1367: @*/
1368: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1369: {
1373:   ksp->guess_knoll = flg;
1374:   return(0);
1375: }

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

1381:    Not Collective

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

1386:    Output Parameter:
1387: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1389:    Level: advanced

1391: .keywords: set, initial guess, nonzero

1393: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1394: @*/
1395: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1396: {
1400:   *flag = ksp->guess_knoll;
1401:   return(0);
1402: }

1404: /*@
1405:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1406:    values will be calculated via a Lanczos or Arnoldi process as the linear
1407:    system is solved.

1409:    Not Collective

1411:    Input Parameter:
1412: .  ksp - iterative context obtained from KSPCreate()

1414:    Output Parameter:
1415: .  flg - PETSC_TRUE or PETSC_FALSE

1417:    Options Database Key:
1418: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1420:    Notes:
1421:    Currently this option is not valid for all iterative methods.

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

1427:    Level: advanced

1429: .keywords: set, compute, singular values

1431: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1432: @*/
1433: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1434: {
1438:   *flg = ksp->calc_sings;
1439:   return(0);
1440: }

1442: /*@
1443:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1444:    values will be calculated via a Lanczos or Arnoldi process as the linear
1445:    system is solved.

1447:    Logically Collective on KSP

1449:    Input Parameters:
1450: +  ksp - iterative context obtained from KSPCreate()
1451: -  flg - PETSC_TRUE or PETSC_FALSE

1453:    Options Database Key:
1454: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1456:    Notes:
1457:    Currently this option is not valid for all iterative methods.

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

1463:    Level: advanced

1465: .keywords: set, compute, singular values

1467: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1468: @*/
1469: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1470: {
1474:   ksp->calc_sings = flg;
1475:   return(0);
1476: }

1478: /*@
1479:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1480:    values will be calculated via a Lanczos or Arnoldi process as the linear
1481:    system is solved.

1483:    Not Collective

1485:    Input Parameter:
1486: .  ksp - iterative context obtained from KSPCreate()

1488:    Output Parameter:
1489: .  flg - PETSC_TRUE or PETSC_FALSE

1491:    Notes:
1492:    Currently this option is not valid for all iterative methods.

1494:    Level: advanced

1496: .keywords: set, compute, eigenvalues

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

1509: /*@
1510:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1511:    values will be calculated via a Lanczos or Arnoldi process as the linear
1512:    system is solved.

1514:    Logically Collective on KSP

1516:    Input Parameters:
1517: +  ksp - iterative context obtained from KSPCreate()
1518: -  flg - PETSC_TRUE or PETSC_FALSE

1520:    Notes:
1521:    Currently this option is not valid for all iterative methods.

1523:    Level: advanced

1525: .keywords: set, compute, eigenvalues

1527: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1528: @*/
1529: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1530: {
1534:   ksp->calc_sings = flg;
1535:   return(0);
1536: }

1538: /*@
1539:    KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1540:    will be calculated via a Lanczos or Arnoldi process as the linear
1541:    system is solved.

1543:    Logically Collective on KSP

1545:    Input Parameters:
1546: +  ksp - iterative context obtained from KSPCreate()
1547: -  flg - PETSC_TRUE or PETSC_FALSE

1549:    Notes:
1550:    Currently this option is only valid for the GMRES method.

1552:    Level: advanced

1554: .keywords: set, compute, ritz

1556: .seealso: KSPComputeRitz(), KSP
1557: @*/
1558: PetscErrorCode  KSPSetComputeRitz(KSP ksp, PetscBool flg)
1559: {
1563:   ksp->calc_ritz = flg;
1564:   return(0);
1565: }

1567: /*@
1568:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1569:    be solved.

1571:    Not Collective

1573:    Input Parameter:
1574: .  ksp - iterative context obtained from KSPCreate()

1576:    Output Parameter:
1577: .  r - right-hand-side vector

1579:    Level: developer

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

1583: .seealso: KSPGetSolution(), KSPSolve(), KSP
1584: @*/
1585: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1586: {
1590:   *r = ksp->vec_rhs;
1591:   return(0);
1592: }

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

1599:    Not Collective

1601:    Input Parameters:
1602: .  ksp - iterative context obtained from KSPCreate()

1604:    Output Parameters:
1605: .  v - solution vector

1607:    Level: developer

1609: .keywords: get, solution

1611: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve(), KSP
1612: @*/
1613: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1614: {
1618:   *v = ksp->vec_sol;
1619:   return(0);
1620: }

1622: /*@
1623:    KSPSetPC - Sets the preconditioner to be used to calculate the
1624:    application of the preconditioner on a vector.

1626:    Collective on KSP

1628:    Input Parameters:
1629: +  ksp - iterative context obtained from KSPCreate()
1630: -  pc   - the preconditioner object

1632:    Notes:
1633:    Use KSPGetPC() to retrieve the preconditioner context (for example,
1634:    to free it at the end of the computations).

1636:    Level: developer

1638: .keywords: set, precondition, Binv

1640: .seealso: KSPGetPC(), KSP
1641: @*/
1642: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1643: {

1650:   PetscObjectReference((PetscObject)pc);
1651:   PCDestroy(&ksp->pc);
1652:   ksp->pc = pc;
1653:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1654:   return(0);
1655: }

1657: /*@
1658:    KSPGetPC - Returns a pointer to the preconditioner context
1659:    set with KSPSetPC().

1661:    Not Collective

1663:    Input Parameters:
1664: .  ksp - iterative context obtained from KSPCreate()

1666:    Output Parameter:
1667: .  pc - preconditioner context

1669:    Level: developer

1671: .keywords: get, preconditioner, Binv

1673: .seealso: KSPSetPC(), KSP
1674: @*/
1675: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1676: {

1682:   if (!ksp->pc) {
1683:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1684:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1685:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1686:   }
1687:   *pc = ksp->pc;
1688:   return(0);
1689: }

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

1694:    Collective on KSP

1696:    Input Parameters:
1697: +  ksp - iterative context obtained from KSPCreate()
1698: .  it - iteration number
1699: -  rnorm - relative norm of the residual

1701:    Notes:
1702:    This routine is called by the KSP implementations.
1703:    It does not typically need to be called by the user.

1705:    Level: developer

1707: .seealso: KSPMonitorSet()
1708: @*/
1709: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1710: {
1711:   PetscInt       i, n = ksp->numbermonitors;

1715:   for (i=0; i<n; i++) {
1716:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1717:   }
1718:   return(0);
1719: }

1721: /*

1723:     Checks if two monitors are identical; if they are then it destroys the new one
1724: */
1725: PetscErrorCode PetscMonitorCompare(PetscErrorCode (*nmon)(void),void *nmctx,PetscErrorCode (*nmdestroy)(void**),PetscErrorCode (*mon)(void),void *mctx,PetscErrorCode (*mdestroy)(void**),PetscBool *identical)
1726: {
1727:   *identical = PETSC_FALSE;
1728:   if (nmon == mon && nmdestroy == mdestroy) {
1729:     if (nmctx == mctx) *identical = PETSC_TRUE;
1730:     else if (nmdestroy == (PetscErrorCode (*)(void**)) PetscViewerAndFormatDestroy) {
1731:       PetscViewerAndFormat *old = (PetscViewerAndFormat*)mctx, *newo = (PetscViewerAndFormat*)nmctx;
1732:       if (old->viewer == newo->viewer && old->format == newo->format) *identical = PETSC_TRUE;
1733:     }
1734:     if (*identical) {
1735:       if (mdestroy) {
1737:         (*mdestroy)(&nmctx);
1738:       }
1739:     }
1740:   }
1741:   return(0);
1742: }

1744: /*@C
1745:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1746:    the residual/error etc.

1748:    Logically Collective on KSP

1750:    Input Parameters:
1751: +  ksp - iterative context obtained from KSPCreate()
1752: .  monitor - pointer to function (if this is NULL, it turns off monitoring
1753: .  mctx    - [optional] context for private data for the
1754:              monitor routine (use NULL if no context is desired)
1755: -  monitordestroy - [optional] routine that frees monitor context
1756:           (may be NULL)

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

1761: +  ksp - iterative context obtained from KSPCreate()
1762: .  it - iteration number
1763: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1764: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

1766:    Options Database Keys:
1767: +    -ksp_monitor        - sets KSPMonitorDefault()
1768: .    -ksp_monitor_true_residual    - sets KSPMonitorTrueResidualNorm()
1769: .    -ksp_monitor_max    - sets KSPMonitorTrueResidualMaxNorm()
1770: .    -ksp_monitor_lg_residualnorm    - sets line graph monitor,
1771:                            uses KSPMonitorLGResidualNormCreate()
1772: .    -ksp_monitor_lg_true_residualnorm   - sets line graph monitor,
1773:                            uses KSPMonitorLGResidualNormCreate()
1774: .    -ksp_monitor_singular_value    - sets KSPMonitorSingularValue()
1775: -    -ksp_monitor_cancel - cancels all monitors that have
1776:                           been hardwired into a code by
1777:                           calls to KSPMonitorSet(), but
1778:                           does not cancel those set via
1779:                           the options database.

1781:    Notes:
1782:    The default is to do nothing.  To print the residual, or preconditioned
1783:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1784:    KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1785:    context.

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

1791:    Fortran Notes:
1792:     Only a single monitor function can be set for each KSP object

1794:    Level: beginner

1796: .keywords: set, monitor

1798: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel(), KSP
1799: @*/
1800: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1801: {
1802:   PetscInt       i;
1804:   PetscBool      identical;

1808:   for (i=0; i<ksp->numbermonitors;i++) {
1809:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
1810:     if (identical) return(0);
1811:   }
1812:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1813:   ksp->monitor[ksp->numbermonitors]          = monitor;
1814:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
1815:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1816:   return(0);
1817: }

1819: /*@
1820:    KSPMonitorCancel - Clears all monitors for a KSP object.

1822:    Logically Collective on KSP

1824:    Input Parameters:
1825: .  ksp - iterative context obtained from KSPCreate()

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

1832:    Level: intermediate

1834: .keywords: set, monitor

1836: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet(), KSP
1837: @*/
1838: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1839: {
1841:   PetscInt       i;

1845:   for (i=0; i<ksp->numbermonitors; i++) {
1846:     if (ksp->monitordestroy[i]) {
1847:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1848:     }
1849:   }
1850:   ksp->numbermonitors = 0;
1851:   return(0);
1852: }

1854: /*@C
1855:    KSPGetMonitorContext - Gets the monitoring context, as set by
1856:    KSPMonitorSet() for the FIRST monitor only.

1858:    Not Collective

1860:    Input Parameter:
1861: .  ksp - iterative context obtained from KSPCreate()

1863:    Output Parameter:
1864: .  ctx - monitoring context

1866:    Level: intermediate

1868: .keywords: get, monitor, context

1870: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSP
1871: @*/
1872: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1873: {
1876:   *ctx =      (ksp->monitorcontext[0]);
1877:   return(0);
1878: }

1880: /*@
1881:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1882:    If set, this array will contain the residual norms computed at each
1883:    iteration of the solver.

1885:    Not Collective

1887:    Input Parameters:
1888: +  ksp - iterative context obtained from KSPCreate()
1889: .  a   - array to hold history
1890: .  na  - size of a
1891: -  reset - PETSC_TRUE indicates the history counter is reset to zero
1892:            for each new linear solve

1894:    Level: advanced

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

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

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

1905: .seealso: KSPGetResidualHistory(), KSP

1907: @*/
1908: PetscErrorCode  KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
1909: {


1915:   PetscFree(ksp->res_hist_alloc);
1916:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
1917:     ksp->res_hist     = a;
1918:     ksp->res_hist_max = na;
1919:   } else {
1920:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
1921:     else                                           ksp->res_hist_max = 10000; /* like default ksp->max_it */
1922:     PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);

1924:     ksp->res_hist = ksp->res_hist_alloc;
1925:   }
1926:   ksp->res_hist_len   = 0;
1927:   ksp->res_hist_reset = reset;
1928:   return(0);
1929: }

1931: /*@C
1932:    KSPGetResidualHistory - Gets the array used to hold the residual history
1933:    and the number of residuals it contains.

1935:    Not Collective

1937:    Input Parameter:
1938: .  ksp - iterative context obtained from KSPCreate()

1940:    Output Parameters:
1941: +  a   - pointer to array to hold history (or NULL)
1942: -  na  - number of used entries in a (or NULL)

1944:    Level: advanced

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

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

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

1957: .seealso: KSPGetResidualHistory(), KSP

1959: @*/
1960: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1961: {
1964:   if (a) *a = ksp->res_hist;
1965:   if (na) *na = ksp->res_hist_len;
1966:   return(0);
1967: }

1969: /*@C
1970:    KSPSetConvergenceTest - Sets the function to be used to determine
1971:    convergence.

1973:    Logically Collective on KSP

1975:    Input Parameters:
1976: +  ksp - iterative context obtained from KSPCreate()
1977: .  converge - pointer to int function
1978: .  cctx    - context for private data for the convergence routine (may be null)
1979: -  destroy - a routine for destroying the context (may be null)

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

1984: +  ksp - iterative context obtained from KSPCreate()
1985: .  it - iteration number
1986: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1987: .  reason - the reason why it has converged or diverged
1988: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


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

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

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

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

2005:    Level: advanced

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

2009: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP
2010: @*/
2011: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2012: {

2017:   if (ksp->convergeddestroy) {
2018:     (*ksp->convergeddestroy)(ksp->cnvP);
2019:   }
2020:   ksp->converged        = converge;
2021:   ksp->convergeddestroy = destroy;
2022:   ksp->cnvP             = (void*)cctx;
2023:   return(0);
2024: }

2026: /*@C
2027:    KSPGetConvergenceContext - Gets the convergence context set with
2028:    KSPSetConvergenceTest().

2030:    Not Collective

2032:    Input Parameter:
2033: .  ksp - iterative context obtained from KSPCreate()

2035:    Output Parameter:
2036: .  ctx - monitoring context

2038:    Level: advanced

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

2042: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
2043: @*/
2044: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
2045: {
2048:   *ctx = ksp->cnvP;
2049:   return(0);
2050: }

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

2056:    Collective on KSP

2058:    Input Parameter:
2059: .  ctx - iterative context obtained from KSPCreate()

2061:    Output Parameter:
2062:    Provide exactly one of
2063: +  v - location to stash solution.
2064: -  V - the solution is returned in this location. This vector is created
2065:        internally. This vector should NOT be destroyed by the user with
2066:        VecDestroy().

2068:    Notes:
2069:    This routine can be used in one of two ways
2070: .vb
2071:       KSPBuildSolution(ksp,NULL,&V);
2072:    or
2073:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2074: .ve
2075:    In the first case an internal vector is allocated to store the solution
2076:    (the user cannot destroy this vector). In the second case the solution
2077:    is generated in the vector that the user provides. Note that for certain
2078:    methods, such as KSPCG, the second case requires a copy of the solution,
2079:    while in the first case the call is essentially free since it simply
2080:    returns the vector where the solution already is stored. For some methods
2081:    like GMRES this is a reasonably expensive operation and should only be
2082:    used in truly needed.

2084:    Level: advanced

2086: .keywords: build, solution

2088: .seealso: KSPGetSolution(), KSPBuildResidual(), KSP
2089: @*/
2090: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2091: {

2096:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2097:   if (!V) V = &v;
2098:   (*ksp->ops->buildsolution)(ksp,v,V);
2099:   return(0);
2100: }

2102: /*@C
2103:    KSPBuildResidual - Builds the residual in a vector provided.

2105:    Collective on KSP

2107:    Input Parameter:
2108: .  ksp - iterative context obtained from KSPCreate()

2110:    Output Parameters:
2111: +  v - optional location to stash residual.  If v is not provided,
2112:        then a location is generated.
2113: .  t - work vector.  If not provided then one is generated.
2114: -  V - the residual

2116:    Notes:
2117:    Regardless of whether or not v is provided, the residual is
2118:    returned in V.

2120:    Level: advanced

2122: .keywords: KSP, build, residual

2124: .seealso: KSPBuildSolution()
2125: @*/
2126: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2127: {
2129:   PetscBool      flag = PETSC_FALSE;
2130:   Vec            w    = v,tt = t;

2134:   if (!w) {
2135:     VecDuplicate(ksp->vec_rhs,&w);
2136:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2137:   }
2138:   if (!tt) {
2139:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2140:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2141:   }
2142:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
2143:   if (flag) {VecDestroy(&tt);}
2144:   return(0);
2145: }

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

2151:    Logically Collective on KSP

2153:    Input Parameter:
2154: +  ksp - the KSP context
2155: -  scale - PETSC_TRUE or PETSC_FALSE

2157:    Options Database Key:
2158: +   -ksp_diagonal_scale -
2159: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


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

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

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

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

2177:    Level: intermediate

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

2181: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2182: @*/
2183: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2184: {
2188:   ksp->dscale = scale;
2189:   return(0);
2190: }

2192: /*@
2193:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2194:                           right hand side

2196:    Not Collective

2198:    Input Parameter:
2199: .  ksp - the KSP context

2201:    Output Parameter:
2202: .  scale - PETSC_TRUE or PETSC_FALSE

2204:    Notes:
2205:     BE CAREFUL with this routine: it actually scales the matrix and right
2206:     hand side that define the system. After the system is solved the matrix
2207:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

2209:    Level: intermediate

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

2213: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2214: @*/
2215: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
2216: {
2220:   *scale = ksp->dscale;
2221:   return(0);
2222: }

2224: /*@
2225:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2226:      back after solving.

2228:    Logically Collective on KSP

2230:    Input Parameter:
2231: +  ksp - the KSP context
2232: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2233:          rescale (default)

2235:    Notes:
2236:      Must be called after KSPSetDiagonalScale()

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

2243:    Level: intermediate

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

2247: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix(), KSP
2248: @*/
2249: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2250: {
2254:   ksp->dscalefix = fix;
2255:   return(0);
2256: }

2258: /*@
2259:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2260:      back after solving.

2262:    Not Collective

2264:    Input Parameter:
2265: .  ksp - the KSP context

2267:    Output Parameter:
2268: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2269:          rescale (default)

2271:    Notes:
2272:      Must be called after KSPSetDiagonalScale()

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

2279:    Level: intermediate

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

2283: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2284: @*/
2285: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2286: {
2290:   *fix = ksp->dscalefix;
2291:   return(0);
2292: }

2294: /*@C
2295:    KSPSetComputeOperators - set routine to compute the linear operators

2297:    Logically Collective

2299:    Input Arguments:
2300: +  ksp - the KSP context
2301: .  func - function to compute the operators
2302: -  ctx - optional context

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

2307: +  ksp - the KSP context
2308: .  A - the linear operator
2309: .  B - preconditioning matrix
2310: -  ctx - optional user-provided context

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

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

2318:    Level: beginner

2320: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2321: @*/
2322: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2323: {
2325:   DM             dm;

2329:   KSPGetDM(ksp,&dm);
2330:   DMKSPSetComputeOperators(dm,func,ctx);
2331:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2332:   return(0);
2333: }

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

2338:    Logically Collective

2340:    Input Arguments:
2341: +  ksp - the KSP context
2342: .  func - function to compute the right hand side
2343: -  ctx - optional context

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

2348: +  ksp - the KSP context
2349: .  b - right hand side of linear system
2350: -  ctx - optional user-provided context

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

2355:    Level: beginner

2357: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2358: @*/
2359: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2360: {
2362:   DM             dm;

2366:   KSPGetDM(ksp,&dm);
2367:   DMKSPSetComputeRHS(dm,func,ctx);
2368:   return(0);
2369: }

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

2374:    Logically Collective

2376:    Input Arguments:
2377: +  ksp - the KSP context
2378: .  func - function to compute the initial guess
2379: -  ctx - optional context

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

2384: +  ksp - the KSP context
2385: .  x - solution vector
2386: -  ctx - optional user-provided context

2388:    Level: beginner

2390: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2391: @*/
2392: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2393: {
2395:   DM             dm;

2399:   KSPGetDM(ksp,&dm);
2400:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2401:   return(0);
2402: }