Actual source code: itfunc.c

petsc-master 2019-10-23
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: .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_compute_eigenvalues - Prints eigenvalues to stdout
 90: -  -ksp_plot_eigenvalues - Plots eigenvalues in an x-window display

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

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

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

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

111:    Level: advanced

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

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

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

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

140:    Not Collective

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

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

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

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

169:    Level: advanced

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

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

188:    Collective on ksp

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

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

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

202:    Level: advanced

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

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

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

226:    Collective on ksp

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

232:    Level: intermediate

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

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

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

251:    Collective on ksp

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

257:    Level: intermediate

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

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

273:    Collective on ksp

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

278:    Level: developer

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

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

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

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

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

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

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

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

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

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

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

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

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

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

432:    Collective on ksp

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


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

442:    Level: beginner

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

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

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

463:   Collective on ksp

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

468:   Level: intermediate

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

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

487:  #include <petscdraw.h>

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

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

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

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

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

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

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

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

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

582: /*@
583:    KSPSolve - Solves linear system.

585:    Collective on ksp

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

592:    Options Database Keys:
593: +  -ksp_compute_eigenvalues - compute preconditioned operators eigenvalues
594: .  -ksp_plot_eigenvalues - plot the computed eigenvalues in an X-window
595: .  -ksp_plot_eigencontours - plot the computed eigenvalues in an X-window with contours
596: .  -ksp_compute_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and using LAPACK
597: .  -ksp_plot_eigenvalues_explicitly - plot the explicitly computing eigenvalues
598: .  -ksp_view_mat binary - save matrix to the default binary viewer
599: .  -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
600: .  -ksp_view_rhs binary - save right hand side vector to the default binary viewer
601: .  -ksp_view_solution binary - save computed solution vector to the default binary viewer
602:            (can be read later with src/ksp/examples/tutorials/ex10.c for testing solvers)
603: .  -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
604: .  -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
605: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
606: .  -ksp_view_final_residual - print 2-norm of true linear system residual at the end of the solution process
607: -  -ksp_view - print the ksp data structure at the end of the system solution

609:    Notes:

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

613:    The operator is specified with KSPSetOperators().

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

618:    If you provide a matrix that has a MatSetNullSpace() and MatSetTransposeNullSpace() this will use that information to solve singular systems
619:    in the least squares sense with a norm minimizing solution.
620: $
621: $                   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()
622: $
623: $    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
624: $    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
625: $    direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
626: $
627: $    We recommend always using GMRES for such singular systems.
628: $    If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
629: $    If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).

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


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

641:    Understanding Convergence:
642:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
643:    KSPComputeEigenvaluesExplicitly() provide information on additional
644:    options to monitor convergence and print eigenvalue information.

646:    Level: beginner

648: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
649:           KSPSolveTranspose(), KSPGetIterationNumber(), MatNullSpaceCreate(), MatSetNullSpace(), MatSetTransposeNullSpace(), KSP
650: @*/
651: PetscErrorCode KSPSolve(KSP ksp,Vec b,Vec x)
652: {
653:   PetscErrorCode    ierr;
654:   PetscBool         flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
655:   Mat               mat,pmat;
656:   MPI_Comm          comm;
657:   MatNullSpace      nullsp;
658:   Vec               btmp,vec_rhs=0;

664:   comm = PetscObjectComm((PetscObject)ksp);
665:   if (x && x == b) {
666:     if (!ksp->guess_zero) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
667:     VecDuplicate(b,&x);
668:     inXisinB = PETSC_TRUE;
669:   }
670:   if (b) {
671:     PetscObjectReference((PetscObject)b);
672:     VecDestroy(&ksp->vec_rhs);
673:     ksp->vec_rhs = b;
674:   }
675:   if (x) {
676:     PetscObjectReference((PetscObject)x);
677:     VecDestroy(&ksp->vec_sol);
678:     ksp->vec_sol = x;
679:   }
680:   if (ksp->viewPre) {ObjectView((PetscObject) ksp, ksp->viewerPre, ksp->formatPre);}

682:   ksp->transpose_solve = PETSC_FALSE;

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

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

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

691:   if (ksp->guess) {
692:     PetscObjectState ostate,state;

694:     KSPGuessSetUp(ksp->guess);
695:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
696:     KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
697:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
698:     if (state != ostate) {
699:       ksp->guess_zero = PETSC_FALSE;
700:     } else {
701:       PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
702:       ksp->guess_zero = PETSC_TRUE;
703:     }
704:   }

706:   /* KSPSetUp() scales the matrix if needed */
707:   KSPSetUp(ksp);
708:   KSPSetUpOnBlocks(ksp);

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

712:   PCGetOperators(ksp->pc,&mat,&pmat);
713:   /* diagonal scale RHS if called for */
714:   if (ksp->dscale) {
715:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
716:     /* second time in, but matrix was scaled back to original */
717:     if (ksp->dscalefix && ksp->dscalefix2) {
718:       Mat mat,pmat;

720:       PCGetOperators(ksp->pc,&mat,&pmat);
721:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
722:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
723:     }

725:     /* scale initial guess */
726:     if (!ksp->guess_zero) {
727:       if (!ksp->truediagonal) {
728:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
729:         VecCopy(ksp->diagonal,ksp->truediagonal);
730:         VecReciprocal(ksp->truediagonal);
731:       }
732:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
733:     }
734:   }
735:   PCPreSolve(ksp->pc,ksp);

737:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
738:   if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
739:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
740:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
741:     ksp->guess_zero = PETSC_FALSE;
742:   }

744:   /* can we mark the initial guess as zero for this solve? */
745:   guess_zero = ksp->guess_zero;
746:   if (!ksp->guess_zero) {
747:     PetscReal norm;

749:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
750:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
751:   }
752:   MatGetTransposeNullSpace(pmat,&nullsp);
753:   if (nullsp) {
754:     VecDuplicate(ksp->vec_rhs,&btmp);
755:     VecCopy(ksp->vec_rhs,btmp);
756:     MatNullSpaceRemove(nullsp,btmp);
757:     vec_rhs      = ksp->vec_rhs;
758:     ksp->vec_rhs = btmp;
759:   }
760:   VecLockReadPush(ksp->vec_rhs);
761:   if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
762:     VecSetInf(ksp->vec_sol);
763:   }
764:   (*ksp->ops->solve)(ksp);

766:   VecLockReadPop(ksp->vec_rhs);
767:   if (nullsp) {
768:     ksp->vec_rhs = vec_rhs;
769:     VecDestroy(&btmp);
770:   }

772:   ksp->guess_zero = guess_zero;

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

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

780:   /* diagonal scale solution if called for */
781:   if (ksp->dscale) {
782:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
783:     /* unscale right hand side and matrix */
784:     if (ksp->dscalefix) {
785:       Mat mat,pmat;

787:       VecReciprocal(ksp->diagonal);
788:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
789:       PCGetOperators(ksp->pc,&mat,&pmat);
790:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
791:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
792:       VecReciprocal(ksp->diagonal);
793:       ksp->dscalefix2 = PETSC_TRUE;
794:     }
795:   }
796:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
797:   if (ksp->guess) {
798:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
799:   }
800:   if (ksp->postsolve) {
801:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
802:   }

804:   PCGetOperators(ksp->pc,&mat,&pmat);
805:   if (ksp->viewEV)       {KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV,    ksp->formatEV);}
806:   if (ksp->viewEVExp)    {KSPViewEigenvalues_Internal(ksp, PETSC_TRUE,  ksp->viewerEVExp, ksp->formatEVExp);}
807:   if (ksp->viewSV)       {KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV);}
808:   if (ksp->viewFinalRes) {KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes);}
809:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,           ksp->viewerMat,    ksp->formatMat);}
810:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,          ksp->viewerPMat,   ksp->formatPMat);}
811:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs,  ksp->viewerRhs,    ksp->formatRhs);}
812:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol,  ksp->viewerSol,    ksp->formatSol);}
813:   if (ksp->view)         {ObjectView((PetscObject) ksp,           ksp->viewer,       ksp->format);}
814:   if (ksp->viewDScale)   {ObjectView((PetscObject) ksp->diagonal, ksp->viewerDScale, ksp->formatDScale);}
815:   if (ksp->viewMatExp)   {
816:     Mat A, B;

818:     PCGetOperators(ksp->pc, &A, NULL);
819:     MatComputeOperator(A, MATAIJ, &B);
820:     ObjectView((PetscObject) B, ksp->viewerMatExp, ksp->formatMatExp);
821:     MatDestroy(&B);
822:   }
823:   if (ksp->viewPOpExp)   {
824:     Mat B;

826:     KSPComputeOperator(ksp, MATAIJ, &B);
827:     ObjectView((PetscObject) B, ksp->viewerPOpExp, ksp->formatPOpExp);
828:     MatDestroy(&B);
829:   }

831:   if (inXisinB) {
832:     VecCopy(x,b);
833:     VecDestroy(&x);
834:   }
835:   PetscObjectSAWsBlock((PetscObject)ksp);
836:   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]);
837:   return(0);
838: }

840: /*@
841:    KSPSolveTranspose - Solves the transpose of a linear system.

843:    Collective on ksp

845:    Input Parameter:
846: +  ksp - iterative context obtained from KSPCreate()
847: .  b - right hand side vector
848: -  x - solution vector

850:    Notes:
851:     For complex numbers this solve the non-Hermitian transpose system.

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

855:    Developer Notes:
856:     We need to implement a KSPSolveHermitianTranspose()

858:    Level: developer

860: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
861:           KSPSolve(), KSP
862: @*/

864: PetscErrorCode  KSPSolveTranspose(KSP ksp,Vec b,Vec x)
865: {
867:   PetscBool      inXisinB=PETSC_FALSE;
868:   Vec            vec_rhs = 0,btmp;
869:   Mat            mat,pmat;
870:   MatNullSpace   nullsp;

876:   if (x == b) {
877:     VecDuplicate(b,&x);
878:     inXisinB = PETSC_TRUE;
879:   }
880:   PetscObjectReference((PetscObject)b);
881:   PetscObjectReference((PetscObject)x);
882:   VecDestroy(&ksp->vec_rhs);
883:   VecDestroy(&ksp->vec_sol);

885:   ksp->vec_rhs         = b;
886:   ksp->vec_sol         = x;
887:   ksp->transpose_solve = PETSC_TRUE;

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

891:   PetscLogEventBegin(KSP_SolveTranspose,ksp,ksp->vec_rhs,ksp->vec_sol,0);
892:   if (ksp->guess) {
893:     PetscObjectState ostate,state;

895:     KSPGuessSetUp(ksp->guess);
896:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
897:     KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
898:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
899:     if (state != ostate) {
900:       ksp->guess_zero = PETSC_FALSE;
901:     } else {
902:       PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
903:       ksp->guess_zero = PETSC_TRUE;
904:     }
905:   }

907:   KSPSetUp(ksp);
908:   KSPSetUpOnBlocks(ksp);
909:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}

911:   PCGetOperators(ksp->pc,&mat,&pmat);
912:   MatGetNullSpace(pmat,&nullsp);
913:   if (nullsp) {
914:     VecDuplicate(ksp->vec_rhs,&btmp);
915:     VecCopy(ksp->vec_rhs,btmp);
916:     MatNullSpaceRemove(nullsp,btmp);
917:     vec_rhs      = ksp->vec_rhs;
918:     ksp->vec_rhs = btmp;
919:   }

921:   (*ksp->ops->solve)(ksp);
922:   ksp->totalits += ksp->its;
923:   if (nullsp) {
924:     ksp->vec_rhs = vec_rhs;
925:     VecDestroy(&btmp);
926:   }
927:   if (!ksp->reason) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason");
928:   if (ksp->viewReason) {KSPReasonView_Internal(ksp, ksp->viewerReason, ksp->formatReason);}
929:   PetscLogEventEnd(KSP_SolveTranspose,ksp,ksp->vec_rhs,ksp->vec_sol,0);
930:   if (ksp->guess) {
931:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
932:   }
933:   if (ksp->postsolve) {
934:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
935:   }

937:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,          ksp->viewerMat,  ksp->formatMat);}
938:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,         ksp->viewerPMat, ksp->formatPMat);}
939:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs, ksp->viewerRhs,  ksp->formatRhs);}
940:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol, ksp->viewerSol,  ksp->formatSol);}
941:   if (ksp->view)         {ObjectView((PetscObject) ksp,          ksp->viewer,     ksp->format);}

943:   if (inXisinB) {
944:     VecCopy(x,b);
945:     VecDestroy(&x);
946:   }
947:   if (ksp->errorifnotconverged && ksp->reason < 0) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_NOT_CONVERGED,"KSPSolve has not converged");
948:   return(0);
949: }

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

954:    Collective on ksp

956:    Input Parameter:
957: .  ksp - iterative context obtained from KSPCreate()

959:    Level: beginner

961: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSPSetFromOptions(), KSP
962: @*/
963: PetscErrorCode  KSPResetViewers(KSP ksp)
964: {

969:   if (!ksp) return(0);
970:   PetscViewerDestroy(&ksp->viewer);
971:   PetscViewerDestroy(&ksp->viewerPre);
972:   PetscViewerDestroy(&ksp->viewerReason);
973:   PetscViewerDestroy(&ksp->viewerMat);
974:   PetscViewerDestroy(&ksp->viewerPMat);
975:   PetscViewerDestroy(&ksp->viewerRhs);
976:   PetscViewerDestroy(&ksp->viewerSol);
977:   PetscViewerDestroy(&ksp->viewerMatExp);
978:   PetscViewerDestroy(&ksp->viewerEV);
979:   PetscViewerDestroy(&ksp->viewerSV);
980:   PetscViewerDestroy(&ksp->viewerEVExp);
981:   PetscViewerDestroy(&ksp->viewerFinalRes);
982:   PetscViewerDestroy(&ksp->viewerPOpExp);
983:   PetscViewerDestroy(&ksp->viewerDScale);
984:   ksp->view         = PETSC_FALSE;
985:   ksp->viewPre      = PETSC_FALSE;
986:   ksp->viewReason   = PETSC_FALSE;
987:   ksp->viewMat      = PETSC_FALSE;
988:   ksp->viewPMat     = PETSC_FALSE;
989:   ksp->viewRhs      = PETSC_FALSE;
990:   ksp->viewSol      = PETSC_FALSE;
991:   ksp->viewMatExp   = PETSC_FALSE;
992:   ksp->viewEV       = PETSC_FALSE;
993:   ksp->viewSV       = PETSC_FALSE;
994:   ksp->viewEVExp    = PETSC_FALSE;
995:   ksp->viewFinalRes = PETSC_FALSE;
996:   ksp->viewPOpExp   = PETSC_FALSE;
997:   ksp->viewDScale   = PETSC_FALSE;
998:   return(0);
999: }

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

1004:    Collective on ksp

1006:    Input Parameter:
1007: .  ksp - iterative context obtained from KSPCreate()

1009:    Level: beginner

1011: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1012: @*/
1013: PetscErrorCode  KSPReset(KSP ksp)
1014: {

1019:   if (!ksp) return(0);
1020:   if (ksp->ops->reset) {
1021:     (*ksp->ops->reset)(ksp);
1022:   }
1023:   if (ksp->pc) {PCReset(ksp->pc);}
1024:   if (ksp->guess) {
1025:     KSPGuess guess = ksp->guess;
1026:     if (guess->ops->reset) { (*guess->ops->reset)(guess); }
1027:   }
1028:   VecDestroyVecs(ksp->nwork,&ksp->work);
1029:   VecDestroy(&ksp->vec_rhs);
1030:   VecDestroy(&ksp->vec_sol);
1031:   VecDestroy(&ksp->diagonal);
1032:   VecDestroy(&ksp->truediagonal);

1034:   KSPResetViewers(ksp);

1036:   ksp->setupstage = KSP_SETUP_NEW;
1037:   return(0);
1038: }

1040: /*@
1041:    KSPDestroy - Destroys KSP context.

1043:    Collective on ksp

1045:    Input Parameter:
1046: .  ksp - iterative context obtained from KSPCreate()

1048:    Level: beginner

1050: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1051: @*/
1052: PetscErrorCode  KSPDestroy(KSP *ksp)
1053: {
1055:   PC             pc;

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

1062:   PetscObjectSAWsViewOff((PetscObject)*ksp);

1064:   /*
1065:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1066:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1067:    refcount (and may be shared, e.g., by other ksps).
1068:    */
1069:   pc         = (*ksp)->pc;
1070:   (*ksp)->pc = NULL;
1071:   KSPReset((*ksp));
1072:   (*ksp)->pc = pc;
1073:   if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

1075:   KSPGuessDestroy(&(*ksp)->guess);
1076:   DMDestroy(&(*ksp)->dm);
1077:   PCDestroy(&(*ksp)->pc);
1078:   PetscFree((*ksp)->res_hist_alloc);
1079:   if ((*ksp)->convergeddestroy) {
1080:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
1081:   }
1082:   KSPMonitorCancel((*ksp));
1083:   PetscViewerDestroy(&(*ksp)->eigviewer);
1084:   PetscHeaderDestroy(ksp);
1085:   return(0);
1086: }

1088: /*@
1089:     KSPSetPCSide - Sets the preconditioning side.

1091:     Logically Collective on ksp

1093:     Input Parameter:
1094: .   ksp - iterative context obtained from KSPCreate()

1096:     Output Parameter:
1097: .   side - the preconditioning side, where side is one of
1098: .vb
1099:       PC_LEFT - left preconditioning (default)
1100:       PC_RIGHT - right preconditioning
1101:       PC_SYMMETRIC - symmetric preconditioning
1102: .ve

1104:     Options Database Keys:
1105: .   -ksp_pc_side <right,left,symmetric>

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

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

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

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

1118:     Level: intermediate

1120: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType(), KSP
1121: @*/
1122: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
1123: {
1127:   ksp->pc_side = ksp->pc_side_set = side;
1128:   return(0);
1129: }

1131: /*@
1132:     KSPGetPCSide - Gets the preconditioning side.

1134:     Not Collective

1136:     Input Parameter:
1137: .   ksp - iterative context obtained from KSPCreate()

1139:     Output Parameter:
1140: .   side - the preconditioning side, where side is one of
1141: .vb
1142:       PC_LEFT - left preconditioning (default)
1143:       PC_RIGHT - right preconditioning
1144:       PC_SYMMETRIC - symmetric preconditioning
1145: .ve

1147:     Level: intermediate

1149: .seealso: KSPSetPCSide(), KSP
1150: @*/
1151: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
1152: {

1158:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
1159:   *side = ksp->pc_side;
1160:   return(0);
1161: }

1163: /*@
1164:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1165:    iteration tolerances used by the default KSP convergence tests.

1167:    Not Collective

1169:    Input Parameter:
1170: .  ksp - the Krylov subspace context

1172:    Output Parameters:
1173: +  rtol - the relative convergence tolerance
1174: .  abstol - the absolute convergence tolerance
1175: .  dtol - the divergence tolerance
1176: -  maxits - maximum number of iterations

1178:    Notes:
1179:    The user can specify NULL for any parameter that is not needed.

1181:    Level: intermediate

1183:            maximum, iterations

1185: .seealso: KSPSetTolerances(), KSP
1186: @*/
1187: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1188: {
1191:   if (abstol) *abstol = ksp->abstol;
1192:   if (rtol) *rtol = ksp->rtol;
1193:   if (dtol) *dtol = ksp->divtol;
1194:   if (maxits) *maxits = ksp->max_it;
1195:   return(0);
1196: }

1198: /*@
1199:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1200:    iteration tolerances used by the default KSP convergence testers.

1202:    Logically Collective on ksp

1204:    Input Parameters:
1205: +  ksp - the Krylov subspace context
1206: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1207: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1208: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1209: -  maxits - maximum number of iterations to use

1211:    Options Database Keys:
1212: +  -ksp_atol <abstol> - Sets abstol
1213: .  -ksp_rtol <rtol> - Sets rtol
1214: .  -ksp_divtol <dtol> - Sets dtol
1215: -  -ksp_max_it <maxits> - Sets maxits

1217:    Notes:
1218:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

1223:    Level: intermediate

1225:            convergence, maximum, iterations

1227: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
1228: @*/
1229: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1230: {

1238:   if (rtol != PETSC_DEFAULT) {
1239:     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);
1240:     ksp->rtol = rtol;
1241:   }
1242:   if (abstol != PETSC_DEFAULT) {
1243:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1244:     ksp->abstol = abstol;
1245:   }
1246:   if (dtol != PETSC_DEFAULT) {
1247:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1248:     ksp->divtol = dtol;
1249:   }
1250:   if (maxits != PETSC_DEFAULT) {
1251:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1252:     ksp->max_it = maxits;
1253:   }
1254:   return(0);
1255: }

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

1262:    Logically Collective on ksp

1264:    Input Parameters:
1265: +  ksp - iterative context obtained from KSPCreate()
1266: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

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

1271:    Level: beginner

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

1276: .seealso: KSPGetInitialGuessNonzero(), KSPSetGuessType(), KSPGuessType, KSP
1277: @*/
1278: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1279: {
1283:   ksp->guess_zero = (PetscBool) !(int)flg;
1284:   return(0);
1285: }

1287: /*@
1288:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1289:    a zero initial guess.

1291:    Not Collective

1293:    Input Parameter:
1294: .  ksp - iterative context obtained from KSPCreate()

1296:    Output Parameter:
1297: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1299:    Level: intermediate

1301: .seealso: KSPSetInitialGuessNonzero(), KSP
1302: @*/
1303: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1304: {
1308:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1309:   else *flag = PETSC_TRUE;
1310:   return(0);
1311: }

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

1316:    Logically Collective on ksp

1318:    Input Parameters:
1319: +  ksp - iterative context obtained from KSPCreate()
1320: -  flg - PETSC_TRUE indicates you want the error generated

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

1325:    Level: intermediate

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


1332: .seealso: KSPGetErrorIfNotConverged(), KSP
1333: @*/
1334: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1335: {
1339:   ksp->errorifnotconverged = flg;
1340:   return(0);
1341: }

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

1346:    Not Collective

1348:    Input Parameter:
1349: .  ksp - iterative context obtained from KSPCreate()

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

1354:    Level: intermediate

1356: .seealso: KSPSetErrorIfNotConverged(), KSP
1357: @*/
1358: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1359: {
1363:   *flag = ksp->errorifnotconverged;
1364:   return(0);
1365: }

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

1370:    Logically Collective on ksp

1372:    Input Parameters:
1373: +  ksp - iterative context obtained from KSPCreate()
1374: -  flg - PETSC_TRUE or PETSC_FALSE

1376:    Level: advanced

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

1380: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1381: @*/
1382: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1383: {
1387:   ksp->guess_knoll = flg;
1388:   return(0);
1389: }

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

1395:    Not Collective

1397:    Input Parameter:
1398: .  ksp - iterative context obtained from KSPCreate()

1400:    Output Parameter:
1401: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1403:    Level: advanced

1405: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1406: @*/
1407: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1408: {
1412:   *flag = ksp->guess_knoll;
1413:   return(0);
1414: }

1416: /*@
1417:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1418:    values will be calculated via a Lanczos or Arnoldi process as the linear
1419:    system is solved.

1421:    Not Collective

1423:    Input Parameter:
1424: .  ksp - iterative context obtained from KSPCreate()

1426:    Output Parameter:
1427: .  flg - PETSC_TRUE or PETSC_FALSE

1429:    Options Database Key:
1430: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1432:    Notes:
1433:    Currently this option is not valid for all iterative methods.

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

1439:    Level: advanced

1441: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1442: @*/
1443: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1444: {
1448:   *flg = ksp->calc_sings;
1449:   return(0);
1450: }

1452: /*@
1453:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1454:    values will be calculated via a Lanczos or Arnoldi process as the linear
1455:    system is solved.

1457:    Logically Collective on ksp

1459:    Input Parameters:
1460: +  ksp - iterative context obtained from KSPCreate()
1461: -  flg - PETSC_TRUE or PETSC_FALSE

1463:    Options Database Key:
1464: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

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

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

1473:    Level: advanced

1475: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1476: @*/
1477: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1478: {
1482:   ksp->calc_sings = flg;
1483:   return(0);
1484: }

1486: /*@
1487:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1488:    values will be calculated via a Lanczos or Arnoldi process as the linear
1489:    system is solved.

1491:    Not Collective

1493:    Input Parameter:
1494: .  ksp - iterative context obtained from KSPCreate()

1496:    Output Parameter:
1497: .  flg - PETSC_TRUE or PETSC_FALSE

1499:    Notes:
1500:    Currently this option is not valid for all iterative methods.

1502:    Level: advanced

1504: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1505: @*/
1506: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1507: {
1511:   *flg = ksp->calc_sings;
1512:   return(0);
1513: }

1515: /*@
1516:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1517:    values will be calculated via a Lanczos or Arnoldi process as the linear
1518:    system is solved.

1520:    Logically Collective on ksp

1522:    Input Parameters:
1523: +  ksp - iterative context obtained from KSPCreate()
1524: -  flg - PETSC_TRUE or PETSC_FALSE

1526:    Notes:
1527:    Currently this option is not valid for all iterative methods.

1529:    Level: advanced

1531: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1532: @*/
1533: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1534: {
1538:   ksp->calc_sings = flg;
1539:   return(0);
1540: }

1542: /*@
1543:    KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1544:    will be calculated via a Lanczos or Arnoldi process as the linear
1545:    system is solved.

1547:    Logically Collective on ksp

1549:    Input Parameters:
1550: +  ksp - iterative context obtained from KSPCreate()
1551: -  flg - PETSC_TRUE or PETSC_FALSE

1553:    Notes:
1554:    Currently this option is only valid for the GMRES method.

1556:    Level: advanced

1558: .seealso: KSPComputeRitz(), KSP
1559: @*/
1560: PetscErrorCode  KSPSetComputeRitz(KSP ksp, PetscBool flg)
1561: {
1565:   ksp->calc_ritz = flg;
1566:   return(0);
1567: }

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

1573:    Not Collective

1575:    Input Parameter:
1576: .  ksp - iterative context obtained from KSPCreate()

1578:    Output Parameter:
1579: .  r - right-hand-side vector

1581:    Level: developer

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: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve(), KSP
1610: @*/
1611: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1612: {
1616:   *v = ksp->vec_sol;
1617:   return(0);
1618: }

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

1624:    Collective on ksp

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

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

1634:    Level: developer

1636: .seealso: KSPGetPC(), KSP
1637: @*/
1638: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1639: {

1646:   PetscObjectReference((PetscObject)pc);
1647:   PCDestroy(&ksp->pc);
1648:   ksp->pc = pc;
1649:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1650:   return(0);
1651: }

1653: /*@
1654:    KSPGetPC - Returns a pointer to the preconditioner context
1655:    set with KSPSetPC().

1657:    Not Collective

1659:    Input Parameters:
1660: .  ksp - iterative context obtained from KSPCreate()

1662:    Output Parameter:
1663: .  pc - preconditioner context

1665:    Level: developer

1667: .seealso: KSPSetPC(), KSP
1668: @*/
1669: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1670: {

1676:   if (!ksp->pc) {
1677:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1678:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1679:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1680:     PetscObjectSetOptions((PetscObject)ksp->pc,((PetscObject)ksp)->options);
1681:   }
1682:   *pc = ksp->pc;
1683:   return(0);
1684: }

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

1689:    Collective on ksp

1691:    Input Parameters:
1692: +  ksp - iterative context obtained from KSPCreate()
1693: .  it - iteration number
1694: -  rnorm - relative norm of the residual

1696:    Notes:
1697:    This routine is called by the KSP implementations.
1698:    It does not typically need to be called by the user.

1700:    Level: developer

1702: .seealso: KSPMonitorSet()
1703: @*/
1704: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1705: {
1706:   PetscInt       i, n = ksp->numbermonitors;

1710:   for (i=0; i<n; i++) {
1711:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1712:   }
1713:   return(0);
1714: }

1716: /*

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

1739: /*@C
1740:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1741:    the residual/error etc.

1743:    Logically Collective on ksp

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

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

1756: +  ksp - iterative context obtained from KSPCreate()
1757: .  it - iteration number
1758: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1759: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

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

1776:    Notes:
1777:    The default is to do nothing.  To print the residual, or preconditioned
1778:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1779:    KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1780:    context.

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

1786:    Fortran Notes:
1787:     Only a single monitor function can be set for each KSP object

1789:    Level: beginner

1791: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel(), KSP
1792: @*/
1793: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1794: {
1795:   PetscInt       i;
1797:   PetscBool      identical;

1801:   for (i=0; i<ksp->numbermonitors;i++) {
1802:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
1803:     if (identical) return(0);
1804:   }
1805:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1806:   ksp->monitor[ksp->numbermonitors]          = monitor;
1807:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
1808:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1809:   return(0);
1810: }

1812: /*@
1813:    KSPMonitorCancel - Clears all monitors for a KSP object.

1815:    Logically Collective on ksp

1817:    Input Parameters:
1818: .  ksp - iterative context obtained from KSPCreate()

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

1825:    Level: intermediate

1827: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet(), KSP
1828: @*/
1829: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1830: {
1832:   PetscInt       i;

1836:   for (i=0; i<ksp->numbermonitors; i++) {
1837:     if (ksp->monitordestroy[i]) {
1838:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1839:     }
1840:   }
1841:   ksp->numbermonitors = 0;
1842:   return(0);
1843: }

1845: /*@C
1846:    KSPGetMonitorContext - Gets the monitoring context, as set by
1847:    KSPMonitorSet() for the FIRST monitor only.

1849:    Not Collective

1851:    Input Parameter:
1852: .  ksp - iterative context obtained from KSPCreate()

1854:    Output Parameter:
1855: .  ctx - monitoring context

1857:    Level: intermediate

1859: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSP
1860: @*/
1861: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1862: {
1865:   *ctx =      (ksp->monitorcontext[0]);
1866:   return(0);
1867: }

1869: /*@
1870:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1871:    If set, this array will contain the residual norms computed at each
1872:    iteration of the solver.

1874:    Not Collective

1876:    Input Parameters:
1877: +  ksp - iterative context obtained from KSPCreate()
1878: .  a   - array to hold history
1879: .  na  - size of a
1880: -  reset - PETSC_TRUE indicates the history counter is reset to zero
1881:            for each new linear solve

1883:    Level: advanced

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

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

1892: .seealso: KSPGetResidualHistory(), KSP

1894: @*/
1895: PetscErrorCode  KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
1896: {


1902:   PetscFree(ksp->res_hist_alloc);
1903:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
1904:     ksp->res_hist     = a;
1905:     ksp->res_hist_max = na;
1906:   } else {
1907:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
1908:     else                                           ksp->res_hist_max = 10000; /* like default ksp->max_it */
1909:     PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);

1911:     ksp->res_hist = ksp->res_hist_alloc;
1912:   }
1913:   ksp->res_hist_len   = 0;
1914:   ksp->res_hist_reset = reset;
1915:   return(0);
1916: }

1918: /*@C
1919:    KSPGetResidualHistory - Gets the array used to hold the residual history
1920:    and the number of residuals it contains.

1922:    Not Collective

1924:    Input Parameter:
1925: .  ksp - iterative context obtained from KSPCreate()

1927:    Output Parameters:
1928: +  a   - pointer to array to hold history (or NULL)
1929: -  na  - number of used entries in a (or NULL)

1931:    Level: advanced

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

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

1942: .seealso: KSPGetResidualHistory(), KSP

1944: @*/
1945: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
1946: {
1949:   if (a) *a = ksp->res_hist;
1950:   if (na) *na = ksp->res_hist_len;
1951:   return(0);
1952: }

1954: /*@C
1955:    KSPSetConvergenceTest - Sets the function to be used to determine
1956:    convergence.

1958:    Logically Collective on ksp

1960:    Input Parameters:
1961: +  ksp - iterative context obtained from KSPCreate()
1962: .  converge - pointer to the function
1963: .  cctx    - context for private data for the convergence routine (may be null)
1964: -  destroy - a routine for destroying the context (may be null)

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

1969: +  ksp - iterative context obtained from KSPCreate()
1970: .  it - iteration number
1971: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1972: .  reason - the reason why it has converged or diverged
1973: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()


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

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

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

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

1990:    Level: advanced

1992: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPGetConvergenceTest(), KSPGetAndClearConvergenceTest()
1993: @*/
1994: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
1995: {

2000:   if (ksp->convergeddestroy) {
2001:     (*ksp->convergeddestroy)(ksp->cnvP);
2002:   }
2003:   ksp->converged        = converge;
2004:   ksp->convergeddestroy = destroy;
2005:   ksp->cnvP             = (void*)cctx;
2006:   return(0);
2007: }

2009: /*@C
2010:    KSPGetConvergenceTest - Gets the function to be used to determine
2011:    convergence.

2013:    Logically Collective on ksp

2015:    Input Parameter:
2016: .   ksp - iterative context obtained from KSPCreate()

2018:    Output Parameter:
2019: +  converge - pointer to convergence test function
2020: .  cctx    - context for private data for the convergence routine (may be null)
2021: -  destroy - a routine for destroying the context (may be null)

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

2026: +  ksp - iterative context obtained from KSPCreate()
2027: .  it - iteration number
2028: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2029: .  reason - the reason why it has converged or diverged
2030: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2032:    Level: advanced

2034: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetAndClearConvergenceTest()
2035: @*/
2036: PetscErrorCode  KSPGetConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2037: {
2040:   if (converge) *converge = ksp->converged;
2041:   if (destroy)  *destroy  = ksp->convergeddestroy;
2042:   if (cctx)     *cctx     = ksp->cnvP;
2043:   return(0);
2044: }

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

2049:    Logically Collective on ksp

2051:    Input Parameter:
2052: .   ksp - iterative context obtained from KSPCreate()

2054:    Output Parameter:
2055: +  converge - pointer to convergence test function
2056: .  cctx    - context for private data for the convergence routine
2057: -  destroy - a routine for destroying the context

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

2062: +  ksp - iterative context obtained from KSPCreate()
2063: .  it - iteration number
2064: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2065: .  reason - the reason why it has converged or diverged
2066: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2068:    Level: advanced

2070:    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
2071:           KSPSetConvergenceTest() on this original KSP. If you just called KSPGetConvergenceTest() followed by KSPSetConvergenceTest() the original context information
2072:           would be destroyed and hence the transferred context would be invalid and trigger a crash on use

2074: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetConvergenceTest()
2075: @*/
2076: PetscErrorCode  KSPGetAndClearConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2077: {
2080:   *converge             = ksp->converged;
2081:   *destroy              = ksp->convergeddestroy;
2082:   *cctx                 = ksp->cnvP;
2083:   ksp->converged        = NULL;
2084:   ksp->cnvP             = NULL;
2085:   ksp->convergeddestroy = NULL;
2086:   return(0);
2087: }

2089: /*@C
2090:    KSPGetConvergenceContext - Gets the convergence context set with
2091:    KSPSetConvergenceTest().

2093:    Not Collective

2095:    Input Parameter:
2096: .  ksp - iterative context obtained from KSPCreate()

2098:    Output Parameter:
2099: .  ctx - monitoring context

2101:    Level: advanced

2103: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
2104: @*/
2105: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
2106: {
2109:   *ctx = ksp->cnvP;
2110:   return(0);
2111: }

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

2117:    Collective on ksp

2119:    Input Parameter:
2120: .  ctx - iterative context obtained from KSPCreate()

2122:    Output Parameter:
2123:    Provide exactly one of
2124: +  v - location to stash solution.
2125: -  V - the solution is returned in this location. This vector is created
2126:        internally. This vector should NOT be destroyed by the user with
2127:        VecDestroy().

2129:    Notes:
2130:    This routine can be used in one of two ways
2131: .vb
2132:       KSPBuildSolution(ksp,NULL,&V);
2133:    or
2134:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2135: .ve
2136:    In the first case an internal vector is allocated to store the solution
2137:    (the user cannot destroy this vector). In the second case the solution
2138:    is generated in the vector that the user provides. Note that for certain
2139:    methods, such as KSPCG, the second case requires a copy of the solution,
2140:    while in the first case the call is essentially free since it simply
2141:    returns the vector where the solution already is stored. For some methods
2142:    like GMRES this is a reasonably expensive operation and should only be
2143:    used in truly needed.

2145:    Level: advanced

2147: .seealso: KSPGetSolution(), KSPBuildResidual(), KSP
2148: @*/
2149: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2150: {

2155:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2156:   if (!V) V = &v;
2157:   (*ksp->ops->buildsolution)(ksp,v,V);
2158:   return(0);
2159: }

2161: /*@C
2162:    KSPBuildResidual - Builds the residual in a vector provided.

2164:    Collective on ksp

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

2169:    Output Parameters:
2170: +  v - optional location to stash residual.  If v is not provided,
2171:        then a location is generated.
2172: .  t - work vector.  If not provided then one is generated.
2173: -  V - the residual

2175:    Notes:
2176:    Regardless of whether or not v is provided, the residual is
2177:    returned in V.

2179:    Level: advanced

2181: .seealso: KSPBuildSolution()
2182: @*/
2183: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2184: {
2186:   PetscBool      flag = PETSC_FALSE;
2187:   Vec            w    = v,tt = t;

2191:   if (!w) {
2192:     VecDuplicate(ksp->vec_rhs,&w);
2193:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2194:   }
2195:   if (!tt) {
2196:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2197:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2198:   }
2199:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
2200:   if (flag) {VecDestroy(&tt);}
2201:   return(0);
2202: }

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

2208:    Logically Collective on ksp

2210:    Input Parameter:
2211: +  ksp - the KSP context
2212: -  scale - PETSC_TRUE or PETSC_FALSE

2214:    Options Database Key:
2215: +   -ksp_diagonal_scale -
2216: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


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

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

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

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

2234:    Level: intermediate

2236: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2237: @*/
2238: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2239: {
2243:   ksp->dscale = scale;
2244:   return(0);
2245: }

2247: /*@
2248:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2249:                           right hand side

2251:    Not Collective

2253:    Input Parameter:
2254: .  ksp - the KSP context

2256:    Output Parameter:
2257: .  scale - PETSC_TRUE or PETSC_FALSE

2259:    Notes:
2260:     BE CAREFUL with this routine: it actually scales the matrix and right
2261:     hand side that define the system. After the system is solved the matrix
2262:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

2264:    Level: intermediate

2266: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2267: @*/
2268: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
2269: {
2273:   *scale = ksp->dscale;
2274:   return(0);
2275: }

2277: /*@
2278:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2279:      back after solving.

2281:    Logically Collective on ksp

2283:    Input Parameter:
2284: +  ksp - the KSP context
2285: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2286:          rescale (default)

2288:    Notes:
2289:      Must be called after KSPSetDiagonalScale()

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

2296:    Level: intermediate

2298: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix(), KSP
2299: @*/
2300: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2301: {
2305:   ksp->dscalefix = fix;
2306:   return(0);
2307: }

2309: /*@
2310:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2311:      back after solving.

2313:    Not Collective

2315:    Input Parameter:
2316: .  ksp - the KSP context

2318:    Output Parameter:
2319: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2320:          rescale (default)

2322:    Notes:
2323:      Must be called after KSPSetDiagonalScale()

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

2330:    Level: intermediate

2332: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2333: @*/
2334: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2335: {
2339:   *fix = ksp->dscalefix;
2340:   return(0);
2341: }

2343: /*@C
2344:    KSPSetComputeOperators - set routine to compute the linear operators

2346:    Logically Collective

2348:    Input Arguments:
2349: +  ksp - the KSP context
2350: .  func - function to compute the operators
2351: -  ctx - optional context

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

2356: +  ksp - the KSP context
2357: .  A - the linear operator
2358: .  B - preconditioning matrix
2359: -  ctx - optional user-provided context

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

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

2367:    Level: beginner

2369: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2370: @*/
2371: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2372: {
2374:   DM             dm;

2378:   KSPGetDM(ksp,&dm);
2379:   DMKSPSetComputeOperators(dm,func,ctx);
2380:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2381:   return(0);
2382: }

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

2387:    Logically Collective

2389:    Input Arguments:
2390: +  ksp - the KSP context
2391: .  func - function to compute the right hand side
2392: -  ctx - optional context

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

2397: +  ksp - the KSP context
2398: .  b - right hand side of linear system
2399: -  ctx - optional user-provided context

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

2404:    Level: beginner

2406: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2407: @*/
2408: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2409: {
2411:   DM             dm;

2415:   KSPGetDM(ksp,&dm);
2416:   DMKSPSetComputeRHS(dm,func,ctx);
2417:   return(0);
2418: }

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

2423:    Logically Collective

2425:    Input Arguments:
2426: +  ksp - the KSP context
2427: .  func - function to compute the initial guess
2428: -  ctx - optional context

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

2433: +  ksp - the KSP context
2434: .  x - solution vector
2435: -  ctx - optional user-provided context

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

2440:    Level: beginner

2442: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2443: @*/
2444: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2445: {
2447:   DM             dm;

2451:   KSPGetDM(ksp,&dm);
2452:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2453:   return(0);
2454: }