Actual source code: itfunc.c

petsc-master 2020-08-04
Report Typos and Errors
  1: /*
  2:       Interface KSP routines that the user calls.
  3: */

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

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

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

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

 23:    Not Collective

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

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

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

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

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

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

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

 49:    Level: advanced

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

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

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

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

 76:    Not Collective

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

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

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

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

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

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

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

110:    Level: advanced

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

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

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

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

139:    Not Collective

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

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

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

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

168:    Level: advanced

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

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

187:    Collective on ksp

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

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

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

201:    Level: advanced

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

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

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

225:    Collective on ksp

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

231:    Level: intermediate

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

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

247: /*@
248:    KSPGetReusePreconditioner - Determines if the KSP reuses the current preconditioner even if the operator in the preconditioner has changed.

250:    Collective on ksp

252:    Input Parameters:
253: .  ksp   - iterative context obtained from KSPCreate()

255:    Output Parameters:
256: .  flag - the boolean flag

258:    Level: intermediate

260: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), KSPSetReusePreconditioner(), KSP
261: @*/
262: PetscErrorCode  KSPGetReusePreconditioner(KSP ksp,PetscBool *flag)
263: {

269:   *flag = PETSC_FALSE;
270:   if (ksp->pc) {
271:     PCGetReusePreconditioner(ksp->pc,flag);
272:   }
273:   return(0);
274: }

276: /*@
277:    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

279:    Collective on ksp

281:    Input Parameters:
282: +  ksp   - iterative context obtained from KSPCreate()
283: -  flag - PETSC_TRUE to skip calling the PCSetFromOptions()

285:    Level: intermediate

287: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), PCSetReusePreconditioner(), KSP
288: @*/
289: PetscErrorCode  KSPSetSkipPCSetFromOptions(KSP ksp,PetscBool flag)
290: {
293:   ksp->skippcsetfromoptions = flag;
294:   return(0);
295: }

297: /*@
298:    KSPSetUp - Sets up the internal data structures for the
299:    later use of an iterative solver.

301:    Collective on ksp

303:    Input Parameter:
304: .  ksp   - iterative context obtained from KSPCreate()

306:    Level: developer

308: .seealso: KSPCreate(), KSPSolve(), KSPDestroy(), KSP
309: @*/
310: PetscErrorCode KSPSetUp(KSP ksp)
311: {
313:   Mat            A,B;
314:   Mat            mat,pmat;
315:   MatNullSpace   nullsp;
316:   PCFailedReason pcreason;


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

324:   if (!((PetscObject)ksp)->type_name) {
325:     KSPSetType(ksp,KSPGMRES);
326:   }
327:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);

329:   if (ksp->dmActive && !ksp->setupstage) {
330:     /* first time in so build matrix and vector data structures using DM */
331:     if (!ksp->vec_rhs) {DMCreateGlobalVector(ksp->dm,&ksp->vec_rhs);}
332:     if (!ksp->vec_sol) {DMCreateGlobalVector(ksp->dm,&ksp->vec_sol);}
333:     DMCreateMatrix(ksp->dm,&A);
334:     KSPSetOperators(ksp,A,A);
335:     PetscObjectDereference((PetscObject)A);
336:   }

338:   if (ksp->dmActive) {
339:     DMKSP kdm;
340:     DMGetDMKSP(ksp->dm,&kdm);

342:     if (kdm->ops->computeinitialguess && ksp->setupstage != KSP_SETUP_NEWRHS) {
343:       /* only computes initial guess the first time through */
344:       (*kdm->ops->computeinitialguess)(ksp,ksp->vec_sol,kdm->initialguessctx);
345:       KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);
346:     }
347:     if (kdm->ops->computerhs) {
348:       (*kdm->ops->computerhs)(ksp,ksp->vec_rhs,kdm->rhsctx);
349:     }

351:     if (ksp->setupstage != KSP_SETUP_NEWRHS) {
352:       if (kdm->ops->computeoperators) {
353:         KSPGetOperators(ksp,&A,&B);
354:         (*kdm->ops->computeoperators)(ksp,A,B,kdm->operatorsctx);
355:       } else SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONGSTATE,"You called KSPSetDM() but did not use DMKSPSetComputeOperators() or KSPSetDMActive(ksp,PETSC_FALSE);");
356:     }
357:   }

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

362:   switch (ksp->setupstage) {
363:   case KSP_SETUP_NEW:
364:     (*ksp->ops->setup)(ksp);
365:     break;
366:   case KSP_SETUP_NEWMATRIX: {   /* This should be replaced with a more general mechanism */
367:     if (ksp->setupnewmatrix) {
368:       (*ksp->ops->setup)(ksp);
369:     }
370:   } break;
371:   default: break;
372:   }

374:   if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
375:   PCGetOperators(ksp->pc,&mat,&pmat);
376:   /* scale the matrix if requested */
377:   if (ksp->dscale) {
378:     PetscScalar *xx;
379:     PetscInt    i,n;
380:     PetscBool   zeroflag = PETSC_FALSE;
381:     if (!ksp->pc) {KSPGetPC(ksp,&ksp->pc);}
382:     if (!ksp->diagonal) { /* allocate vector to hold diagonal */
383:       MatCreateVecs(pmat,&ksp->diagonal,NULL);
384:     }
385:     MatGetDiagonal(pmat,ksp->diagonal);
386:     VecGetLocalSize(ksp->diagonal,&n);
387:     VecGetArray(ksp->diagonal,&xx);
388:     for (i=0; i<n; i++) {
389:       if (xx[i] != 0.0) xx[i] = 1.0/PetscSqrtReal(PetscAbsScalar(xx[i]));
390:       else {
391:         xx[i]    = 1.0;
392:         zeroflag = PETSC_TRUE;
393:       }
394:     }
395:     VecRestoreArray(ksp->diagonal,&xx);
396:     if (zeroflag) {
397:       PetscInfo(ksp,"Zero detected in diagonal of matrix, using 1 at those locations\n");
398:     }
399:     MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
400:     if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
401:     ksp->dscalefix2 = PETSC_FALSE;
402:   }
403:   PetscLogEventEnd(KSP_SetUp,ksp,ksp->vec_rhs,ksp->vec_sol,0);
404:   PCSetErrorIfFailure(ksp->pc,ksp->errorifnotconverged);
405:   PCSetUp(ksp->pc);
406:   PCGetFailedReason(ksp->pc,&pcreason);
407:   if (pcreason) {
408:     ksp->reason = KSP_DIVERGED_PC_FAILED;
409:   }

411:   MatGetNullSpace(mat,&nullsp);
412:   if (nullsp) {
413:     PetscBool test = PETSC_FALSE;
414:     PetscOptionsGetBool(((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_test_null_space",&test,NULL);
415:     if (test) {
416:       MatNullSpaceTest(nullsp,mat,NULL);
417:     }
418:   }
419:   ksp->setupstage = KSP_SETUP_NEWRHS;
420:   return(0);
421: }

423: static PetscErrorCode KSPReasonView_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
424: {
426:   PetscBool      isAscii;

429:   if (format != PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,format);}
430:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isAscii);
431:   if (isAscii) {
432:     PetscViewerASCIIAddTab(viewer,((PetscObject)ksp)->tablevel);
433:     if (ksp->reason > 0) {
434:       if (((PetscObject) ksp)->prefix) {
435:         PetscViewerASCIIPrintf(viewer,"Linear %s solve converged due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
436:       } else {
437:         PetscViewerASCIIPrintf(viewer,"Linear solve converged due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
438:       }
439:     } else {
440:       if (((PetscObject) ksp)->prefix) {
441:         PetscViewerASCIIPrintf(viewer,"Linear %s solve did not converge due to %s iterations %D\n",((PetscObject) ksp)->prefix,KSPConvergedReasons[ksp->reason],ksp->its);
442:       } else {
443:         PetscViewerASCIIPrintf(viewer,"Linear solve did not converge due to %s iterations %D\n",KSPConvergedReasons[ksp->reason],ksp->its);
444:       }
445:       if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
446:         PCFailedReason reason;
447:         PCGetFailedReason(ksp->pc,&reason);
448:         PetscViewerASCIIPrintf(viewer,"               PC_FAILED due to %s \n",PCFailedReasons[reason]);
449:       }
450:     }
451:     PetscViewerASCIISubtractTab(viewer,((PetscObject)ksp)->tablevel);
452:   }
453:   if (format != PETSC_VIEWER_DEFAULT) {PetscViewerPopFormat(viewer);}
454:   return(0);
455: }

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

460:    Collective on ksp

462:    Parameter:
463: +  ksp - iterative context obtained from KSPCreate()
464: -  viewer - the viewer to display the reason


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

470:    Level: beginner

472: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
473:           KSPSolveTranspose(), KSPGetIterationNumber(), KSP
474: @*/
475: PetscErrorCode KSPReasonView(KSP ksp,PetscViewer viewer)
476: {

480:   KSPReasonView_Internal(ksp, viewer, PETSC_VIEWER_DEFAULT);
481:   return(0);
482: }

484: #if defined(PETSC_HAVE_THREADSAFETY)
485: #define KSPReasonViewFromOptions KSPReasonViewFromOptionsUnsafe
486: #else
487: #endif
488: /*@C
489:   KSPReasonViewFromOptions - Processes command line options to determine if/how a KSPReason is to be viewed.

491:   Collective on ksp

493:   Input Parameters:
494: . ksp   - the KSP object

496:   Level: intermediate

498: @*/
499: PetscErrorCode KSPReasonViewFromOptions(KSP ksp)
500: {
501:   PetscViewer       viewer;
502:   PetscBool         flg;
503:   PetscViewerFormat format;
504:   PetscErrorCode    ierr;

507:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)ksp),((PetscObject)ksp)->options,((PetscObject)ksp)->prefix,"-ksp_converged_reason",&viewer,&format,&flg);
508:   if (flg) {
509:     KSPReasonView_Internal(ksp, viewer, format);
510:     PetscViewerDestroy(&viewer);
511:   }
512:   return(0);
513: }

515:  #include <petscdraw.h>

517: static PetscErrorCode KSPViewEigenvalues_Internal(KSP ksp, PetscBool isExplicit, PetscViewer viewer, PetscViewerFormat format)
518: {
519:   PetscReal     *r, *c;
520:   PetscInt       n, i, neig;
521:   PetscBool      isascii, isdraw;
522:   PetscMPIInt    rank;

526:   MPI_Comm_rank(PetscObjectComm((PetscObject) ksp), &rank);
527:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
528:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW,  &isdraw);
529:   if (isExplicit) {
530:     VecGetSize(ksp->vec_sol,&n);
531:     PetscMalloc2(n, &r, n, &c);
532:     KSPComputeEigenvaluesExplicitly(ksp, n, r, c);
533:     neig = n;
534:   } else {
535:     PetscInt nits;

537:     KSPGetIterationNumber(ksp, &nits);
538:     n    = nits+2;
539:     if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any eigenvalues\n");return(0);}
540:     PetscMalloc2(n, &r, n, &c);
541:     KSPComputeEigenvalues(ksp, n, r, c, &neig);
542:   }
543:   if (isascii) {
544:     PetscViewerASCIIPrintf(viewer, "%s computed eigenvalues\n", isExplicit ? "Explicitly" : "Iteratively");
545:     for (i = 0; i < neig; ++i) {
546:       if (c[i] >= 0.0) {PetscViewerASCIIPrintf(viewer, "%g + %gi\n", (double) r[i],  (double) c[i]);}
547:       else             {PetscViewerASCIIPrintf(viewer, "%g - %gi\n", (double) r[i], -(double) c[i]);}
548:     }
549:   } else if (isdraw && !rank) {
550:     PetscDraw   draw;
551:     PetscDrawSP drawsp;

553:     if (format == PETSC_VIEWER_DRAW_CONTOUR) {
554:       KSPPlotEigenContours_Private(ksp,neig,r,c);
555:     } else {
556:       if (!ksp->eigviewer) {PetscViewerDrawOpen(PETSC_COMM_SELF,NULL,isExplicit ? "Explicitly Computed Eigenvalues" : "Iteratively Computed Eigenvalues",PETSC_DECIDE,PETSC_DECIDE,400,400,&ksp->eigviewer);}
557:       PetscViewerDrawGetDraw(ksp->eigviewer,0,&draw);
558:       PetscDrawSPCreate(draw,1,&drawsp);
559:       PetscDrawSPReset(drawsp);
560:       for (i = 0; i < neig; ++i) {PetscDrawSPAddPoint(drawsp,r+i,c+i);}
561:       PetscDrawSPDraw(drawsp,PETSC_TRUE);
562:       PetscDrawSPSave(drawsp);
563:       PetscDrawSPDestroy(&drawsp);
564:     }
565:   }
566:   PetscFree2(r, c);
567:   return(0);
568: }

570: static PetscErrorCode KSPViewSingularvalues_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
571: {
572:   PetscReal      smax, smin;
573:   PetscInt       nits;
574:   PetscBool      isascii;

578:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
579:   KSPGetIterationNumber(ksp, &nits);
580:   if (!nits) {PetscViewerASCIIPrintf(viewer, "Zero iterations in solver, cannot approximate any singular values\n");return(0);}
581:   KSPComputeExtremeSingularValues(ksp, &smax, &smin);
582:   if (isascii) {PetscViewerASCIIPrintf(viewer, "Iteratively computed extreme singular values: max %g min %g max/min %g\n",(double)smax,(double)smin,(double)(smax/smin));}
583:   return(0);
584: }

586: static PetscErrorCode KSPViewFinalResidual_Internal(KSP ksp, PetscViewer viewer, PetscViewerFormat format)
587: {
588:   PetscBool      isascii;

592:   PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
593:   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");
594:   if (isascii) {
595:     Mat       A;
596:     Vec       t;
597:     PetscReal norm;

599:     PCGetOperators(ksp->pc, &A, NULL);
600:     VecDuplicate(ksp->vec_rhs, &t);
601:     KSP_MatMult(ksp, A, ksp->vec_sol, t);
602:     VecAYPX(t, -1.0, ksp->vec_rhs);
603:     VecNorm(t, NORM_2, &norm);
604:     VecDestroy(&t);
605:     PetscViewerASCIIPrintf(viewer, "KSP final norm of residual %g\n", (double) norm);
606:   }
607:   return(0);
608: }

610: static PetscErrorCode KSPSolve_Private(KSP ksp,Vec b,Vec x)
611: {
613:   PetscBool      flg = PETSC_FALSE,inXisinB=PETSC_FALSE,guess_zero;
614:   Mat            mat,pmat;
615:   MPI_Comm       comm;
616:   MatNullSpace   nullsp;
617:   Vec            btmp,vec_rhs=NULL;

620:   comm = PetscObjectComm((PetscObject)ksp);
621:   if (x && x == b) {
622:     if (!ksp->guess_zero) SETERRQ(comm,PETSC_ERR_ARG_INCOMP,"Cannot use x == b with nonzero initial guess");
623:     VecDuplicate(b,&x);
624:     inXisinB = PETSC_TRUE;
625:   }
626:   if (b) {
627:     PetscObjectReference((PetscObject)b);
628:     VecDestroy(&ksp->vec_rhs);
629:     ksp->vec_rhs = b;
630:   }
631:   if (x) {
632:     PetscObjectReference((PetscObject)x);
633:     VecDestroy(&ksp->vec_sol);
634:     ksp->vec_sol = x;
635:   }

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

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

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

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

646:   if (ksp->guess) {
647:     PetscObjectState ostate,state;

649:     KSPGuessSetUp(ksp->guess);
650:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&ostate);
651:     KSPGuessFormGuess(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
652:     PetscObjectStateGet((PetscObject)ksp->vec_sol,&state);
653:     if (state != ostate) {
654:       ksp->guess_zero = PETSC_FALSE;
655:     } else {
656:       PetscInfo(ksp,"Using zero initial guess since the KSPGuess object did not change the vector\n");
657:       ksp->guess_zero = PETSC_TRUE;
658:     }
659:   }

661:   /* KSPSetUp() scales the matrix if needed */
662:   KSPSetUp(ksp);
663:   KSPSetUpOnBlocks(ksp);

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

667:   PCGetOperators(ksp->pc,&mat,&pmat);
668:   /* diagonal scale RHS if called for */
669:   if (ksp->dscale) {
670:     VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
671:     /* second time in, but matrix was scaled back to original */
672:     if (ksp->dscalefix && ksp->dscalefix2) {
673:       Mat mat,pmat;

675:       PCGetOperators(ksp->pc,&mat,&pmat);
676:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
677:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
678:     }

680:     /* scale initial guess */
681:     if (!ksp->guess_zero) {
682:       if (!ksp->truediagonal) {
683:         VecDuplicate(ksp->diagonal,&ksp->truediagonal);
684:         VecCopy(ksp->diagonal,ksp->truediagonal);
685:         VecReciprocal(ksp->truediagonal);
686:       }
687:       VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->truediagonal);
688:     }
689:   }
690:   PCPreSolve(ksp->pc,ksp);

692:   if (ksp->guess_zero) { VecSet(ksp->vec_sol,0.0);}
693:   if (ksp->guess_knoll) { /* The Knoll trick is independent on the KSPGuess specified */
694:     PCApply(ksp->pc,ksp->vec_rhs,ksp->vec_sol);
695:     KSP_RemoveNullSpace(ksp,ksp->vec_sol);
696:     ksp->guess_zero = PETSC_FALSE;
697:   }

699:   /* can we mark the initial guess as zero for this solve? */
700:   guess_zero = ksp->guess_zero;
701:   if (!ksp->guess_zero) {
702:     PetscReal norm;

704:     VecNormAvailable(ksp->vec_sol,NORM_2,&flg,&norm);
705:     if (flg && !norm) ksp->guess_zero = PETSC_TRUE;
706:   }
707:   if (ksp->transpose_solve) {
708:     MatGetNullSpace(pmat,&nullsp);
709:   } else {
710:     MatGetTransposeNullSpace(pmat,&nullsp);
711:   }
712:   if (nullsp) {
713:     VecDuplicate(ksp->vec_rhs,&btmp);
714:     VecCopy(ksp->vec_rhs,btmp);
715:     MatNullSpaceRemove(nullsp,btmp);
716:     vec_rhs      = ksp->vec_rhs;
717:     ksp->vec_rhs = btmp;
718:   }
719:   VecLockReadPush(ksp->vec_rhs);
720:   if (ksp->reason == KSP_DIVERGED_PC_FAILED) {
721:     VecSetInf(ksp->vec_sol);
722:   }
723:   (*ksp->ops->solve)(ksp);

725:   VecLockReadPop(ksp->vec_rhs);
726:   if (nullsp) {
727:     ksp->vec_rhs = vec_rhs;
728:     VecDestroy(&btmp);
729:   }

731:   ksp->guess_zero = guess_zero;

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

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

739:   /* diagonal scale solution if called for */
740:   if (ksp->dscale) {
741:     VecPointwiseMult(ksp->vec_sol,ksp->vec_sol,ksp->diagonal);
742:     /* unscale right hand side and matrix */
743:     if (ksp->dscalefix) {
744:       Mat mat,pmat;

746:       VecReciprocal(ksp->diagonal);
747:       VecPointwiseMult(ksp->vec_rhs,ksp->vec_rhs,ksp->diagonal);
748:       PCGetOperators(ksp->pc,&mat,&pmat);
749:       MatDiagonalScale(pmat,ksp->diagonal,ksp->diagonal);
750:       if (mat != pmat) {MatDiagonalScale(mat,ksp->diagonal,ksp->diagonal);}
751:       VecReciprocal(ksp->diagonal);
752:       ksp->dscalefix2 = PETSC_TRUE;
753:     }
754:   }
755:   PetscLogEventEnd(KSP_Solve,ksp,ksp->vec_rhs,ksp->vec_sol,0);
756:   if (ksp->guess) {
757:     KSPGuessUpdate(ksp->guess,ksp->vec_rhs,ksp->vec_sol);
758:   }
759:   if (ksp->postsolve) {
760:     (*ksp->postsolve)(ksp,ksp->vec_rhs,ksp->vec_sol,ksp->postctx);
761:   }

763:   PCGetOperators(ksp->pc,&mat,&pmat);
764:   if (ksp->viewEV)       {KSPViewEigenvalues_Internal(ksp, PETSC_FALSE, ksp->viewerEV,    ksp->formatEV);}
765:   if (ksp->viewEVExp)    {KSPViewEigenvalues_Internal(ksp, PETSC_TRUE,  ksp->viewerEVExp, ksp->formatEVExp);}
766:   if (ksp->viewSV)       {KSPViewSingularvalues_Internal(ksp, ksp->viewerSV, ksp->formatSV);}
767:   if (ksp->viewFinalRes) {KSPViewFinalResidual_Internal(ksp, ksp->viewerFinalRes, ksp->formatFinalRes);}
768:   if (ksp->viewMat)      {ObjectView((PetscObject) mat,           ksp->viewerMat,    ksp->formatMat);}
769:   if (ksp->viewPMat)     {ObjectView((PetscObject) pmat,          ksp->viewerPMat,   ksp->formatPMat);}
770:   if (ksp->viewRhs)      {ObjectView((PetscObject) ksp->vec_rhs,  ksp->viewerRhs,    ksp->formatRhs);}
771:   if (ksp->viewSol)      {ObjectView((PetscObject) ksp->vec_sol,  ksp->viewerSol,    ksp->formatSol);}
772:   if (ksp->view)         {ObjectView((PetscObject) ksp,           ksp->viewer,       ksp->format);}
773:   if (ksp->viewDScale)   {ObjectView((PetscObject) ksp->diagonal, ksp->viewerDScale, ksp->formatDScale);}
774:   if (ksp->viewMatExp)   {
775:     Mat A, B;

777:     PCGetOperators(ksp->pc, &A, NULL);
778:     if (ksp->transpose_solve) {
779:       Mat AT;

781:       MatCreateTranspose(A, &AT);
782:       MatComputeOperator(AT, MATAIJ, &B);
783:       MatDestroy(&AT);
784:     } else {
785:       MatComputeOperator(A, MATAIJ, &B);
786:     }
787:     ObjectView((PetscObject) B, ksp->viewerMatExp, ksp->formatMatExp);
788:     MatDestroy(&B);
789:   }
790:   if (ksp->viewPOpExp)   {
791:     Mat B;

793:     KSPComputeOperator(ksp, MATAIJ, &B);
794:     ObjectView((PetscObject) B, ksp->viewerPOpExp, ksp->formatPOpExp);
795:     MatDestroy(&B);
796:   }

798:   if (inXisinB) {
799:     VecCopy(x,b);
800:     VecDestroy(&x);
801:   }
802:   PetscObjectSAWsBlock((PetscObject)ksp);
803:   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]);
804:   return(0);
805: }

807: /*@
808:    KSPSolve - Solves linear system.

810:    Collective on ksp

812:    Parameters:
813: +  ksp - iterative context obtained from KSPCreate()
814: .  b - the right hand side vector
815: -  x - the solution (this may be the same vector as b, then b will be overwritten with answer)

817:    Options Database Keys:
818: +  -ksp_view_eigenvalues - compute preconditioned operators eigenvalues
819: .  -ksp_view_eigenvalues_explicitly - compute the eigenvalues by forming the dense operator and using LAPACK
820: .  -ksp_view_mat binary - save matrix to the default binary viewer
821: .  -ksp_view_pmat binary - save matrix used to build preconditioner to the default binary viewer
822: .  -ksp_view_rhs binary - save right hand side vector to the default binary viewer
823: .  -ksp_view_solution binary - save computed solution vector to the default binary viewer
824:            (can be read later with src/ksp/tutorials/ex10.c for testing solvers)
825: .  -ksp_view_mat_explicit - for matrix-free operators, computes the matrix entries and views them
826: .  -ksp_view_preconditioned_operator_explicit - computes the product of the preconditioner and matrix as an explicit matrix and views it
827: .  -ksp_converged_reason - print reason for converged or diverged, also prints number of iterations
828: .  -ksp_view_final_residual - print 2-norm of true linear system residual at the end of the solution process
829: -  -ksp_view - print the ksp data structure at the end of the system solution

831:    Notes:

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

835:    The operator is specified with KSPSetOperators().

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

840:    If you provide a matrix that has a MatSetNullSpace() and MatSetTransposeNullSpace() this will use that information to solve singular systems
841:    in the least squares sense with a norm minimizing solution.
842: $
843: $                   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()
844: $
845: $    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
846: $    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
847: $    direction thus the solution which is a linear combination of the search directions has no component in the nullspace(A).
848: $
849: $    We recommend always using GMRES for such singular systems.
850: $    If nullspace(A) = nullspace(A') (note symmetric matrices always satisfy this property) then both left and right preconditioning will work
851: $    If nullspace(A) != nullspace(A') then left preconditioning will work but right preconditioning may not work (or it may).

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


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

863:    Understanding Convergence:
864:    The routines KSPMonitorSet(), KSPComputeEigenvalues(), and
865:    KSPComputeEigenvaluesExplicitly() provide information on additional
866:    options to monitor convergence and print eigenvalue information.

868:    Level: beginner

870: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
871:           KSPSolveTranspose(), KSPGetIterationNumber(), MatNullSpaceCreate(), MatSetNullSpace(), MatSetTransposeNullSpace(), KSP
872: @*/
873: PetscErrorCode KSPSolve(KSP ksp,Vec b,Vec x)
874: {

881:   ksp->transpose_solve = PETSC_FALSE;
882:   KSPSolve_Private(ksp,b,x);
883:   return(0);
884: }

886: /*@
887:    KSPSolveTranspose - Solves the transpose of a linear system.

889:    Collective on ksp

891:    Input Parameters:
892: +  ksp - iterative context obtained from KSPCreate()
893: .  b - right hand side vector
894: -  x - solution vector

896:    Notes:
897:     For complex numbers this solve the non-Hermitian transpose system.

899:    Developer Notes:
900:     We need to implement a KSPSolveHermitianTranspose()

902:    Level: developer

904: .seealso: KSPCreate(), KSPSetUp(), KSPDestroy(), KSPSetTolerances(), KSPConvergedDefault(),
905:           KSPSolve(), KSP
906: @*/
907: PetscErrorCode KSPSolveTranspose(KSP ksp,Vec b,Vec x)
908: {

915:   ksp->transpose_solve = PETSC_TRUE;
916:   KSPSolve_Private(ksp,b,x);
917:   return(0);
918: }

920: static PetscErrorCode KSPViewFinalMatResidual_Internal(KSP ksp, Mat B, Mat X, PetscViewer viewer, PetscViewerFormat format, PetscInt shift)
921: {
922:   Mat            A, R;
923:   PetscReal      *norms;
924:   PetscInt       i, N;
925:   PetscBool      flg;

929:   PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &flg);
930:   if (flg) {
931:     PCGetOperators(ksp->pc, &A, NULL);
932:     MatMatMult(A, X, MAT_INITIAL_MATRIX, PETSC_DEFAULT, &R);
933:     MatAYPX(R, -1.0, B, SAME_NONZERO_PATTERN);
934:     MatGetSize(R, NULL, &N);
935:     PetscMalloc1(N, &norms);
936:     MatGetColumnNorms(R, NORM_2, norms);
937:     MatDestroy(&R);
938:     for (i = 0; i < N; ++i) {
939:       PetscViewerASCIIPrintf(viewer, "%s #%D %g\n", i == 0 ? "KSP final norm of residual" : "                          ", shift + i, (double)norms[i]);
940:     }
941:     PetscFree(norms);
942:   }
943:   return(0);
944: }

946: /*@
947:      KSPMatSolve - Solves a linear system with multiple right-hand sides stored as a MATDENSE. Unlike KSPSolve(), B and X must be different matrices.

949:    Input Parameters:
950: +     ksp - iterative context
951: -     B - block of right-hand sides

953:    Output Parameter:
954: .     X - block of solutions

956:    Notes:
957:      This is a stripped-down version of KSPSolve(), which only handles -ksp_view, -ksp_converged_reason, and -ksp_view_final_residual.

959:    Level: intermediate

961: .seealso:  KSPSolve(), MatMatSolve(), MATDENSE, KSPHPDDM, PCBJACOBI, PCASM
962: @*/
963: PetscErrorCode KSPMatSolve(KSP ksp, Mat B, Mat X)
964: {
965:   Mat            A, vB, vX;
966:   Vec            cb, cx;
967:   PetscInt       m1, M1, m2, M2, n1, N1, n2, N2, Bbn = PETSC_DECIDE;
968:   PetscBool      match;

977:   if (!B->assembled) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Not for unassembled matrix");
978:   MatCheckPreallocated(X, 3);
979:   if (!X->assembled) {
980:     MatSetOption(X, MAT_NO_OFF_PROC_ENTRIES, PETSC_TRUE);
981:     MatAssemblyBegin(X, MAT_FINAL_ASSEMBLY);
982:     MatAssemblyEnd(X, MAT_FINAL_ASSEMBLY);
983:   }
984:   if (B == X) SETERRQ(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_IDN, "B and X must be different matrices");
985:   KSPGetOperators(ksp, &A, NULL);
986:   MatGetLocalSize(A, &m1, NULL);
987:   MatGetLocalSize(B, &m2, &n2);
988:   MatGetSize(A, &M1, NULL);
989:   MatGetSize(B, &M2, &N2);
990:   if (m1 != m2 || M1 != M2) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Cannot use a block of right-hand sides with (m2,M2) = (%D,%D) for a linear system with (m1,M1) = (%D,%D)", m2, M2, m1, M1);
991:   MatGetLocalSize(X, &m1, &n1);
992:   MatGetSize(X, &M1, &N1);
993:   if (m1 != m2 || M1 != M2 || n1 != n2 || N1 != N2) SETERRQ8(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Incompatible block of right-hand sides (m2,M2)x(n2,N2) = (%D,%D)x(%D,%D) and solutions (m1,M1)x(n1,N1) = (%D,%D)x(%D,%D)", m2, M2, n2, N2, m1, M1, n1, N1);
994:   PetscObjectBaseTypeCompareAny((PetscObject)B, &match, MATSEQDENSE, MATMPIDENSE, "");
995:   if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of right-hand sides not stored in a dense Mat");
996:   PetscObjectBaseTypeCompareAny((PetscObject)X, &match, MATSEQDENSE, MATMPIDENSE, "");
997:   if (!match) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Provided block of solutions not stored in a dense Mat");
998:   KSPSetUp(ksp);
999:   KSPSetUpOnBlocks(ksp);
1000:   if (ksp->ops->matsolve) {
1001:     if (ksp->guess_zero) {
1002:       MatZeroEntries(X);
1003:     }
1004:     PetscLogEventBegin(KSP_MatSolve, ksp, B, X, 0);
1005:     KSPGetMatSolveBlockSize(ksp, &Bbn);
1006:     /* by default, do a single solve with all columns */
1007:     if (Bbn == PETSC_DECIDE) Bbn = N2;
1008:     else if (Bbn < 1)        SETERRQ1(PetscObjectComm((PetscObject)ksp), PETSC_ERR_ARG_OUTOFRANGE, "KSPMatSolve() block size %D must be positive", Bbn);
1009:     PetscInfo2(ksp, "KSP type %s solving using blocks of width at most %D\n", ((PetscObject)ksp)->type_name, Bbn);
1010:     /* if -ksp_matsolve_block_size is greater than the actual number of columns, do a single solve with all columns */
1011:     if (Bbn >= N2) {
1012:       (*ksp->ops->matsolve)(ksp, B, X);
1013:       if (ksp->viewFinalRes) {
1014:         KSPViewFinalMatResidual_Internal(ksp, B, X, ksp->viewerFinalRes, ksp->formatFinalRes, 0);
1015:       }
1016:       if (ksp->viewReason) {
1017:         KSPReasonView(ksp, ksp->viewerReason);
1018:       }
1019:     } else {
1020:       for (n2 = 0; n2 < N2; n2 += Bbn) {
1021:         MatDenseGetSubMatrix(B, n2, PetscMin(n2+Bbn, N2), &vB);
1022:         MatDenseGetSubMatrix(X, n2, PetscMin(n2+Bbn, N2), &vX);
1023:         (*ksp->ops->matsolve)(ksp, vB, vX);
1024:         if (ksp->viewFinalRes) {
1025:           KSPViewFinalMatResidual_Internal(ksp, vB, vX, ksp->viewerFinalRes, ksp->formatFinalRes, n2);
1026:         }
1027:         if (ksp->viewReason) {
1028:           KSPReasonView(ksp, ksp->viewerReason);
1029:         }
1030:         MatDenseRestoreSubMatrix(B, &vB);
1031:         MatDenseRestoreSubMatrix(X, &vX);
1032:       }
1033:     }
1034:     if (ksp->view) {
1035:       KSPView(ksp, ksp->viewer);
1036:     }
1037:     PetscLogEventEnd(KSP_MatSolve, ksp, B, X, 0);
1038:   } else {
1039:     PetscInfo1(ksp, "KSP type %s solving column by column\n", ((PetscObject)ksp)->type_name);
1040:     for (n2 = 0; n2 < N2; ++n2) {
1041:       MatDenseGetColumnVecRead(B, n2, &cb);
1042:       MatDenseGetColumnVecWrite(X, n2, &cx);
1043:       KSPSolve(ksp, cb, cx);
1044:       MatDenseRestoreColumnVecWrite(X, n2, &cx);
1045:       MatDenseRestoreColumnVecRead(B, n2, &cb);
1046:     }
1047:   }
1048:   return(0);
1049: }

1051: /*@
1052:      KSPSetMatSolveBlockSize - Sets the maximum number of columns treated simultaneously in KSPMatSolve().

1054:     Logically collective

1056:    Input Parameters:
1057: +     ksp - iterative context
1058: -     bs - block size

1060:    Level: advanced

1062: .seealso:  KSPMatSolve(), KSPGetMatSolveBlockSize(), -mat_mumps_icntl_27, -matmatmult_Bbn
1063: @*/
1064: PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt bs)
1065: {

1071:   PetscTryMethod(ksp, "KSPSetMatSolveBlockSize_C", (KSP, PetscInt), (ksp, bs));
1072:   return(0);
1073: }

1075: /*@
1076:      KSPGetMatSolveBlockSize - Gets the maximum number of columns treated simultaneously in KSPMatSolve().

1078:    Input Parameter:
1079: .     ksp - iterative context

1081:    Output Parameter:
1082: .     bs - block size

1084:    Level: advanced

1086: .seealso:  KSPMatSolve(), KSPSetMatSolveBlockSize(), -mat_mumps_icntl_27, -matmatmult_Bbn
1087: @*/
1088: PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *bs)
1089: {

1094:   *bs = PETSC_DECIDE;
1095:   PetscTryMethod(ksp, "KSPGetMatSolveBlockSize_C", (KSP, PetscInt*), (ksp, bs));
1096:   return(0);
1097: }

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

1102:    Collective on ksp

1104:    Input Parameter:
1105: .  ksp - iterative context obtained from KSPCreate()

1107:    Level: beginner

1109: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSPSetFromOptions(), KSP
1110: @*/
1111: PetscErrorCode  KSPResetViewers(KSP ksp)
1112: {

1117:   if (!ksp) return(0);
1118:   PetscViewerDestroy(&ksp->viewer);
1119:   PetscViewerDestroy(&ksp->viewerPre);
1120:   PetscViewerDestroy(&ksp->viewerReason);
1121:   PetscViewerDestroy(&ksp->viewerMat);
1122:   PetscViewerDestroy(&ksp->viewerPMat);
1123:   PetscViewerDestroy(&ksp->viewerRhs);
1124:   PetscViewerDestroy(&ksp->viewerSol);
1125:   PetscViewerDestroy(&ksp->viewerMatExp);
1126:   PetscViewerDestroy(&ksp->viewerEV);
1127:   PetscViewerDestroy(&ksp->viewerSV);
1128:   PetscViewerDestroy(&ksp->viewerEVExp);
1129:   PetscViewerDestroy(&ksp->viewerFinalRes);
1130:   PetscViewerDestroy(&ksp->viewerPOpExp);
1131:   PetscViewerDestroy(&ksp->viewerDScale);
1132:   ksp->view         = PETSC_FALSE;
1133:   ksp->viewPre      = PETSC_FALSE;
1134:   ksp->viewReason   = PETSC_FALSE;
1135:   ksp->viewMat      = PETSC_FALSE;
1136:   ksp->viewPMat     = PETSC_FALSE;
1137:   ksp->viewRhs      = PETSC_FALSE;
1138:   ksp->viewSol      = PETSC_FALSE;
1139:   ksp->viewMatExp   = PETSC_FALSE;
1140:   ksp->viewEV       = PETSC_FALSE;
1141:   ksp->viewSV       = PETSC_FALSE;
1142:   ksp->viewEVExp    = PETSC_FALSE;
1143:   ksp->viewFinalRes = PETSC_FALSE;
1144:   ksp->viewPOpExp   = PETSC_FALSE;
1145:   ksp->viewDScale   = PETSC_FALSE;
1146:   return(0);
1147: }

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

1152:    Collective on ksp

1154:    Input Parameter:
1155: .  ksp - iterative context obtained from KSPCreate()

1157:    Level: beginner

1159: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1160: @*/
1161: PetscErrorCode  KSPReset(KSP ksp)
1162: {

1167:   if (!ksp) return(0);
1168:   if (ksp->ops->reset) {
1169:     (*ksp->ops->reset)(ksp);
1170:   }
1171:   if (ksp->pc) {PCReset(ksp->pc);}
1172:   if (ksp->guess) {
1173:     KSPGuess guess = ksp->guess;
1174:     if (guess->ops->reset) { (*guess->ops->reset)(guess); }
1175:   }
1176:   VecDestroyVecs(ksp->nwork,&ksp->work);
1177:   VecDestroy(&ksp->vec_rhs);
1178:   VecDestroy(&ksp->vec_sol);
1179:   VecDestroy(&ksp->diagonal);
1180:   VecDestroy(&ksp->truediagonal);

1182:   KSPResetViewers(ksp);

1184:   ksp->setupstage = KSP_SETUP_NEW;
1185:   return(0);
1186: }

1188: /*@
1189:    KSPDestroy - Destroys KSP context.

1191:    Collective on ksp

1193:    Input Parameter:
1194: .  ksp - iterative context obtained from KSPCreate()

1196:    Level: beginner

1198: .seealso: KSPCreate(), KSPSetUp(), KSPSolve(), KSP
1199: @*/
1200: PetscErrorCode  KSPDestroy(KSP *ksp)
1201: {
1203:   PC             pc;

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

1210:   PetscObjectSAWsViewOff((PetscObject)*ksp);

1212:   /*
1213:    Avoid a cascading call to PCReset(ksp->pc) from the following call:
1214:    PCReset() shouldn't be called from KSPDestroy() as it is unprotected by pc's
1215:    refcount (and may be shared, e.g., by other ksps).
1216:    */
1217:   pc         = (*ksp)->pc;
1218:   (*ksp)->pc = NULL;
1219:   KSPReset((*ksp));
1220:   (*ksp)->pc = pc;
1221:   if ((*ksp)->ops->destroy) {(*(*ksp)->ops->destroy)(*ksp);}

1223:   KSPGuessDestroy(&(*ksp)->guess);
1224:   DMDestroy(&(*ksp)->dm);
1225:   PCDestroy(&(*ksp)->pc);
1226:   PetscFree((*ksp)->res_hist_alloc);
1227:   if ((*ksp)->convergeddestroy) {
1228:     (*(*ksp)->convergeddestroy)((*ksp)->cnvP);
1229:   }
1230:   KSPMonitorCancel((*ksp));
1231:   PetscViewerDestroy(&(*ksp)->eigviewer);
1232:   PetscHeaderDestroy(ksp);
1233:   return(0);
1234: }

1236: /*@
1237:     KSPSetPCSide - Sets the preconditioning side.

1239:     Logically Collective on ksp

1241:     Input Parameter:
1242: .   ksp - iterative context obtained from KSPCreate()

1244:     Output Parameter:
1245: .   side - the preconditioning side, where side is one of
1246: .vb
1247:       PC_LEFT - left preconditioning (default)
1248:       PC_RIGHT - right preconditioning
1249:       PC_SYMMETRIC - symmetric preconditioning
1250: .ve

1252:     Options Database Keys:
1253: .   -ksp_pc_side <right,left,symmetric>

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

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

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

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

1266:     Level: intermediate

1268: .seealso: KSPGetPCSide(), KSPSetNormType(), KSPGetNormType(), KSP
1269: @*/
1270: PetscErrorCode  KSPSetPCSide(KSP ksp,PCSide side)
1271: {
1275:   ksp->pc_side = ksp->pc_side_set = side;
1276:   return(0);
1277: }

1279: /*@
1280:     KSPGetPCSide - Gets the preconditioning side.

1282:     Not Collective

1284:     Input Parameter:
1285: .   ksp - iterative context obtained from KSPCreate()

1287:     Output Parameter:
1288: .   side - the preconditioning side, where side is one of
1289: .vb
1290:       PC_LEFT - left preconditioning (default)
1291:       PC_RIGHT - right preconditioning
1292:       PC_SYMMETRIC - symmetric preconditioning
1293: .ve

1295:     Level: intermediate

1297: .seealso: KSPSetPCSide(), KSP
1298: @*/
1299: PetscErrorCode  KSPGetPCSide(KSP ksp,PCSide *side)
1300: {

1306:   KSPSetUpNorms_Private(ksp,PETSC_TRUE,&ksp->normtype,&ksp->pc_side);
1307:   *side = ksp->pc_side;
1308:   return(0);
1309: }

1311: /*@
1312:    KSPGetTolerances - Gets the relative, absolute, divergence, and maximum
1313:    iteration tolerances used by the default KSP convergence tests.

1315:    Not Collective

1317:    Input Parameter:
1318: .  ksp - the Krylov subspace context

1320:    Output Parameters:
1321: +  rtol - the relative convergence tolerance
1322: .  abstol - the absolute convergence tolerance
1323: .  dtol - the divergence tolerance
1324: -  maxits - maximum number of iterations

1326:    Notes:
1327:    The user can specify NULL for any parameter that is not needed.

1329:    Level: intermediate

1331:            maximum, iterations

1333: .seealso: KSPSetTolerances(), KSP
1334: @*/
1335: PetscErrorCode  KSPGetTolerances(KSP ksp,PetscReal *rtol,PetscReal *abstol,PetscReal *dtol,PetscInt *maxits)
1336: {
1339:   if (abstol) *abstol = ksp->abstol;
1340:   if (rtol) *rtol = ksp->rtol;
1341:   if (dtol) *dtol = ksp->divtol;
1342:   if (maxits) *maxits = ksp->max_it;
1343:   return(0);
1344: }

1346: /*@
1347:    KSPSetTolerances - Sets the relative, absolute, divergence, and maximum
1348:    iteration tolerances used by the default KSP convergence testers.

1350:    Logically Collective on ksp

1352:    Input Parameters:
1353: +  ksp - the Krylov subspace context
1354: .  rtol - the relative convergence tolerance, relative decrease in the (possibly preconditioned) residual norm
1355: .  abstol - the absolute convergence tolerance   absolute size of the (possibly preconditioned) residual norm
1356: .  dtol - the divergence tolerance,   amount (possibly preconditioned) residual norm can increase before KSPConvergedDefault() concludes that the method is diverging
1357: -  maxits - maximum number of iterations to use

1359:    Options Database Keys:
1360: +  -ksp_atol <abstol> - Sets abstol
1361: .  -ksp_rtol <rtol> - Sets rtol
1362: .  -ksp_divtol <dtol> - Sets dtol
1363: -  -ksp_max_it <maxits> - Sets maxits

1365:    Notes:
1366:    Use PETSC_DEFAULT to retain the default value of any of the tolerances.

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

1371:    Level: intermediate

1373:            convergence, maximum, iterations

1375: .seealso: KSPGetTolerances(), KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
1376: @*/
1377: PetscErrorCode  KSPSetTolerances(KSP ksp,PetscReal rtol,PetscReal abstol,PetscReal dtol,PetscInt maxits)
1378: {

1386:   if (rtol != PETSC_DEFAULT) {
1387:     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);
1388:     ksp->rtol = rtol;
1389:   }
1390:   if (abstol != PETSC_DEFAULT) {
1391:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
1392:     ksp->abstol = abstol;
1393:   }
1394:   if (dtol != PETSC_DEFAULT) {
1395:     if (dtol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Divergence tolerance %g must be larger than 1.0",(double)dtol);
1396:     ksp->divtol = dtol;
1397:   }
1398:   if (maxits != PETSC_DEFAULT) {
1399:     if (maxits < 0) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxits);
1400:     ksp->max_it = maxits;
1401:   }
1402:   return(0);
1403: }

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

1410:    Logically Collective on ksp

1412:    Input Parameters:
1413: +  ksp - iterative context obtained from KSPCreate()
1414: -  flg - PETSC_TRUE indicates the guess is non-zero, PETSC_FALSE indicates the guess is zero

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

1419:    Level: beginner

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

1424: .seealso: KSPGetInitialGuessNonzero(), KSPSetGuessType(), KSPGuessType, KSP
1425: @*/
1426: PetscErrorCode  KSPSetInitialGuessNonzero(KSP ksp,PetscBool flg)
1427: {
1431:   ksp->guess_zero = (PetscBool) !(int)flg;
1432:   return(0);
1433: }

1435: /*@
1436:    KSPGetInitialGuessNonzero - Determines whether the KSP solver is using
1437:    a zero initial guess.

1439:    Not Collective

1441:    Input Parameter:
1442: .  ksp - iterative context obtained from KSPCreate()

1444:    Output Parameter:
1445: .  flag - PETSC_TRUE if guess is nonzero, else PETSC_FALSE

1447:    Level: intermediate

1449: .seealso: KSPSetInitialGuessNonzero(), KSP
1450: @*/
1451: PetscErrorCode  KSPGetInitialGuessNonzero(KSP ksp,PetscBool  *flag)
1452: {
1456:   if (ksp->guess_zero) *flag = PETSC_FALSE;
1457:   else *flag = PETSC_TRUE;
1458:   return(0);
1459: }

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

1464:    Logically Collective on ksp

1466:    Input Parameters:
1467: +  ksp - iterative context obtained from KSPCreate()
1468: -  flg - PETSC_TRUE indicates you want the error generated

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

1473:    Level: intermediate

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


1480: .seealso: KSPGetErrorIfNotConverged(), KSP
1481: @*/
1482: PetscErrorCode  KSPSetErrorIfNotConverged(KSP ksp,PetscBool flg)
1483: {
1487:   ksp->errorifnotconverged = flg;
1488:   return(0);
1489: }

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

1494:    Not Collective

1496:    Input Parameter:
1497: .  ksp - iterative context obtained from KSPCreate()

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

1502:    Level: intermediate

1504: .seealso: KSPSetErrorIfNotConverged(), KSP
1505: @*/
1506: PetscErrorCode  KSPGetErrorIfNotConverged(KSP ksp,PetscBool  *flag)
1507: {
1511:   *flag = ksp->errorifnotconverged;
1512:   return(0);
1513: }

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

1518:    Logically Collective on ksp

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

1524:    Level: advanced

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

1528: .seealso: KSPGetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1529: @*/
1530: PetscErrorCode  KSPSetInitialGuessKnoll(KSP ksp,PetscBool flg)
1531: {
1535:   ksp->guess_knoll = flg;
1536:   return(0);
1537: }

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

1543:    Not Collective

1545:    Input Parameter:
1546: .  ksp - iterative context obtained from KSPCreate()

1548:    Output Parameter:
1549: .  flag - PETSC_TRUE if using Knoll trick, else PETSC_FALSE

1551:    Level: advanced

1553: .seealso: KSPSetInitialGuessKnoll(), KSPSetInitialGuessNonzero(), KSPGetInitialGuessNonzero(), KSP
1554: @*/
1555: PetscErrorCode  KSPGetInitialGuessKnoll(KSP ksp,PetscBool  *flag)
1556: {
1560:   *flag = ksp->guess_knoll;
1561:   return(0);
1562: }

1564: /*@
1565:    KSPGetComputeSingularValues - Gets the flag indicating whether the extreme singular
1566:    values will be calculated via a Lanczos or Arnoldi process as the linear
1567:    system is solved.

1569:    Not Collective

1571:    Input Parameter:
1572: .  ksp - iterative context obtained from KSPCreate()

1574:    Output Parameter:
1575: .  flg - PETSC_TRUE or PETSC_FALSE

1577:    Options Database Key:
1578: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1580:    Notes:
1581:    Currently this option is not valid for all iterative methods.

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

1587:    Level: advanced

1589: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1590: @*/
1591: PetscErrorCode  KSPGetComputeSingularValues(KSP ksp,PetscBool  *flg)
1592: {
1596:   *flg = ksp->calc_sings;
1597:   return(0);
1598: }

1600: /*@
1601:    KSPSetComputeSingularValues - Sets a flag so that the extreme singular
1602:    values will be calculated via a Lanczos or Arnoldi process as the linear
1603:    system is solved.

1605:    Logically Collective on ksp

1607:    Input Parameters:
1608: +  ksp - iterative context obtained from KSPCreate()
1609: -  flg - PETSC_TRUE or PETSC_FALSE

1611:    Options Database Key:
1612: .  -ksp_monitor_singular_value - Activates KSPSetComputeSingularValues()

1614:    Notes:
1615:    Currently this option is not valid for all iterative methods.

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

1621:    Level: advanced

1623: .seealso: KSPComputeExtremeSingularValues(), KSPMonitorSingularValue(), KSP
1624: @*/
1625: PetscErrorCode  KSPSetComputeSingularValues(KSP ksp,PetscBool flg)
1626: {
1630:   ksp->calc_sings = flg;
1631:   return(0);
1632: }

1634: /*@
1635:    KSPGetComputeEigenvalues - Gets the flag indicating that the extreme eigenvalues
1636:    values will be calculated via a Lanczos or Arnoldi process as the linear
1637:    system is solved.

1639:    Not Collective

1641:    Input Parameter:
1642: .  ksp - iterative context obtained from KSPCreate()

1644:    Output Parameter:
1645: .  flg - PETSC_TRUE or PETSC_FALSE

1647:    Notes:
1648:    Currently this option is not valid for all iterative methods.

1650:    Level: advanced

1652: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1653: @*/
1654: PetscErrorCode  KSPGetComputeEigenvalues(KSP ksp,PetscBool  *flg)
1655: {
1659:   *flg = ksp->calc_sings;
1660:   return(0);
1661: }

1663: /*@
1664:    KSPSetComputeEigenvalues - Sets a flag so that the extreme eigenvalues
1665:    values will be calculated via a Lanczos or Arnoldi process as the linear
1666:    system is solved.

1668:    Logically Collective on ksp

1670:    Input Parameters:
1671: +  ksp - iterative context obtained from KSPCreate()
1672: -  flg - PETSC_TRUE or PETSC_FALSE

1674:    Notes:
1675:    Currently this option is not valid for all iterative methods.

1677:    Level: advanced

1679: .seealso: KSPComputeEigenvalues(), KSPComputeEigenvaluesExplicitly(), KSP
1680: @*/
1681: PetscErrorCode  KSPSetComputeEigenvalues(KSP ksp,PetscBool flg)
1682: {
1686:   ksp->calc_sings = flg;
1687:   return(0);
1688: }

1690: /*@
1691:    KSPSetComputeRitz - Sets a flag so that the Ritz or harmonic Ritz pairs
1692:    will be calculated via a Lanczos or Arnoldi process as the linear
1693:    system is solved.

1695:    Logically Collective on ksp

1697:    Input Parameters:
1698: +  ksp - iterative context obtained from KSPCreate()
1699: -  flg - PETSC_TRUE or PETSC_FALSE

1701:    Notes:
1702:    Currently this option is only valid for the GMRES method.

1704:    Level: advanced

1706: .seealso: KSPComputeRitz(), KSP
1707: @*/
1708: PetscErrorCode  KSPSetComputeRitz(KSP ksp, PetscBool flg)
1709: {
1713:   ksp->calc_ritz = flg;
1714:   return(0);
1715: }

1717: /*@
1718:    KSPGetRhs - Gets the right-hand-side vector for the linear system to
1719:    be solved.

1721:    Not Collective

1723:    Input Parameter:
1724: .  ksp - iterative context obtained from KSPCreate()

1726:    Output Parameter:
1727: .  r - right-hand-side vector

1729:    Level: developer

1731: .seealso: KSPGetSolution(), KSPSolve(), KSP
1732: @*/
1733: PetscErrorCode  KSPGetRhs(KSP ksp,Vec *r)
1734: {
1738:   *r = ksp->vec_rhs;
1739:   return(0);
1740: }

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

1747:    Not Collective

1749:    Input Parameters:
1750: .  ksp - iterative context obtained from KSPCreate()

1752:    Output Parameters:
1753: .  v - solution vector

1755:    Level: developer

1757: .seealso: KSPGetRhs(),  KSPBuildSolution(), KSPSolve(), KSP
1758: @*/
1759: PetscErrorCode  KSPGetSolution(KSP ksp,Vec *v)
1760: {
1764:   *v = ksp->vec_sol;
1765:   return(0);
1766: }

1768: /*@
1769:    KSPSetPC - Sets the preconditioner to be used to calculate the
1770:    Section 1.5 Writing Application Codes with PETSc of the preconditioner on a vector.

1772:    Collective on ksp

1774:    Input Parameters:
1775: +  ksp - iterative context obtained from KSPCreate()
1776: -  pc   - the preconditioner object (can be NULL)

1778:    Notes:
1779:    Use KSPGetPC() to retrieve the preconditioner context.

1781:    Level: developer

1783: .seealso: KSPGetPC(), KSP
1784: @*/
1785: PetscErrorCode  KSPSetPC(KSP ksp,PC pc)
1786: {

1791:   if (pc) {
1794:   }
1795:   PetscObjectReference((PetscObject)pc);
1796:   PCDestroy(&ksp->pc);
1797:   ksp->pc = pc;
1798:   PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1799:   return(0);
1800: }

1802: /*@
1803:    KSPGetPC - Returns a pointer to the preconditioner context
1804:    set with KSPSetPC().

1806:    Not Collective

1808:    Input Parameters:
1809: .  ksp - iterative context obtained from KSPCreate()

1811:    Output Parameter:
1812: .  pc - preconditioner context

1814:    Level: developer

1816: .seealso: KSPSetPC(), KSP
1817: @*/
1818: PetscErrorCode  KSPGetPC(KSP ksp,PC *pc)
1819: {

1825:   if (!ksp->pc) {
1826:     PCCreate(PetscObjectComm((PetscObject)ksp),&ksp->pc);
1827:     PetscObjectIncrementTabLevel((PetscObject)ksp->pc,(PetscObject)ksp,0);
1828:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)ksp->pc);
1829:     PetscObjectSetOptions((PetscObject)ksp->pc,((PetscObject)ksp)->options);
1830:   }
1831:   *pc = ksp->pc;
1832:   return(0);
1833: }

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

1838:    Collective on ksp

1840:    Input Parameters:
1841: +  ksp - iterative context obtained from KSPCreate()
1842: .  it - iteration number
1843: -  rnorm - relative norm of the residual

1845:    Notes:
1846:    This routine is called by the KSP implementations.
1847:    It does not typically need to be called by the user.

1849:    Level: developer

1851: .seealso: KSPMonitorSet()
1852: @*/
1853: PetscErrorCode KSPMonitor(KSP ksp,PetscInt it,PetscReal rnorm)
1854: {
1855:   PetscInt       i, n = ksp->numbermonitors;

1859:   for (i=0; i<n; i++) {
1860:     (*ksp->monitor[i])(ksp,it,rnorm,ksp->monitorcontext[i]);
1861:   }
1862:   return(0);
1863: }

1865: /*@C
1866:    KSPMonitorSet - Sets an ADDITIONAL function to be called at every iteration to monitor
1867:    the residual/error etc.

1869:    Logically Collective on ksp

1871:    Input Parameters:
1872: +  ksp - iterative context obtained from KSPCreate()
1873: .  monitor - pointer to function (if this is NULL, it turns off monitoring
1874: .  mctx    - [optional] context for private data for the
1875:              monitor routine (use NULL if no context is desired)
1876: -  monitordestroy - [optional] routine that frees monitor context
1877:           (may be NULL)

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

1882: +  ksp - iterative context obtained from KSPCreate()
1883: .  it - iteration number
1884: .  rnorm - (estimated) 2-norm of (preconditioned) residual
1885: -  mctx  - optional monitoring context, as set by KSPMonitorSet()

1887:    Options Database Keys:
1888: +    -ksp_monitor        - sets KSPMonitorDefault()
1889: .    -ksp_monitor_true_residual    - sets KSPMonitorTrueResidualNorm()
1890: .    -ksp_monitor_max    - sets KSPMonitorTrueResidualMaxNorm()
1891: .    -ksp_monitor_lg_residualnorm    - sets line graph monitor,
1892:                            uses KSPMonitorLGResidualNormCreate()
1893: .    -ksp_monitor_lg_true_residualnorm   - sets line graph monitor,
1894:                            uses KSPMonitorLGResidualNormCreate()
1895: .    -ksp_monitor_singular_value    - sets KSPMonitorSingularValue()
1896: -    -ksp_monitor_cancel - cancels all monitors that have
1897:                           been hardwired into a code by
1898:                           calls to KSPMonitorSet(), but
1899:                           does not cancel those set via
1900:                           the options database.

1902:    Notes:
1903:    The default is to do nothing.  To print the residual, or preconditioned
1904:    residual if KSPSetNormType(ksp,KSP_NORM_PRECONDITIONED) was called, use
1905:    KSPMonitorDefault() as the monitoring routine, with a ASCII viewer as the
1906:    context.

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

1912:    Fortran Notes:
1913:     Only a single monitor function can be set for each KSP object

1915:    Level: beginner

1917: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorCancel(), KSP
1918: @*/
1919: PetscErrorCode  KSPMonitorSet(KSP ksp,PetscErrorCode (*monitor)(KSP,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void**))
1920: {
1921:   PetscInt       i;
1923:   PetscBool      identical;

1927:   for (i=0; i<ksp->numbermonitors;i++) {
1928:     PetscMonitorCompare((PetscErrorCode (*)(void))monitor,mctx,monitordestroy,(PetscErrorCode (*)(void))ksp->monitor[i],ksp->monitorcontext[i],ksp->monitordestroy[i],&identical);
1929:     if (identical) return(0);
1930:   }
1931:   if (ksp->numbermonitors >= MAXKSPMONITORS) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_OUTOFRANGE,"Too many KSP monitors set");
1932:   ksp->monitor[ksp->numbermonitors]          = monitor;
1933:   ksp->monitordestroy[ksp->numbermonitors]   = monitordestroy;
1934:   ksp->monitorcontext[ksp->numbermonitors++] = (void*)mctx;
1935:   return(0);
1936: }

1938: /*@
1939:    KSPMonitorCancel - Clears all monitors for a KSP object.

1941:    Logically Collective on ksp

1943:    Input Parameters:
1944: .  ksp - iterative context obtained from KSPCreate()

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

1951:    Level: intermediate

1953: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSPMonitorSet(), KSP
1954: @*/
1955: PetscErrorCode  KSPMonitorCancel(KSP ksp)
1956: {
1958:   PetscInt       i;

1962:   for (i=0; i<ksp->numbermonitors; i++) {
1963:     if (ksp->monitordestroy[i]) {
1964:       (*ksp->monitordestroy[i])(&ksp->monitorcontext[i]);
1965:     }
1966:   }
1967:   ksp->numbermonitors = 0;
1968:   return(0);
1969: }

1971: /*@C
1972:    KSPGetMonitorContext - Gets the monitoring context, as set by
1973:    KSPMonitorSet() for the FIRST monitor only.

1975:    Not Collective

1977:    Input Parameter:
1978: .  ksp - iterative context obtained from KSPCreate()

1980:    Output Parameter:
1981: .  ctx - monitoring context

1983:    Level: intermediate

1985: .seealso: KSPMonitorDefault(), KSPMonitorLGResidualNormCreate(), KSP
1986: @*/
1987: PetscErrorCode  KSPGetMonitorContext(KSP ksp,void **ctx)
1988: {
1991:   *ctx =      (ksp->monitorcontext[0]);
1992:   return(0);
1993: }

1995: /*@
1996:    KSPSetResidualHistory - Sets the array used to hold the residual history.
1997:    If set, this array will contain the residual norms computed at each
1998:    iteration of the solver.

2000:    Not Collective

2002:    Input Parameters:
2003: +  ksp - iterative context obtained from KSPCreate()
2004: .  a   - array to hold history
2005: .  na  - size of a
2006: -  reset - PETSC_TRUE indicates the history counter is reset to zero
2007:            for each new linear solve

2009:    Level: advanced

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

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

2018: .seealso: KSPGetResidualHistory(), KSP

2020: @*/
2021: PetscErrorCode  KSPSetResidualHistory(KSP ksp,PetscReal a[],PetscInt na,PetscBool reset)
2022: {


2028:   PetscFree(ksp->res_hist_alloc);
2029:   if (na != PETSC_DECIDE && na != PETSC_DEFAULT && a) {
2030:     ksp->res_hist     = a;
2031:     ksp->res_hist_max = na;
2032:   } else {
2033:     if (na != PETSC_DECIDE && na != PETSC_DEFAULT) ksp->res_hist_max = na;
2034:     else                                           ksp->res_hist_max = 10000; /* like default ksp->max_it */
2035:     PetscCalloc1(ksp->res_hist_max,&ksp->res_hist_alloc);

2037:     ksp->res_hist = ksp->res_hist_alloc;
2038:   }
2039:   ksp->res_hist_len   = 0;
2040:   ksp->res_hist_reset = reset;
2041:   return(0);
2042: }

2044: /*@C
2045:    KSPGetResidualHistory - Gets the array used to hold the residual history
2046:    and the number of residuals it contains.

2048:    Not Collective

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

2053:    Output Parameters:
2054: +  a   - pointer to array to hold history (or NULL)
2055: -  na  - number of used entries in a (or NULL)

2057:    Level: advanced

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

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

2068: .seealso: KSPGetResidualHistory(), KSP

2070: @*/
2071: PetscErrorCode  KSPGetResidualHistory(KSP ksp,PetscReal *a[],PetscInt *na)
2072: {
2075:   if (a) *a = ksp->res_hist;
2076:   if (na) *na = ksp->res_hist_len;
2077:   return(0);
2078: }

2080: /*@C
2081:    KSPSetConvergenceTest - Sets the function to be used to determine
2082:    convergence.

2084:    Logically Collective on ksp

2086:    Input Parameters:
2087: +  ksp - iterative context obtained from KSPCreate()
2088: .  converge - pointer to the function
2089: .  cctx    - context for private data for the convergence routine (may be null)
2090: -  destroy - a routine for destroying the context (may be null)

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

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


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

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

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

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

2116:    Level: advanced

2118: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPGetConvergenceTest(), KSPGetAndClearConvergenceTest()
2119: @*/
2120: PetscErrorCode  KSPSetConvergenceTest(KSP ksp,PetscErrorCode (*converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*))
2121: {

2126:   if (ksp->convergeddestroy) {
2127:     (*ksp->convergeddestroy)(ksp->cnvP);
2128:   }
2129:   ksp->converged        = converge;
2130:   ksp->convergeddestroy = destroy;
2131:   ksp->cnvP             = (void*)cctx;
2132:   return(0);
2133: }

2135: /*@C
2136:    KSPGetConvergenceTest - Gets the function to be used to determine
2137:    convergence.

2139:    Logically Collective on ksp

2141:    Input Parameter:
2142: .   ksp - iterative context obtained from KSPCreate()

2144:    Output Parameter:
2145: +  converge - pointer to convergence test function
2146: .  cctx    - context for private data for the convergence routine (may be null)
2147: -  destroy - a routine for destroying the context (may be null)

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

2152: +  ksp - iterative context obtained from KSPCreate()
2153: .  it - iteration number
2154: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2155: .  reason - the reason why it has converged or diverged
2156: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2158:    Level: advanced

2160: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetAndClearConvergenceTest()
2161: @*/
2162: PetscErrorCode  KSPGetConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2163: {
2166:   if (converge) *converge = ksp->converged;
2167:   if (destroy)  *destroy  = ksp->convergeddestroy;
2168:   if (cctx)     *cctx     = ksp->cnvP;
2169:   return(0);
2170: }

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

2175:    Logically Collective on ksp

2177:    Input Parameter:
2178: .   ksp - iterative context obtained from KSPCreate()

2180:    Output Parameter:
2181: +  converge - pointer to convergence test function
2182: .  cctx    - context for private data for the convergence routine
2183: -  destroy - a routine for destroying the context

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

2188: +  ksp - iterative context obtained from KSPCreate()
2189: .  it - iteration number
2190: .  rnorm - (estimated) 2-norm of (preconditioned) residual
2191: .  reason - the reason why it has converged or diverged
2192: -  cctx  - optional convergence context, as set by KSPSetConvergenceTest()

2194:    Level: advanced

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

2200: .seealso: KSPConvergedDefault(), KSPGetConvergenceContext(), KSPSetTolerances(), KSP, KSPSetConvergenceTest(), KSPGetConvergenceTest()
2201: @*/
2202: PetscErrorCode  KSPGetAndClearConvergenceTest(KSP ksp,PetscErrorCode (**converge)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void **cctx,PetscErrorCode (**destroy)(void*))
2203: {
2206:   *converge             = ksp->converged;
2207:   *destroy              = ksp->convergeddestroy;
2208:   *cctx                 = ksp->cnvP;
2209:   ksp->converged        = NULL;
2210:   ksp->cnvP             = NULL;
2211:   ksp->convergeddestroy = NULL;
2212:   return(0);
2213: }

2215: /*@C
2216:    KSPGetConvergenceContext - Gets the convergence context set with
2217:    KSPSetConvergenceTest().

2219:    Not Collective

2221:    Input Parameter:
2222: .  ksp - iterative context obtained from KSPCreate()

2224:    Output Parameter:
2225: .  ctx - monitoring context

2227:    Level: advanced

2229: .seealso: KSPConvergedDefault(), KSPSetConvergenceTest(), KSP
2230: @*/
2231: PetscErrorCode  KSPGetConvergenceContext(KSP ksp,void **ctx)
2232: {
2235:   *ctx = ksp->cnvP;
2236:   return(0);
2237: }

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

2243:    Collective on ksp

2245:    Input Parameter:
2246: .  ctx - iterative context obtained from KSPCreate()

2248:    Output Parameter:
2249:    Provide exactly one of
2250: +  v - location to stash solution.
2251: -  V - the solution is returned in this location. This vector is created
2252:        internally. This vector should NOT be destroyed by the user with
2253:        VecDestroy().

2255:    Notes:
2256:    This routine can be used in one of two ways
2257: .vb
2258:       KSPBuildSolution(ksp,NULL,&V);
2259:    or
2260:       KSPBuildSolution(ksp,v,NULL); or KSPBuildSolution(ksp,v,&v);
2261: .ve
2262:    In the first case an internal vector is allocated to store the solution
2263:    (the user cannot destroy this vector). In the second case the solution
2264:    is generated in the vector that the user provides. Note that for certain
2265:    methods, such as KSPCG, the second case requires a copy of the solution,
2266:    while in the first case the call is essentially free since it simply
2267:    returns the vector where the solution already is stored. For some methods
2268:    like GMRES this is a reasonably expensive operation and should only be
2269:    used in truly needed.

2271:    Level: advanced

2273: .seealso: KSPGetSolution(), KSPBuildResidual(), KSP
2274: @*/
2275: PetscErrorCode  KSPBuildSolution(KSP ksp,Vec v,Vec *V)
2276: {

2281:   if (!V && !v) SETERRQ(PetscObjectComm((PetscObject)ksp),PETSC_ERR_ARG_WRONG,"Must provide either v or V");
2282:   if (!V) V = &v;
2283:   (*ksp->ops->buildsolution)(ksp,v,V);
2284:   return(0);
2285: }

2287: /*@C
2288:    KSPBuildResidual - Builds the residual in a vector provided.

2290:    Collective on ksp

2292:    Input Parameter:
2293: .  ksp - iterative context obtained from KSPCreate()

2295:    Output Parameters:
2296: +  v - optional location to stash residual.  If v is not provided,
2297:        then a location is generated.
2298: .  t - work vector.  If not provided then one is generated.
2299: -  V - the residual

2301:    Notes:
2302:    Regardless of whether or not v is provided, the residual is
2303:    returned in V.

2305:    Level: advanced

2307: .seealso: KSPBuildSolution()
2308: @*/
2309: PetscErrorCode  KSPBuildResidual(KSP ksp,Vec t,Vec v,Vec *V)
2310: {
2312:   PetscBool      flag = PETSC_FALSE;
2313:   Vec            w    = v,tt = t;

2317:   if (!w) {
2318:     VecDuplicate(ksp->vec_rhs,&w);
2319:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)w);
2320:   }
2321:   if (!tt) {
2322:     VecDuplicate(ksp->vec_sol,&tt); flag = PETSC_TRUE;
2323:     PetscLogObjectParent((PetscObject)ksp,(PetscObject)tt);
2324:   }
2325:   (*ksp->ops->buildresidual)(ksp,tt,w,V);
2326:   if (flag) {VecDestroy(&tt);}
2327:   return(0);
2328: }

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

2334:    Logically Collective on ksp

2336:    Input Parameter:
2337: +  ksp - the KSP context
2338: -  scale - PETSC_TRUE or PETSC_FALSE

2340:    Options Database Key:
2341: +   -ksp_diagonal_scale -
2342: -   -ksp_diagonal_scale_fix - scale the matrix back AFTER the solve


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

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

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

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

2360:    Level: intermediate

2362: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2363: @*/
2364: PetscErrorCode  KSPSetDiagonalScale(KSP ksp,PetscBool scale)
2365: {
2369:   ksp->dscale = scale;
2370:   return(0);
2371: }

2373: /*@
2374:    KSPGetDiagonalScale - Checks if KSP solver scales the matrix and
2375:                           right hand side

2377:    Not Collective

2379:    Input Parameter:
2380: .  ksp - the KSP context

2382:    Output Parameter:
2383: .  scale - PETSC_TRUE or PETSC_FALSE

2385:    Notes:
2386:     BE CAREFUL with this routine: it actually scales the matrix and right
2387:     hand side that define the system. After the system is solved the matrix
2388:     and right hand side remain scaled  unless you use KSPSetDiagonalScaleFix()

2390:    Level: intermediate

2392: .seealso: KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2393: @*/
2394: PetscErrorCode  KSPGetDiagonalScale(KSP ksp,PetscBool  *scale)
2395: {
2399:   *scale = ksp->dscale;
2400:   return(0);
2401: }

2403: /*@
2404:    KSPSetDiagonalScaleFix - Tells KSP to diagonally scale the system
2405:      back after solving.

2407:    Logically Collective on ksp

2409:    Input Parameter:
2410: +  ksp - the KSP context
2411: -  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2412:          rescale (default)

2414:    Notes:
2415:      Must be called after KSPSetDiagonalScale()

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

2422:    Level: intermediate

2424: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPGetDiagonalScaleFix(), KSP
2425: @*/
2426: PetscErrorCode  KSPSetDiagonalScaleFix(KSP ksp,PetscBool fix)
2427: {
2431:   ksp->dscalefix = fix;
2432:   return(0);
2433: }

2435: /*@
2436:    KSPGetDiagonalScaleFix - Determines if KSP diagonally scales the system
2437:      back after solving.

2439:    Not Collective

2441:    Input Parameter:
2442: .  ksp - the KSP context

2444:    Output Parameter:
2445: .  fix - PETSC_TRUE to scale back after the system solve, PETSC_FALSE to not
2446:          rescale (default)

2448:    Notes:
2449:      Must be called after KSPSetDiagonalScale()

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

2456:    Level: intermediate

2458: .seealso: KSPGetDiagonalScale(), KSPSetDiagonalScale(), KSPSetDiagonalScaleFix(), KSP
2459: @*/
2460: PetscErrorCode  KSPGetDiagonalScaleFix(KSP ksp,PetscBool  *fix)
2461: {
2465:   *fix = ksp->dscalefix;
2466:   return(0);
2467: }

2469: /*@C
2470:    KSPSetComputeOperators - set routine to compute the linear operators

2472:    Logically Collective

2474:    Input Arguments:
2475: +  ksp - the KSP context
2476: .  func - function to compute the operators
2477: -  ctx - optional context

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

2482: +  ksp - the KSP context
2483: .  A - the linear operator
2484: .  B - preconditioning matrix
2485: -  ctx - optional user-provided context

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

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

2493:    Level: beginner

2495: .seealso: KSPSetOperators(), KSPSetComputeRHS(), DMKSPSetComputeOperators(), KSPSetComputeInitialGuess()
2496: @*/
2497: PetscErrorCode KSPSetComputeOperators(KSP ksp,PetscErrorCode (*func)(KSP,Mat,Mat,void*),void *ctx)
2498: {
2500:   DM             dm;

2504:   KSPGetDM(ksp,&dm);
2505:   DMKSPSetComputeOperators(dm,func,ctx);
2506:   if (ksp->setupstage == KSP_SETUP_NEWRHS) ksp->setupstage = KSP_SETUP_NEWMATRIX;
2507:   return(0);
2508: }

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

2513:    Logically Collective

2515:    Input Arguments:
2516: +  ksp - the KSP context
2517: .  func - function to compute the right hand side
2518: -  ctx - optional context

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

2523: +  ksp - the KSP context
2524: .  b - right hand side of linear system
2525: -  ctx - optional user-provided context

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

2530:    Level: beginner

2532: .seealso: KSPSolve(), DMKSPSetComputeRHS(), KSPSetComputeOperators()
2533: @*/
2534: PetscErrorCode KSPSetComputeRHS(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2535: {
2537:   DM             dm;

2541:   KSPGetDM(ksp,&dm);
2542:   DMKSPSetComputeRHS(dm,func,ctx);
2543:   return(0);
2544: }

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

2549:    Logically Collective

2551:    Input Arguments:
2552: +  ksp - the KSP context
2553: .  func - function to compute the initial guess
2554: -  ctx - optional context

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

2559: +  ksp - the KSP context
2560: .  x - solution vector
2561: -  ctx - optional user-provided context

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

2566:    Level: beginner

2568: .seealso: KSPSolve(), KSPSetComputeRHS(), KSPSetComputeOperators(), DMKSPSetComputeInitialGuess()
2569: @*/
2570: PetscErrorCode KSPSetComputeInitialGuess(KSP ksp,PetscErrorCode (*func)(KSP,Vec,void*),void *ctx)
2571: {
2573:   DM             dm;

2577:   KSPGetDM(ksp,&dm);
2578:   DMKSPSetComputeInitialGuess(dm,func,ctx);
2579:   return(0);
2580: }