Actual source code: itfunc.c

petsc-master 2020-08-10
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 && format != PETSC_VIEWER_FAILED) {
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 if (ksp->reason <= 0) {
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
469: .  -ksp_converged_reason ::failed - only print reason and number of iterations when diverged

471:    Level: beginner

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

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

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

492:   Collective on ksp

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

497:   Level: intermediate

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

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

516:  #include <petscdraw.h>

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

732:   ksp->guess_zero = guess_zero;

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

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

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

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

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

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

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

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

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

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

811:    Collective on ksp

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

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

832:    Notes:

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

836:    The operator is specified with KSPSetOperators().

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

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

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


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

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

869:    Level: beginner

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

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

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

890:    Collective on ksp

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

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

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

903:    Level: developer

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

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

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

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

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

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

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

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

960:    Level: intermediate

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

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

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

1055:     Logically collective

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

1061:    Level: advanced

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

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

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

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

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

1085:    Level: advanced

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

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

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

1103:    Collective on ksp

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

1108:    Level: beginner

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

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

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

1153:    Collective on ksp

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

1158:    Level: beginner

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

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

1183:   KSPResetViewers(ksp);

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

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

1192:    Collective on ksp

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

1197:    Level: beginner

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

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

1211:   PetscObjectSAWsViewOff((PetscObject)*ksp);

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

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

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

1240:     Logically Collective on ksp

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

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

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

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

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

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

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

1267:     Level: intermediate

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

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

1283:     Not Collective

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

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

1296:     Level: intermediate

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

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

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

1316:    Not Collective

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

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

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

1330:    Level: intermediate

1332:            maximum, iterations

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

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

1351:    Logically Collective on ksp

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

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

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

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

1372:    Level: intermediate

1374:            convergence, maximum, iterations

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

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

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

1411:    Logically Collective on ksp

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

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

1420:    Level: beginner

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

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

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

1440:    Not Collective

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

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

1448:    Level: intermediate

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

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

1465:    Logically Collective on ksp

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

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

1474:    Level: intermediate

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


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

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

1495:    Not Collective

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

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

1503:    Level: intermediate

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

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

1519:    Logically Collective on ksp

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

1525:    Level: advanced

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

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

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

1544:    Not Collective

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

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

1552:    Level: advanced

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

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

1570:    Not Collective

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

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

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

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

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

1588:    Level: advanced

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

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

1606:    Logically Collective on ksp

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

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

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

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

1622:    Level: advanced

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

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

1640:    Not Collective

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

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

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

1651:    Level: advanced

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

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

1669:    Logically Collective on ksp

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

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

1678:    Level: advanced

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

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

1696:    Logically Collective on ksp

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

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

1705:    Level: advanced

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

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

1722:    Not Collective

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

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

1730:    Level: developer

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

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

1748:    Not Collective

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

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

1756:    Level: developer

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

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

1773:    Collective on ksp

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

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

1782:    Level: developer

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

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

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

1807:    Not Collective

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

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

1815:    Level: developer

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

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

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

1839:    Collective on ksp

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

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

1850:    Level: developer

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

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

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

1870:    Logically Collective on ksp

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

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

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

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

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

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

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

1916:    Level: beginner

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

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

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

1942:    Logically Collective on ksp

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

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

1952:    Level: intermediate

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

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

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

1976:    Not Collective

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

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

1984:    Level: intermediate

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

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

2001:    Not Collective

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

2010:    Level: advanced

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

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

2019: .seealso: KSPGetResidualHistory(), KSP

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


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

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

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

2049:    Not Collective

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

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

2058:    Level: advanced

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

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

2069: .seealso: KSPGetResidualHistory(), KSP

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

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

2085:    Logically Collective on ksp

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

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

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


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

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

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

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

2117:    Level: advanced

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

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

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

2140:    Logically Collective on ksp

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

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

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

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

2159:    Level: advanced

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

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

2176:    Logically Collective on ksp

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

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

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

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

2195:    Level: advanced

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

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

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

2220:    Not Collective

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

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

2228:    Level: advanced

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

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

2244:    Collective on ksp

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

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

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

2272:    Level: advanced

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

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

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

2291:    Collective on ksp

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

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

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

2306:    Level: advanced

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

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

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

2335:    Logically Collective on ksp

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

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


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

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

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

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

2361:    Level: intermediate

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

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

2378:    Not Collective

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

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

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

2391:    Level: intermediate

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

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

2408:    Logically Collective on ksp

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

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

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

2423:    Level: intermediate

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

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

2440:    Not Collective

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

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

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

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

2457:    Level: intermediate

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

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

2473:    Logically Collective

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

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

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

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

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

2494:    Level: beginner

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

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

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

2514:    Logically Collective

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

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

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

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

2531:    Level: beginner

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

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

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

2550:    Logically Collective

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

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

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

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

2567:    Level: beginner

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

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